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