1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_COV_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_COV_LOG_HPP
20 template <
bool propto,
21 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
23 boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
24 lkj_cov_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
25 const Eigen::Matrix<T_loc, Eigen::Dynamic, 1>& mu,
26 const Eigen::Matrix<T_scale, Eigen::Dynamic, 1>& sigma,
28 static const char*
function(
"stan::math::lkj_cov_log");
34 using boost::math::tools::promote_args;
36 typename promote_args<T_y, T_loc, T_scale, T_shape>::type lp(0.0);
38 "Rows of location parameter", mu.rows(),
39 "columns of scale parameter", sigma.rows());
42 "Rows of random variable", y.rows(),
43 "rows of location parameter", mu.rows());
48 for (
int m = 0; m < y.rows(); ++m)
49 for (
int n = 0; n < y.cols(); ++n)
52 const unsigned int K = y.rows();
53 const Eigen::Array<T_y, Eigen::Dynamic, 1> sds
54 = y.diagonal().array().sqrt();
55 for (
unsigned int k = 0; k < K; k++) {
56 lp += lognormal_log<propto>(sds(k), mu(k), sigma(k));
61 lp += lkj_corr_log<propto, T_y, T_shape>(y, eta);
64 Eigen::DiagonalMatrix<T_y, Eigen::Dynamic> D(K);
65 D.diagonal() = sds.inverse();
66 lp += lkj_corr_log<propto, T_y, T_shape>(D * y * D, eta);
70 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
73 boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
74 lkj_cov_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
75 const Eigen::Matrix<T_loc, Eigen::Dynamic, 1>& mu,
76 const Eigen::Matrix<T_scale, Eigen::Dynamic, 1>& sigma,
78 return lkj_cov_log<false>(y, mu, sigma, eta);
83 template <
bool propto,
84 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
86 boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
87 lkj_cov_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
91 static const char*
function(
"stan::math::lkj_cov_log");
95 using boost::math::tools::promote_args;
97 typename promote_args<T_y, T_loc, T_scale, T_shape>::type lp(0.0);
102 const unsigned int K = y.rows();
103 const Eigen::Array<T_y, Eigen::Dynamic, 1> sds
104 = y.diagonal().array().sqrt();
105 for (
unsigned int k = 0; k < K; k++) {
106 lp += lognormal_log<propto>(sds(k), mu, sigma);
111 lp += lkj_corr_log<propto>(y, eta);
114 Eigen::DiagonalMatrix<T_y, Eigen::Dynamic> D(K);
115 D.diagonal() = sds.inverse();
116 lp += lkj_corr_log<propto, T_y, T_shape>(D * y * D, eta);
120 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
122 typename boost::math::tools::promote_args
123 <T_y, T_loc, T_scale, T_shape>::type
124 lkj_cov_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
126 const T_scale& sigma,
127 const T_shape& eta) {
128 return lkj_cov_log<false>(y, mu, sigma, eta);
boost::math::tools::promote_args< T_y, T_loc, T_scale, T_shape >::type lkj_cov_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const Eigen::Matrix< T_loc, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< T_scale, Eigen::Dynamic, 1 > &sigma, const T_shape &eta)
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.
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
bool check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Return true if the provided sizes match.
bool check_finite(const char *function, const char *name, const T_y &y)
Return true if y is finite.
bool check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is square.