1 #ifndef STAN_MATH_PRIM_FUN_READ_CORR_MATRIX_HPP
2 #define STAN_MATH_PRIM_FUN_READ_CORR_MATRIX_HPP
3 
4 #include <stan/math/prim/fun/Eigen.hpp>
5 #include <stan/math/prim/fun/read_corr_L.hpp>
6 #include <stan/math/prim/fun/multiply_lower_tri_self_transpose.hpp>
7 
8 namespace stan {
9 namespace math {
10 
11 /**
12  * Return the correlation matrix of the specified dimensionality
13  * corresponding to the specified canonical partial correlations.
14  *
15  * <p>See <code>read_corr_matrix(Array, size_t, T)</code>
16  * for more information.
17  *
18  * @tparam T_CPCs type of the array (must be derived from \c Eigen::ArrayBase
19  * and have one compile-time dimension equal to 1)
20  * @param CPCs The (K choose 2) canonical partial correlations in (-1, 1).
21  * @param K Dimensionality of correlation matrix.
22  * @return Cholesky factor of correlation matrix for specified
23  * canonical partial correlations.
24  */
25 template <typename T_CPCs, require_eigen_vector_t<T_CPCs>* = nullptr>
26 Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
read_corr_matrix(const T_CPCs & CPCs,size_t K)27 read_corr_matrix(const T_CPCs& CPCs, size_t K) {
28   if (K == 0) {
29     return {};
30   }
31 
32   return multiply_lower_tri_self_transpose(read_corr_L(CPCs, K));
33 }
34 
35 /**
36  * Return the correlation matrix of the specified dimensionality
37  * corresponding to the specified canonical partial correlations,
38  * incrementing the specified scalar reference with the log
39  * absolute determinant of the Jacobian of the transformation.
40  *
41  * It is usually preferable to utilize the version that returns
42  * the Cholesky factor of the correlation matrix rather than the
43  * correlation matrix itself in statistical calculations.
44  *
45  * @tparam T_CPCs type of the array (must be derived from \c Eigen::ArrayBase
46  * and have one compile-time dimension equal to 1)
47  * @param CPCs The (K choose 2) canonical partial correlations in
48  * (-1, 1).
49  * @param K Dimensionality of correlation matrix.
50  * @param log_prob Reference to variable to increment with the log
51  * Jacobian determinant.
52  * @return Correlation matrix for specified partial correlations.
53  */
54 template <typename T_CPCs, require_eigen_vector_t<T_CPCs>* = nullptr>
55 Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
read_corr_matrix(const T_CPCs & CPCs,size_t K,value_type_t<T_CPCs> & log_prob)56 read_corr_matrix(const T_CPCs& CPCs, size_t K, value_type_t<T_CPCs>& log_prob) {
57   if (K == 0) {
58     return {};
59   }
60 
61   return multiply_lower_tri_self_transpose(read_corr_L(CPCs, K, log_prob));
62 }
63 
64 }  // namespace math
65 }  // namespace stan
66 
67 #endif
68