Automatic Differentiation
 
Loading...
Searching...
No Matches
read_corr_matrix.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_FUN_READ_CORR_MATRIX_HPP
2#define STAN_MATH_PRIM_FUN_READ_CORR_MATRIX_HPP
3
7
8namespace stan {
9namespace math {
10
25template <typename T_CPCs, require_eigen_vector_t<T_CPCs>* = nullptr>
26Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
27read_corr_matrix(const T_CPCs& CPCs, size_t K) {
28 if (K == 0) {
29 return {};
30 }
31
33}
34
54template <typename T_CPCs, require_eigen_vector_t<T_CPCs>* = nullptr>
55Eigen::Matrix<value_type_t<T_CPCs>, Eigen::Dynamic, Eigen::Dynamic>
56read_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
typename value_type< T >::type value_type_t
Helper function for accessing underlying type.
Eigen::Matrix< value_type_t< EigMat >, EigMat::RowsAtCompileTime, EigMat::RowsAtCompileTime > multiply_lower_tri_self_transpose(const EigMat &m)
Eigen::Matrix< value_type_t< T_CPCs >, Eigen::Dynamic, Eigen::Dynamic > read_corr_matrix(const T_CPCs &CPCs, size_t K)
Return the correlation matrix of the specified dimensionality corresponding to the specified canonica...
Eigen::Matrix< value_type_t< T >, Eigen::Dynamic, Eigen::Dynamic > read_corr_L(const T &CPCs, size_t K)
Return the Cholesky factor of the correlation matrix of the specified dimensionality corresponding to...
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...