Stan Math Library  2.14.0
reverse mode automatic differentiation
lkj_corr_lpdf.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LPDF_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LPDF_HPP
3 
48 
49 namespace stan {
50  namespace math {
51 
52  template <typename T_shape>
53  T_shape do_lkj_constant(const T_shape& eta, const unsigned int& K) {
54  // Lewandowski, Kurowicka, and Joe (2009) theorem 5
55  T_shape constant;
56  const int Km1 = K - 1;
57  if (eta == 1.0) {
58  // C++ integer division is appropriate in this block
59  Eigen::VectorXd numerator(Km1 / 2);
60  for (int k = 1; k <= numerator.rows(); k++)
61  numerator(k - 1) = lgamma(2.0 * k);
62  constant = sum(numerator);
63  if ((K % 2) == 1)
64  constant += 0.25 * (K * K - 1) * LOG_PI
65  - 0.25 * (Km1 * Km1) * LOG_TWO - Km1 * lgamma(0.5 * (K + 1));
66  else
67  constant += 0.25 * K * (K - 2) * LOG_PI
68  + 0.25 * (3 * K * K - 4 * K) * LOG_TWO
69  + K * lgamma(0.5 * K) - Km1 * lgamma(static_cast<double>(K));
70  } else {
71  constant = -Km1 * lgamma(eta + 0.5 * Km1);
72  for (int k = 1; k <= Km1; k++)
73  constant += 0.5 * k * LOG_PI + lgamma(eta + 0.5 * (Km1 - k));
74  }
75  return constant;
76  }
77 
78  // LKJ_Corr(y|eta) [ y correlation matrix (not covariance matrix)
79  // eta > 0; eta == 1 <-> uniform]
80  template <bool propto,
81  typename T_y, typename T_shape>
82  typename boost::math::tools::promote_args<T_y, T_shape>::type
83  lkj_corr_lpdf(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
84  const T_shape& eta) {
85  static const char* function("lkj_corr_lpdf");
86 
87  using boost::math::tools::promote_args;
88 
89  typename promote_args<T_y, T_shape>::type lp(0.0);
90  check_positive(function, "Shape parameter", eta);
91  check_corr_matrix(function, "Correlation matrix", y);
92 
93  const unsigned int K = y.rows();
94  if (K == 0)
95  return 0.0;
96 
98  lp += do_lkj_constant(eta, K);
99 
100  if ( (eta == 1.0) &&
102  return lp;
103 
105  return lp;
106 
107  Eigen::Matrix<T_y, Eigen::Dynamic, 1> values =
108  y.ldlt().vectorD().array().log().matrix();
109  lp += (eta - 1.0) * sum(values);
110  return lp;
111  }
112 
113  template <typename T_y, typename T_shape>
114  inline
115  typename boost::math::tools::promote_args<T_y, T_shape>::type
116  lkj_corr_lpdf(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
117  const T_shape& eta) {
118  return lkj_corr_lpdf<false>(y, eta);
119  }
120 
121  }
122 }
123 #endif
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition: sum.hpp:20
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Definition: is_constant.hpp:22
Metaprogram structure to determine the base scalar type of a template argument.
Definition: scalar_type.hpp:33
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition: lgamma.hpp:20
boost::math::tools::promote_args< T_y, T_shape >::type lkj_corr_lpdf(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const T_shape &eta)
const double LOG_PI
Definition: constants.hpp:167
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
void check_corr_matrix(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Check if the specified matrix is a valid correlation matrix.
const double LOG_TWO
Definition: constants.hpp:174
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)

     [ Stan Home Page ] © 2011–2016, Stan Development Team.