1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_HPP
54 template <
typename T_shape>
61 const int Km1 = K - 1;
64 Eigen::VectorXd numerator(Km1 / 2);
65 for (
int k = 1; k <= numerator.rows(); k++)
66 numerator(k-1) =
lgamma(2 * k);
67 constant =
sum(numerator);
69 constant += 0.25 * (K * K - 1) *
LOG_PI
72 constant += 0.25 * K * (K - 2) *
LOG_PI
73 + 0.25 * (3 * K * K - 4 * K) *
LOG_TWO
76 constant = -Km1 *
lgamma(eta + 0.5 * Km1);
77 for (
int k = 1; k <= Km1; k++)
78 constant += 0.5 * k *
LOG_PI +
lgamma(eta + 0.5 * (Km1 - k));
85 template <
bool propto,
86 typename T_y,
typename T_shape>
87 typename boost::math::tools::promote_args<T_y, T_shape>::type
88 lkj_corr_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
90 static const char*
function(
"stan::math::lkj_corr_log");
95 using boost::math::tools::promote_args;
97 typename promote_args<T_y, T_shape>::type lp(0.0);
101 const unsigned int K = y.rows();
115 Eigen::Matrix<T_y, Eigen::Dynamic, 1> values =
116 y.ldlt().vectorD().array().log().matrix();
117 lp += (eta - 1.0) *
sum(values);
121 template <
typename T_y,
typename T_shape>
123 typename boost::math::tools::promote_args<T_y, T_shape>::type
124 lkj_corr_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
125 const T_shape& eta) {
126 return lkj_corr_log<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.
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)
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)
bool check_corr_matrix(const char *function, const char *name, const Eigen::Matrix< T_y, Dynamic, Dynamic > &y)
Return true if the specified matrix is a valid correlation matrix.