Stan Math Library  2.20.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 
9 
10 namespace stan {
11 namespace math {
12 
13 template <typename T_shape>
14 typename boost::math::tools::promote_args<double, T_shape>::type
15 do_lkj_constant(const T_shape& eta, const unsigned int& K) {
16  // Lewandowski, Kurowicka, and Joe (2009) theorem 5
17  typename boost::math::tools::promote_args<double, T_shape>::type constant;
18  const int Km1 = K - 1;
19  using stan::math::lgamma;
20  if (eta == 1.0) {
21  // C++ integer division is appropriate in this block
22  Eigen::VectorXd denominator(Km1 / 2);
23  for (int k = 1; k <= denominator.rows(); k++)
24  denominator(k - 1) = lgamma(2.0 * k);
25  constant = -denominator.sum();
26  if ((K % 2) == 1)
27  constant -= 0.25 * (K * K - 1) * LOG_PI - 0.25 * (Km1 * Km1) * LOG_2
28  - Km1 * lgamma(0.5 * (K + 1));
29  else
30  constant -= 0.25 * K * (K - 2) * LOG_PI
31  + 0.25 * (3 * K * K - 4 * K) * LOG_2 + K * lgamma(0.5 * K)
32  - Km1 * lgamma(static_cast<double>(K));
33  } else {
34  constant = Km1 * lgamma(eta + 0.5 * Km1);
35  for (int k = 1; k <= Km1; k++)
36  constant -= 0.5 * k * LOG_PI + lgamma(eta + 0.5 * (Km1 - k));
37  }
38  return constant;
39 }
40 
41 // LKJ_Corr(y|eta) [ y correlation matrix (not covariance matrix)
42 // eta > 0; eta == 1 <-> uniform]
43 template <bool propto, typename T_y, typename T_shape>
44 typename boost::math::tools::promote_args<T_y, T_shape>::type lkj_corr_lpdf(
45  const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
46  const T_shape& eta) {
47  static const char* function = "lkj_corr_lpdf";
48 
49  using boost::math::tools::promote_args;
50 
51  typename promote_args<T_y, T_shape>::type lp(0.0);
52  check_positive(function, "Shape parameter", eta);
53  check_corr_matrix(function, "Correlation matrix", y);
54 
55  const unsigned int K = y.rows();
56  if (K == 0)
57  return 0.0;
58 
60  lp += do_lkj_constant(eta, K);
61 
62  if ((eta == 1.0)
64  return lp;
65 
67  return lp;
68 
69  Eigen::Matrix<T_y, Eigen::Dynamic, 1> values
70  = y.ldlt().vectorD().array().log().matrix();
71  lp += (eta - 1.0) * sum(values);
72  return lp;
73 }
74 
75 template <typename T_y, typename T_shape>
76 inline typename boost::math::tools::promote_args<T_y, T_shape>::type
77 lkj_corr_lpdf(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
78  const T_shape& eta) {
79  return lkj_corr_lpdf<false>(y, eta);
80 }
81 
82 } // namespace math
83 } // namespace stan
84 #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
const double LOG_2
The natural logarithm of 2, .
Definition: constants.hpp:37
Metaprogram structure to determine the base scalar type of a template argument.
Definition: scalar_type.hpp:16
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition: lgamma.hpp:21
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:144
Extends std::true_type when instantiated with zero or more template parameters, all of which extend t...
Definition: conjunction.hpp:14
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.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
boost::math::tools::promote_args< double, T_shape >::type do_lkj_constant(const T_shape &eta, const unsigned int &K)

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