1 #ifndef STAN_MATH_PRIM_FUN_COV_MATRIX_FREE_LKJ_HPP
2 #define STAN_MATH_PRIM_FUN_COV_MATRIX_FREE_LKJ_HPP
3 
4 #include <stan/math/prim/meta.hpp>
5 #include <stan/math/prim/err.hpp>
6 #include <stan/math/prim/fun/Eigen.hpp>
7 #include <stan/math/prim/fun/factor_cov_matrix.hpp>
8 
9 namespace stan {
10 namespace math {
11 
12 /**
13  * Return the vector of unconstrained partial correlations and
14  * deviations that transform to the specified covariance matrix.
15  *
16  * <p>The constraining transform is defined as for
17  * <code>cov_matrix_constrain(Matrix, size_t)</code>.  The
18  * inverse first factors out the deviations, then applies the
19  * freeing transform of <code>corr_matrix_free(Matrix&)</code>.
20  *
21  * @tparam T type of the matrix (must be derived from \c Eigen::MatrixBase)
22  * @param y Covariance matrix to free.
23  * @return Vector of unconstrained values that transforms to the
24  * specified covariance matrix.
25  * @throw std::domain_error if the correlation matrix has no
26  *    elements or is not a square matrix.
27  * @throw std::runtime_error if the correlation matrix cannot be
28  *    factorized by factor_cov_matrix()
29  */
30 template <typename T, require_eigen_t<T>* = nullptr>
cov_matrix_free_lkj(const T & y)31 Eigen::Matrix<value_type_t<T>, Eigen::Dynamic, 1> cov_matrix_free_lkj(
32     const T& y) {
33   using Eigen::Array;
34   using Eigen::Dynamic;
35   using Eigen::Matrix;
36   using T_scalar = value_type_t<T>;
37 
38   check_nonzero_size("cov_matrix_free_lkj", "y", y);
39   check_square("cov_matrix_free_lkj", "y", y);
40   Eigen::Index k = y.rows();
41   Eigen::Index k_choose_2 = (k * (k - 1)) / 2;
42   Matrix<T_scalar, Dynamic, 1> x(k_choose_2 + k);
43   bool successful
44       = factor_cov_matrix(y, x.head(k_choose_2).array(), x.tail(k).array());
45   if (!successful) {
46     throw_domain_error("cov_matrix_free_lkj", "factor_cov_matrix failed on y",
47                        "", "");
48   }
49   return x;
50 }
51 
52 /**
53  * Overload of `cov_matrix_free_lkj()` to untransform each matrix
54  * in a standard vector.
55  * @tparam T A standard vector with with a `value_type` which inherits from
56  *  `Eigen::MatrixBase`.
57  * @param x The standard vector to untransform.
58  */
59 template <typename T, require_std_vector_t<T>* = nullptr>
cov_matrix_free_lkj(const T & x)60 auto cov_matrix_free_lkj(const T& x) {
61   return apply_vector_unary<T>::apply(
62       x, [](auto&& v) { return cov_matrix_free_lkj(v); });
63 }
64 
65 }  // namespace math
66 }  // namespace stan
67 
68 #endif
69