1 #ifndef __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP__
2 #define __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP__
21 template <
bool propto,
22 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape,
24 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
25 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
26 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
27 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
30 static const char*
function =
"stan::prob::lkj_cov_log(%1%)";
35 using boost::math::tools::promote_args;
37 typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
44 if (!
check_positive(
function, eta,
"Shape parameter", &lp, Policy()))
46 if (!
check_finite(
function, mu,
"Location parameter", &lp, Policy()))
48 if (!
check_finite(
function, sigma,
"Scale parameter", &lp, Policy()))
51 for (
int m = 0; m < y.rows(); ++m)
52 for (
int n = 0; n < y.cols(); ++n)
53 if (!
check_finite(
function, y(m,n),
"Covariance matrix", &lp, Policy()))
56 const unsigned int K = y.rows();
57 const Eigen::Array<T_y,Eigen::Dynamic,1> sds
58 = y.diagonal().array().sqrt();
59 for (
unsigned int k = 0; k < K; k++) {
60 lp += lognormal_log<propto>(sds(k), mu(k), sigma(k), Policy());
65 lp += lkj_corr_log<propto,T_y,T_shape,Policy>(y, eta, Policy());
68 Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
69 D.diagonal() = sds.inverse();
70 lp += lkj_corr_log<propto,T_y,T_shape,Policy>(D * y * D, eta, Policy());
74 template <
bool propto,
75 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
77 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
78 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
79 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
80 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
86 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape,
89 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
90 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
91 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
92 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
95 return lkj_cov_log<false>(y,mu,sigma,eta,Policy());
99 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
101 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
102 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
103 const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
104 const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
105 const T_shape& eta) {
113 template <
bool propto,
114 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape,
116 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
117 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
119 const T_scale& sigma,
122 static const char*
function =
"stan::prob::lkj_cov_log(%1%)";
126 using boost::math::tools::promote_args;
128 typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
129 if (!
check_positive(
function, eta,
"Shape parameter", &lp, Policy()))
131 if (!
check_finite(
function, mu,
"Location parameter", &lp, Policy()))
137 const unsigned int K = y.rows();
138 const Eigen::Array<T_y,Eigen::Dynamic,1> sds
139 = y.diagonal().array().sqrt();
140 for (
unsigned int k = 0; k < K; k++) {
141 lp += lognormal_log<propto>(sds(k), mu, sigma, Policy());
146 lp += lkj_corr_log<propto>(y,eta,Policy());
149 Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
150 D.diagonal() = sds.inverse();
151 lp += lkj_corr_log<propto,T_y,T_shape,Policy>(D * y * D, eta, Policy());
155 template <
bool propto,
156 typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
158 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
159 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
161 const T_scale& sigma,
162 const T_shape& eta) {
166 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape,
169 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
170 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
172 const T_scale& sigma,
175 return lkj_cov_log<false>(y,mu,sigma,eta,Policy());
178 template <
typename T_y,
typename T_loc,
typename T_scale,
typename T_shape>
180 typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
181 lkj_cov_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
183 const T_scale& sigma,
184 const T_shape& eta) {
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
Checks if the variable y is finite.
boost::math::policies::policy default_policy
Default error-handling policy from Boost.
bool check_size_match(const char *function, T_size1 i, T_size2 j, T_result *result, const Policy &)
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, const Policy &)
Probability, optimization and sampling library.
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.