Stan Math Library  2.10.0
reverse mode automatic differentiation
lkj_corr_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_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  using stan::math::sum;
55  using stan::math::lgamma;
56 
57  // Lewandowski, Kurowicka, and Joe (2009) theorem 5
58  T_shape constant;
59  const int Km1 = K - 1;
60  if (eta == 1.0) {
61  // C++ integer division is appropriate in this block
62  Eigen::VectorXd numerator(Km1 / 2);
63  for (int k = 1; k <= numerator.rows(); k++)
64  numerator(k-1) = lgamma(2 * k);
65  constant = sum(numerator);
66  if ( (K % 2) == 1 )
67  constant += 0.25 * (K * K - 1) * LOG_PI
68  - 0.25 * (Km1 * Km1) * LOG_TWO - Km1 * lgamma((K + 1) / 2);
69  else
70  constant += 0.25 * K * (K - 2) * LOG_PI
71  + 0.25 * (3 * K * K - 4 * K) * LOG_TWO
72  + K * lgamma(K / 2) - Km1 * lgamma(K);
73  } else {
74  constant = -Km1 * lgamma(eta + 0.5 * Km1);
75  for (int k = 1; k <= Km1; k++)
76  constant += 0.5 * k * LOG_PI + lgamma(eta + 0.5 * (Km1 - k));
77  }
78  return constant;
79  }
80 
81  // LKJ_Corr(y|eta) [ y correlation matrix (not covariance matrix)
82  // eta > 0; eta == 1 <-> uniform]
83  template <bool propto,
84  typename T_y, typename T_shape>
85  typename boost::math::tools::promote_args<T_y, T_shape>::type
86  lkj_corr_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
87  const T_shape& eta) {
88  static const char* function("stan::math::lkj_corr_log");
89 
92  using stan::math::sum;
93  using boost::math::tools::promote_args;
94 
95  typename promote_args<T_y, T_shape>::type lp(0.0);
96  check_positive(function, "Shape parameter", eta);
97  check_corr_matrix(function, "Correlation matrix", y);
98 
99  const unsigned int K = y.rows();
100  if (K == 0)
101  return 0.0;
102 
104  lp += do_lkj_constant(eta, K);
105 
106  if ( (eta == 1.0) &&
108  return lp;
109 
111  return lp;
112 
113  Eigen::Matrix<T_y, Eigen::Dynamic, 1> values =
114  y.ldlt().vectorD().array().log().matrix();
115  lp += (eta - 1.0) * sum(values);
116  return lp;
117  }
118 
119  template <typename T_y, typename T_shape>
120  inline
121  typename boost::math::tools::promote_args<T_y, T_shape>::type
122  lkj_corr_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
123  const T_shape& eta) {
124  return lkj_corr_log<false>(y, eta);
125  }
126 
127  }
128 }
129 #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:34
boost::math::tools::promote_args< T_y, T_shape >::type lkj_corr_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const T_shape &eta)
fvar< T > lgamma(const fvar< T > &x)
Definition: lgamma.hpp:15
const double LOG_PI
Definition: constants.hpp:170
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
const double LOG_TWO
Definition: constants.hpp:177
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
bool check_corr_matrix(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is a valid correlation matrix.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)

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