1#ifndef STAN_MATH_REV_FUN_MDIVIDE_LEFT_SPD_HPP
2#define STAN_MATH_REV_FUN_MDIVIDE_LEFT_SPD_HPP
18template <
int R1,
int C1,
int R2,
int C2>
23 Eigen::LLT<Eigen::Matrix<double, R1, C1>>
llt_;
24 Eigen::Matrix<double, R2, C2>
C_;
27template <
int R1,
int C1,
int R2,
int C2>
38 const Eigen::Matrix<var, R2, C2> &B)
55 alloc_->llt_ = A.val().llt();
60 =
alloc_->C_.unaryExpr([](
double x) {
return new vari(x,
false); });
65 alloc_->llt_.solveInPlace(adjB);
67 -= adjB *
alloc_->C_.transpose();
72template <
int R1,
int C1,
int R2,
int C2>
82 const Eigen::Matrix<var, R2, C2> &B)
100 =
alloc_->C_.unaryExpr([](
double x) {
return new vari(x,
false); });
105 alloc_->llt_.solveInPlace(adjB);
110template <
int R1,
int C1,
int R2,
int C2>
120 const Eigen::Matrix<double, R2, C2> &B)
132 alloc_->llt_ = A.val().llt();
137 =
alloc_->C_.unaryExpr([](
double x) {
return new vari(x,
false); });
149 typename EigMat1,
typename EigMat2,
151inline Eigen::Matrix<
var, EigMat1::RowsAtCompileTime,
152 EigMat2::ColsAtCompileTime>
154 constexpr int R1 = EigMat1::RowsAtCompileTime;
155 constexpr int C1 = EigMat1::ColsAtCompileTime;
156 constexpr int R2 = EigMat2::RowsAtCompileTime;
157 constexpr int C2 = EigMat2::ColsAtCompileTime;
158 static constexpr const char *function =
"mdivide_left_spd";
160 const auto &A_ref =
to_ref(A);
164 return {0, b.cols()};
174 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
175 res.vi() = Eigen::Map<matrix_vi>(&baseVari->
variRefC_[0], b.rows(), b.cols());
179template <
typename EigMat1,
typename EigMat2,
182inline Eigen::Matrix<
var, EigMat1::RowsAtCompileTime,
183 EigMat2::ColsAtCompileTime>
185 constexpr int R1 = EigMat1::RowsAtCompileTime;
186 constexpr int C1 = EigMat1::ColsAtCompileTime;
187 constexpr int R2 = EigMat2::RowsAtCompileTime;
188 constexpr int C2 = EigMat2::ColsAtCompileTime;
189 static constexpr const char *function =
"mdivide_left_spd";
191 const auto &A_ref =
to_ref(A);
195 return {0, b.cols()};
205 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
206 res.vi() = Eigen::Map<matrix_vi>(&baseVari->
variRefC_[0], b.rows(), b.cols());
210template <
typename EigMat1,
typename EigMat2,
213inline Eigen::Matrix<
var, EigMat1::RowsAtCompileTime,
214 EigMat2::ColsAtCompileTime>
216 constexpr int R1 = EigMat1::RowsAtCompileTime;
217 constexpr int C1 = EigMat1::ColsAtCompileTime;
218 constexpr int R2 = EigMat2::RowsAtCompileTime;
219 constexpr int C2 = EigMat2::ColsAtCompileTime;
220 static constexpr const char *function =
"mdivide_left_spd";
222 const auto &A_ref =
to_ref(A);
226 return {0, b.cols()};
233 internal::mdivide_left_spd_dv_vari<R1, C1, R2, C2> *baseVari
234 =
new internal::mdivide_left_spd_dv_vari<R1, C1, R2, C2>(A_ref, b);
236 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
237 res.vi() = Eigen::Map<matrix_vi>(&baseVari->variRefC_[0], b.rows(), b.cols());
260template <
typename T1,
typename T2, require_all_matrix_t<T1, T2> * =
nullptr,
261 require_any_var_matrix_t<T1, T2> * =
nullptr>
267 return ret_type(ret_val_type(0, B.cols()));
279 auto A_llt = arena_A.val().llt();
289 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
290 arena_A_llt.template triangularView<Eigen::Lower>()
294 arena_A.adj() -= adjB * res.val_op().transpose();
295 arena_B.adj() += adjB;
298 return ret_type(res);
305 auto A_llt = arena_A.val().llt();
315 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
316 arena_A_llt.template triangularView<Eigen::Lower>()
320 arena_A.adj() -= adjB * res.val().transpose().eval();
323 return ret_type(res);
331 auto A_llt = A_ref.llt();
341 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
342 arena_A_llt.template triangularView<Eigen::Lower>()
346 arena_B.adj() += adjB;
349 return ret_type(res);
A chainable_alloc is an object which is constructed and destructed normally but the memory lifespan i...
virtual ~mdivide_left_spd_alloc()
Eigen::Matrix< double, R2, C2 > C_
Eigen::LLT< Eigen::Matrix< double, R1, C1 > > llt_
mdivide_left_spd_dv_vari(const Eigen::Matrix< double, R1, C1 > &A, const Eigen::Matrix< var, R2, C2 > &B)
mdivide_left_spd_alloc< R1, C1, R2, C2 > * alloc_
mdivide_left_spd_alloc< R1, C1, R2, C2 > * alloc_
mdivide_left_spd_vd_vari(const Eigen::Matrix< var, R1, C1 > &A, const Eigen::Matrix< double, R2, C2 > &B)
mdivide_left_spd_alloc< R1, C1, R2, C2 > * alloc_
mdivide_left_spd_vv_vari(const Eigen::Matrix< var, R1, C1 > &A, const Eigen::Matrix< var, R2, C2 > &B)
require_t< container_type_check_base< is_eigen_matrix_base, value_type_t, TypeCheck, Check... > > require_eigen_matrix_base_vt
Require type satisfies is_eigen_matrix_base.
require_all_t< container_type_check_base< is_eigen_matrix_base, value_type_t, TypeCheck, Check >... > require_all_eigen_matrix_base_vt
Require all of the types satisfy is_eigen_matrix_base.
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
int64_t cols(const T_x &x)
Returns the number of columns in the specified kernel generator expression.
int64_t rows(const T_x &x)
Returns the number of rows in the specified kernel generator expression.
void check_multiplicable(const char *function, const char *name1, const T1 &y1, const char *name2, const T2 &y2)
Check if the matrices can be multiplied.
typename promote_scalar_type< std::decay_t< T >, std::decay_t< S > >::type promote_scalar_t
Eigen::Matrix< return_type_t< EigMat1, EigMat2 >, EigMat1::RowsAtCompileTime, EigMat2::ColsAtCompileTime > mdivide_left_spd(const EigMat1 &A, const EigMat2 &b)
Returns the solution of the system Ax=b where A is symmetric positive definite.
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
T value_of(const fvar< T > &v)
Return the value of the specified variable.
void check_pos_definite(const char *function, const char *name, const EigMat &y)
Check if the specified square, symmetric matrix is positive definite.
vari_value< double > vari
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > matrix_d
Type for matrix of double values.
typename plain_type< T >::type plain_type_t
typename internal::arena_type_impl< std::decay_t< T > >::type arena_t
Determines a type that can be used in place of T that does any dynamic allocations on the AD stack.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
This struct always provides access to the autodiff stack using the singleton pattern.