Automatic Differentiation
 
Loading...
Searching...
No Matches
read_cov_matrix.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_READ_COV_MATRIX_HPP
2#define STAN_MATH_REV_FUN_READ_COV_MATRIX_HPP
3
11
12namespace stan {
13namespace math {
14
28template <typename T_CPCs, typename T_sds,
29 require_all_var_vector_t<T_CPCs, T_sds>* = nullptr>
30var_value<Eigen::MatrixXd> read_cov_matrix(const T_CPCs& CPCs, const T_sds& sds,
31 scalar_type_t<T_CPCs>& log_prob) {
32 return multiply_lower_tri_self_transpose(read_cov_L(CPCs, sds, log_prob));
33}
34
45template <typename T_CPCs, typename T_sds,
48 const T_sds& sds) {
49 size_t K = sds.rows();
50
53 = sds.val().matrix().asDiagonal() * corr_L.val();
54
55 reverse_pass_callback([sds, corr_L, prod]() mutable {
56 corr_L.adj() += sds.val().matrix().asDiagonal() * prod.adj();
57 sds.adj() += (prod.adj().cwiseProduct(corr_L.val())).rowwise().sum();
58 });
59
61}
62
63} // namespace math
64} // namespace stan
65
66#endif
require_all_t< is_var_vector< std::decay_t< Types > >... > require_all_var_vector_t
Require all of the types satisfy is_var_vector.
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_cov_matrix(const T_CPCs &CPCs, const T_sds &sds, value_type_t< T_CPCs > &log_prob)
A generally worse alternative to call prior to evaluating the density of an elliptical distribution.
value_type_t< T > prod(const T &m)
Calculates product of given kernel generator expression elements.
Definition prod.hpp:21
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
Eigen::Matrix< value_type_t< T_CPCs >, Eigen::Dynamic, Eigen::Dynamic > read_cov_L(const T_CPCs &CPCs, const T_sds &sds, value_type_t< T_CPCs > &log_prob)
This is the function that should be called prior to evaluating the density of any elliptical distribu...
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...
typename scalar_type< T >::type scalar_type_t
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...