1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_CORR_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_CORR_HPP
14 template <
typename T_shape>
22 const int Km1 = K - 1;
25 Eigen::VectorXd numerator( Km1 / 2 );
26 for(
size_t k = 1; k <= numerator.rows(); k++)
27 numerator(k-1) =
lgamma(2 * k);
28 constant =
sum(numerator);
29 if ( (K % 2) == 1 ) constant += 0.25 * (K * K - 1) * LOG_PI -
30 0.25 * (Km1 * Km1) * LOG_TWO - Km1 *
lgamma( (K + 1) / 2);
31 else constant += 0.25 * K * (K - 2) * LOG_PI +
32 0.25 * (3 * K * K - 4 * K) * LOG_TWO +
36 constant = -Km1 *
lgamma(eta + 0.5 * Km1);
37 for (
size_t k = 1; k <= Km1; k++)
38 constant += 0.5 * k * LOG_PI +
lgamma(eta + 0.5 * (Km1 - k));
45 template <
bool propto,
46 typename T_covar,
typename T_shape>
47 typename boost::math::tools::promote_args<T_covar, T_shape>::type
49 const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L,
52 static const char*
function
53 =
"stan::prob::lkj_corr_cholesky_log(%1%)";
55 using boost::math::tools::promote_args;
60 typename promote_args<T_covar,T_shape>::type lp(0.0);
64 const unsigned int K = L.rows();
71 const int Km1 = K - 1;
72 Eigen::Matrix<T_covar,Eigen::Dynamic,1> log_diagonals =
73 L.diagonal().tail(Km1).array().log();
74 Eigen::Matrix<T_covar,Eigen::Dynamic,1> values(Km1);
75 for (
size_t k = 0; k < Km1; k++)
76 values(k) = (Km1 - k - 1) * log_diagonals(k);
82 values += (2.0 * eta - 2.0) * log_diagonals;
89 template <
typename T_covar,
typename T_shape>
91 typename boost::math::tools::promote_args<T_covar, T_shape>::type
93 const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L,
95 return lkj_corr_cholesky_log<false>(L,eta);
100 template <
bool propto,
101 typename T_y,
typename T_shape>
102 typename boost::math::tools::promote_args<T_y, T_shape>::type
103 lkj_corr_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
104 const T_shape& eta) {
105 static const char*
function =
"stan::prob::lkj_corr_log(%1%)";
112 using boost::math::tools::promote_args;
114 typename promote_args<T_y,T_shape>::type lp(0.0);
117 y.rows(),
"Rows of correlation matrix",
118 y.cols(),
"columns of correlation matrix",
123 const unsigned int K = y.rows();
137 Eigen::Matrix<T_y,Eigen::Dynamic,1> values =
138 y.ldlt().vectorD().array().log().matrix();
139 lp += (eta - 1.0) *
sum(values);
143 template <
typename T_y,
typename T_shape>
145 typename boost::math::tools::promote_args<T_y, T_shape>::type
146 lkj_corr_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
147 const T_shape& eta) {
148 return lkj_corr_log<false>(y,eta);
152 inline Eigen::MatrixXd
156 static const char*
function
157 =
"stan::prob::lkj_corr_cholesky_rng(%1%)";
163 Eigen::ArrayXd CPCs( (K * (K - 1)) / 2 );
164 double alpha = eta + 0.5 * (K - 1);
165 unsigned int count = 0;
166 for (
size_t i = 0; i < (K - 1); i++) {
168 for (
size_t j = i + 1; j < K; j++) {
177 inline Eigen::MatrixXd
182 static const char*
function
183 =
"stan::prob::lkj_corr_rng(%1%)";
Eigen::MatrixXd lkj_corr_cholesky_rng(const size_t K, const double eta, RNG &rng)
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_covar, T_shape >::type lkj_corr_cholesky_log(const Eigen::Matrix< T_covar, Eigen::Dynamic, Eigen::Dynamic > &L, const T_shape &eta)
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)
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 > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
fvar< T > lgamma(const fvar< T > &x)
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Eigen::MatrixXd lkj_corr_rng(const size_t K, const double eta, RNG &rng)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_corr_L(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const size_t K)
Return the Cholesky factor of the correlation matrix of the specified dimensionality corresponding to...
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
bool check_corr_matrix(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is a valid correlation matrix.
double beta_rng(const double alpha, const double beta, RNG &rng)
double sum(std::vector< double > &x)
bool check_not_nan(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is nan.
Eigen::Matrix< fvar< T >, R, R > multiply_lower_tri_self_transpose(const Eigen::Matrix< fvar< T >, R, C > &m)
matrix_d multiply_lower_tri_self_transpose(const matrix_d &L)
Returns the result of multiplying the lower triangular portion of the input matrix by its own transpo...
bool check_lower_triangular(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is lower triangular.