1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LPDF_HPP 2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LPDF_HPP 52 template <
typename T_shape>
56 const int Km1 = K - 1;
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);
64 constant += 0.25 * (K * K - 1) *
LOG_PI 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));
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));
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,
85 static const char*
function(
"lkj_corr_lpdf");
87 using boost::math::tools::promote_args;
89 typename promote_args<T_y, T_shape>::type lp(0.0);
93 const unsigned int K = y.rows();
107 Eigen::Matrix<T_y, Eigen::Dynamic, 1> values =
108 y.ldlt().vectorD().array().log().matrix();
109 lp += (eta - 1.0) *
sum(values);
113 template <
typename T_y,
typename T_shape>
115 typename boost::math::tools::promote_args<T_y, T_shape>::type
117 const T_shape& eta) {
118 return lkj_corr_lpdf<false>(y, eta);
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Metaprogram structure to determine the base scalar type of a template argument.
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
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)
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.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)