Automatic Differentiation
 
Loading...
Searching...
No Matches
mdivide_left_spd.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_MDIVIDE_LEFT_SPD_HPP
2#define STAN_MATH_REV_FUN_MDIVIDE_LEFT_SPD_HPP
3
12#include <vector>
13
14namespace stan {
15namespace math {
16namespace internal {
17
18template <int R1, int C1, int R2, int C2>
20 public:
22
23 Eigen::LLT<Eigen::Matrix<double, R1, C1>> llt_;
24 Eigen::Matrix<double, R2, C2> C_;
25};
26
27template <int R1, int C1, int R2, int C2>
29 public:
30 int M_; // A.rows() = A.cols() = B.rows()
31 int N_; // B.cols()
36
37 mdivide_left_spd_vv_vari(const Eigen::Matrix<var, R1, C1> &A,
38 const Eigen::Matrix<var, R2, C2> &B)
39 : vari(0.0),
40 M_(A.rows()),
41 N_(B.cols()),
42 variRefA_(reinterpret_cast<vari **>(
43 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * A.rows()
44 * A.cols()))),
45 variRefB_(reinterpret_cast<vari **>(
46 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
47 * B.cols()))),
48 variRefC_(reinterpret_cast<vari **>(
49 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
50 * B.cols()))),
51 alloc_(new mdivide_left_spd_alloc<R1, C1, R2, C2>()) {
52 Eigen::Map<matrix_vi>(variRefA_, M_, M_) = A.vi();
53 Eigen::Map<matrix_vi>(variRefB_, M_, N_) = B.vi();
54 alloc_->C_ = B.val();
55 alloc_->llt_ = A.val().llt();
56 check_pos_definite("mdivide_left_spd", "A", alloc_->llt_);
57 alloc_->llt_.solveInPlace(alloc_->C_);
58
59 Eigen::Map<matrix_vi>(variRefC_, M_, N_)
60 = alloc_->C_.unaryExpr([](double x) { return new vari(x, false); });
61 }
62
63 virtual void chain() {
64 matrix_d adjB = Eigen::Map<matrix_vi>(variRefC_, M_, N_).adj();
65 alloc_->llt_.solveInPlace(adjB);
66 Eigen::Map<matrix_vi>(variRefA_, M_, M_).adj()
67 -= adjB * alloc_->C_.transpose();
68 Eigen::Map<matrix_vi>(variRefB_, M_, N_).adj() += adjB;
69 }
70};
71
72template <int R1, int C1, int R2, int C2>
74 public:
75 int M_; // A.rows() = A.cols() = B.rows()
76 int N_; // B.cols()
80
81 mdivide_left_spd_dv_vari(const Eigen::Matrix<double, R1, C1> &A,
82 const Eigen::Matrix<var, R2, C2> &B)
83 : vari(0.0),
84 M_(A.rows()),
85 N_(B.cols()),
86 variRefB_(reinterpret_cast<vari **>(
87 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
88 * B.cols()))),
89 variRefC_(reinterpret_cast<vari **>(
90 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
91 * B.cols()))),
92 alloc_(new mdivide_left_spd_alloc<R1, C1, R2, C2>()) {
93 alloc_->C_ = B.val();
94 Eigen::Map<matrix_vi>(variRefB_, M_, N_) = B.vi();
95 alloc_->llt_ = A.llt();
96 check_pos_definite("mdivide_left_spd", "A", alloc_->llt_);
97 alloc_->llt_.solveInPlace(alloc_->C_);
98
99 Eigen::Map<matrix_vi>(variRefC_, M_, N_)
100 = alloc_->C_.unaryExpr([](double x) { return new vari(x, false); });
101 }
102
103 virtual void chain() {
104 matrix_d adjB = Eigen::Map<matrix_vi>(variRefC_, M_, N_).adj();
105 alloc_->llt_.solveInPlace(adjB);
106 Eigen::Map<matrix_vi>(variRefB_, M_, N_).adj() += adjB;
107 }
108};
109
110template <int R1, int C1, int R2, int C2>
112 public:
113 int M_; // A.rows() = A.cols() = B.rows()
114 int N_; // B.cols()
118
119 mdivide_left_spd_vd_vari(const Eigen::Matrix<var, R1, C1> &A,
120 const Eigen::Matrix<double, R2, C2> &B)
121 : vari(0.0),
122 M_(A.rows()),
123 N_(B.cols()),
124 variRefA_(reinterpret_cast<vari **>(
125 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * A.rows()
126 * A.cols()))),
127 variRefC_(reinterpret_cast<vari **>(
128 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
129 * B.cols()))),
130 alloc_(new mdivide_left_spd_alloc<R1, C1, R2, C2>()) {
131 Eigen::Map<matrix_vi>(variRefA_, M_, M_) = A.vi();
132 alloc_->llt_ = A.val().llt();
133 check_pos_definite("mdivide_left_spd", "A", alloc_->llt_);
134 alloc_->C_ = alloc_->llt_.solve(B);
135
136 Eigen::Map<matrix_vi>(variRefC_, M_, N_)
137 = alloc_->C_.unaryExpr([](double x) { return new vari(x, false); });
138 }
139
140 virtual void chain() {
141 matrix_d adjC = Eigen::Map<matrix_vi>(variRefC_, M_, N_).adj();
142 Eigen::Map<matrix_vi>(variRefA_, M_, M_).adj()
143 -= alloc_->llt_.solve(adjC * alloc_->C_.transpose());
144 }
145};
146} // namespace internal
147
148template <
149 typename EigMat1, typename EigMat2,
151inline Eigen::Matrix<var, EigMat1::RowsAtCompileTime,
152 EigMat2::ColsAtCompileTime>
153mdivide_left_spd(const EigMat1 &A, const EigMat2 &b) {
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";
159 check_multiplicable(function, "A", A, "b", b);
160 const auto &A_ref = to_ref(A);
161 check_symmetric(function, "A", A_ref);
162 check_not_nan(function, "A", A_ref);
163 if (A.size() == 0) {
164 return {0, b.cols()};
165 }
166
167 // NOTE: this is not a memory leak, this vari is used in the
168 // expression graph to evaluate the adjoint, but is not needed
169 // for the returned matrix. Memory will be cleaned up with the
170 // arena allocator.
173
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());
176 return res;
177}
178
179template <typename EigMat1, typename EigMat2,
182inline Eigen::Matrix<var, EigMat1::RowsAtCompileTime,
183 EigMat2::ColsAtCompileTime>
184mdivide_left_spd(const EigMat1 &A, const EigMat2 &b) {
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";
190 check_multiplicable(function, "A", A, "b", b);
191 const auto &A_ref = to_ref(A);
192 check_symmetric(function, "A", A_ref);
193 check_not_nan(function, "A", A_ref);
194 if (A.size() == 0) {
195 return {0, b.cols()};
196 }
197
198 // NOTE: this is not a memory leak, this vari is used in the
199 // expression graph to evaluate the adjoint, but is not needed
200 // for the returned matrix. Memory will be cleaned up with the
201 // arena allocator.
204
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());
207 return res;
208}
209
210template <typename EigMat1, typename EigMat2,
213inline Eigen::Matrix<var, EigMat1::RowsAtCompileTime,
214 EigMat2::ColsAtCompileTime>
215mdivide_left_spd(const EigMat1 &A, const EigMat2 &b) {
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";
221 check_multiplicable(function, "A", A, "b", b);
222 const auto &A_ref = to_ref(A);
223 check_symmetric(function, "A", A_ref);
224 check_not_nan(function, "A", A_ref);
225 if (A.size() == 0) {
226 return {0, b.cols()};
227 }
228
229 // NOTE: this is not a memory leak, this vari is used in the
230 // expression graph to evaluate the adjoint, but is not needed
231 // for the returned matrix. Memory will be cleaned up with the
232 // arena allocator.
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);
235
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());
238
239 return res;
240}
241
260template <typename T1, typename T2, require_all_matrix_t<T1, T2> * = nullptr,
261 require_any_var_matrix_t<T1, T2> * = nullptr>
262inline auto mdivide_left_spd(const T1 &A, const T2 &B) {
263 using ret_val_type = plain_type_t<decltype(value_of(A) * value_of(B))>;
264 using ret_type = var_value<ret_val_type>;
265
266 if (A.size() == 0) {
267 return ret_type(ret_val_type(0, B.cols()));
268 }
269
270 check_multiplicable("mdivide_left_spd", "A", A, "B", B);
271
275
276 check_symmetric("mdivide_left_spd", "A", arena_A.val());
277 check_not_nan("mdivide_left_spd", "A", arena_A.val());
278
279 auto A_llt = arena_A.val().llt();
280
281 check_pos_definite("mdivide_left_spd", "A", A_llt);
282
283 arena_t<Eigen::MatrixXd> arena_A_llt = A_llt.matrixL();
284 arena_t<ret_type> res = A_llt.solve(arena_B.val());
285
286 reverse_pass_callback([arena_A, arena_B, arena_A_llt, res]() mutable {
287 promote_scalar_t<double, T2> adjB = res.adj();
288
289 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
290 arena_A_llt.template triangularView<Eigen::Lower>()
291 .transpose()
292 .solveInPlace(adjB);
293
294 arena_A.adj() -= adjB * res.val_op().transpose();
295 arena_B.adj() += adjB;
296 });
297
298 return ret_type(res);
299 } else if (!is_constant<T1>::value) {
301
302 check_symmetric("mdivide_left_spd", "A", arena_A.val());
303 check_not_nan("mdivide_left_spd", "A", arena_A.val());
304
305 auto A_llt = arena_A.val().llt();
306
307 check_pos_definite("mdivide_left_spd", "A", A_llt);
308
309 arena_t<Eigen::MatrixXd> arena_A_llt = A_llt.matrixL();
310 arena_t<ret_type> res = A_llt.solve(value_of(B));
311
312 reverse_pass_callback([arena_A, arena_A_llt, res]() mutable {
313 promote_scalar_t<double, T2> adjB = res.adj();
314
315 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
316 arena_A_llt.template triangularView<Eigen::Lower>()
317 .transpose()
318 .solveInPlace(adjB);
319
320 arena_A.adj() -= adjB * res.val().transpose().eval();
321 });
322
323 return ret_type(res);
324 } else {
325 const auto &A_ref = to_ref(value_of(A));
327
328 check_symmetric("mdivide_left_spd", "A", A_ref);
329 check_not_nan("mdivide_left_spd", "A", A_ref);
330
331 auto A_llt = A_ref.llt();
332
333 check_pos_definite("mdivide_left_spd", "A", A_llt);
334
335 arena_t<Eigen::MatrixXd> arena_A_llt = A_llt.matrixL();
336 arena_t<ret_type> res = A_llt.solve(arena_B.val());
337
338 reverse_pass_callback([arena_B, arena_A_llt, res]() mutable {
339 promote_scalar_t<double, T2> adjB = res.adj();
340
341 arena_A_llt.template triangularView<Eigen::Lower>().solveInPlace(adjB);
342 arena_A_llt.template triangularView<Eigen::Lower>()
343 .transpose()
344 .solveInPlace(adjB);
345
346 arena_B.adj() += adjB;
347 });
348
349 return ret_type(res);
350 }
351}
352
353} // namespace math
354} // namespace stan
355#endif
A chainable_alloc is an object which is constructed and destructed normally but the memory lifespan i...
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.
Definition cols.hpp:21
int64_t rows(const T_x &x)
Returns the number of rows in the specified kernel generator expression.
Definition rows.hpp:22
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.
Definition value_of.hpp:18
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
Definition vari.hpp:197
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:17
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
var_value< double > var
Definition var.hpp:1187
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > matrix_d
Type for matrix of double values.
Definition typedefs.hpp:19
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.