1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP
20 template <
bool propto,
21 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
22 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
23 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
24 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
25 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
27 static const char*
function =
"stan::prob::lkj_cov_log(%1%)";
32 using boost::math::tools::promote_args;
34 typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
36 mu.rows(),
"Rows of location parameter",
37 sigma.rows(),
"columns of scale parameter",
40 y.rows(),
"Rows of random variable",
41 y.cols(),
"columns of random variable",
44 y.rows(),
"Rows of random variable",
45 mu.rows(),
"rows of location parameter",
51 for (
int m = 0; m < y.rows(); ++m)
52 for (
int n = 0; n < y.cols(); ++n)
53 check_finite(
function, y(m,n),
"Covariance matrix", &lp);
55 const unsigned int K = y.rows();
56 const Eigen::Array<T_y,Eigen::Dynamic,1> sds
57 = y.diagonal().array().sqrt();
58 for (
unsigned int k = 0; k < K; k++) {
59 lp += lognormal_log<propto>(sds(k), mu(k), sigma(k));
64 lp += lkj_corr_log<propto,T_y,T_shape>(y, eta);
67 Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
68 D.diagonal() = sds.inverse();
69 lp += lkj_corr_log<propto,T_y,T_shape>(D * y * D, eta);
73 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
75 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
76 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
77 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
78 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
80 return lkj_cov_log<false>(y,mu,sigma,eta);
85 template <
bool propto,
86 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
87 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
88 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
92 static const char*
function =
"stan::prob::lkj_cov_log(%1%)";
96 using boost::math::tools::promote_args;
98 typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
103 const unsigned int K = y.rows();
104 const Eigen::Array<T_y,Eigen::Dynamic,1> sds
105 = y.diagonal().array().sqrt();
106 for (
unsigned int k = 0; k < K; k++) {
107 lp += lognormal_log<propto>(sds(k), mu, sigma);
112 lp += lkj_corr_log<propto>(y,eta);
115 Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
116 D.diagonal() = sds.inverse();
117 lp += lkj_corr_log<propto,T_y,T_shape>(D * y * D, eta);
121 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
123 typename boost::math::tools::promote_args<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);
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_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)
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)