1 #ifndef STAN_MATH_REV_FUN_READ_COV_L_HPP
2 #define STAN_MATH_REV_FUN_READ_COV_L_HPP
3
4 #include <stan/math/prim/fun/Eigen.hpp>
5 #include <stan/math/rev/core.hpp>
6 #include <stan/math/prim/fun/log.hpp>
7 #include <stan/math/prim/fun/read_cov_L.hpp>
8 #include <stan/math/rev/fun/read_corr_L.hpp>
9 #include <stan/math/prim/fun/sum.hpp>
10 #include <stan/math/prim/fun/constants.hpp>
11
12 namespace stan {
13 namespace math {
14
15 /**
16 * This is the function that should be called prior to evaluating
17 * the density of any elliptical distribution
18 *
19 * @tparam T_CPCs type of CPCs vector (must be a `var_value<T>` where `T`
20 * inherits from EigenBase)
21 * @tparam T_sds type of sds vector (must be a `var_value<T>` where `T`
22 * inherits from EigenBase)
23 * @param CPCs on (-1, 1)
24 * @param sds on (0, inf)
25 * @param log_prob the log probability value to increment with the Jacobian
26 * @return Cholesky factor of covariance matrix for specified
27 * partial correlations.
28 */
29 template <typename T_CPCs, typename T_sds,
30 require_any_var_vector_t<T_CPCs, T_sds>* = nullptr,
31 require_vt_same<T_CPCs, T_sds>* = nullptr>
read_cov_L(const T_CPCs & CPCs,const T_sds & sds,scalar_type_t<T_CPCs> & log_prob)32 inline auto read_cov_L(const T_CPCs& CPCs, const T_sds& sds,
33 scalar_type_t<T_CPCs>& log_prob) {
34 size_t K = sds.rows();
35 // adjust due to transformation from correlations to covariances
36 log_prob += (sum(log(sds.val())) + LOG_TWO) * K;
37
38 auto corr_L = read_corr_L(CPCs, K, log_prob);
39 var_value<Eigen::MatrixXd> res
40 = sds.val().matrix().asDiagonal() * corr_L.val();
41
42 reverse_pass_callback([CPCs, sds, corr_L, log_prob, res]() mutable {
43 size_t K = sds.size();
44
45 corr_L.adj() += sds.val().matrix().asDiagonal() * res.adj();
46 sds.adj() += (res.adj().cwiseProduct(corr_L.val())).rowwise().sum();
47
48 sds.adj() += (K * log_prob.adj() / sds.val().array()).matrix();
49 });
50
51 return res;
52 }
53
54 } // namespace math
55 } // namespace stan
56
57 #endif
58