1 #ifndef STAN_MATH_PRIM_FUN_READ_COV_MATRIX_HPP
2 #define STAN_MATH_PRIM_FUN_READ_COV_MATRIX_HPP
3 
4 #include <stan/math/prim/fun/read_corr_L.hpp>
5 #include <stan/math/prim/fun/read_cov_L.hpp>
6 #include <stan/math/prim/fun/multiply_lower_tri_self_transpose.hpp>
7 #include <stan/math/prim/fun/Eigen.hpp>
8 
9 namespace stan {
10 namespace math {
11 
12 /**
13  * A generally worse alternative to call prior to evaluating the
14  * density of an elliptical distribution
15  *
16  * @tparam T_CPCs type of \c T_CPCs vector (must be derived from \c
17  * Eigen::ArrayBase and have one compile-time dimension equal to 1)
18  * @tparam T_sds type of \c sds vector (must be derived from \c Eigen::ArrayBase
19  * and have one compile-time dimension equal to 1)
20  * @param CPCs on (-1, 1)
21  * @param sds on (0, inf)
22  * @param log_prob the log probability value to increment with the Jacobian
23  * @return Covariance matrix for specified partial correlations.
24  */
25 template <typename T_CPCs, typename T_sds,
26           require_all_eigen_vector_t<T_CPCs, T_sds>* = nullptr,
27           require_vt_same<T_CPCs, T_sds>* = nullptr>
28 Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
read_cov_matrix(const T_CPCs & CPCs,const T_sds & sds,value_type_t<T_CPCs> & log_prob)29 read_cov_matrix(const T_CPCs& CPCs, const T_sds& sds,
30                 value_type_t<T_CPCs>& log_prob) {
31   Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic> L
32       = read_cov_L(CPCs, sds, log_prob);
33   return multiply_lower_tri_self_transpose(L);
34 }
35 
36 /**
37  * Builds a covariance matrix from CPCs and standard deviations
38  *
39  * @tparam T_CPCs type of \c T_CPCs vector (must be derived from \c
40  * Eigen::ArrayBase and have one compile-time dimension equal to 1)
41  * @tparam T_sds type of \c sds vector (must be derived from \c Eigen::ArrayBase
42  * and have one compile-time dimension equal to 1)
43  * @param CPCs in (-1, 1)
44  * @param sds in (0, inf)
45  */
46 template <typename T_CPCs, typename T_sds,
47           require_all_eigen_vector_t<T_CPCs, T_sds>* = nullptr,
48           require_vt_same<T_CPCs, T_sds>* = nullptr>
49 Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
read_cov_matrix(const T_CPCs & CPCs,const T_sds & sds)50 read_cov_matrix(const T_CPCs& CPCs, const T_sds& sds) {
51   size_t K = sds.rows();
52   Eigen::DiagonalMatrix<value_type_t<T_CPCs>, Eigen::Dynamic> D(K);
53   D.diagonal() = sds;
54   Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic> L
55       = D * read_corr_L(CPCs, K);
56   return multiply_lower_tri_self_transpose(L);
57 }
58 
59 }  // namespace math
60 }  // namespace stan
61 
62 #endif
63