Automatic Differentiation
 
Loading...
Searching...
No Matches
mdivide_left_tri.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_MDIVIDE_LEFT_TRI_HPP
2#define STAN_MATH_REV_FUN_MDIVIDE_LEFT_TRI_HPP
3
11
12namespace stan {
13namespace math {
14
15namespace internal {
16
17template <Eigen::UpLoType TriView, int R1, int C1, int R2, int C2>
19 public:
20 int M_; // A.rows() = A.cols() = B.rows()
21 int N_; // B.cols()
22 double *A_;
23 double *C_;
27
28 mdivide_left_tri_vv_vari(const Eigen::Matrix<var, R1, C1> &A,
29 const Eigen::Matrix<var, R2, C2> &B)
30 : vari(0.0),
31 M_(A.rows()),
32 N_(B.cols()),
33 A_(reinterpret_cast<double *>(
34 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * A.rows()
35 * A.cols()))),
36 C_(reinterpret_cast<double *>(
37 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * B.rows()
38 * B.cols()))),
39 variRefA_(reinterpret_cast<vari **>(
40 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * A.rows()
41 * (A.rows() + 1) / 2))),
42 variRefB_(reinterpret_cast<vari **>(
43 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
44 * B.cols()))),
45 variRefC_(reinterpret_cast<vari **>(
46 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
47 * B.cols()))) {
48 using Eigen::Map;
49
50 size_t pos = 0;
51 if (TriView == Eigen::Lower) {
52 for (size_type j = 0; j < M_; j++) {
53 for (size_type i = j; i < M_; i++) {
54 variRefA_[pos++] = A(i, j).vi_;
55 }
56 }
57 } else if (TriView == Eigen::Upper) {
58 for (size_type j = 0; j < M_; j++) {
59 for (size_type i = 0; i < j + 1; i++) {
60 variRefA_[pos++] = A(i, j).vi_;
61 }
62 }
63 }
64
65 Map<matrix_d> c_map(C_, M_, N_);
66 Map<matrix_d> a_map(A_, M_, M_);
67 a_map = A.val();
68 c_map = B.val();
69 Map<matrix_vi>(variRefB_, M_, N_) = B.vi();
70
71 c_map = a_map.template triangularView<TriView>().solve(c_map);
72
73 Map<matrix_vi>(variRefC_, M_, N_)
74 = c_map.unaryExpr([](double x) { return new vari(x, false); });
75 }
76
77 virtual void chain() {
78 using Eigen::Map;
79 matrix_d adjA;
80 matrix_d adjB;
81
82 adjB = Map<matrix_d>(A_, M_, M_)
83 .template triangularView<TriView>()
84 .transpose()
85 .solve(Map<matrix_vi>(variRefC_, M_, N_).adj());
86 adjA = -adjB * Map<matrix_d>(C_, M_, N_).transpose();
87
88 size_t pos = 0;
89 if (TriView == Eigen::Lower) {
90 for (size_type j = 0; j < adjA.cols(); j++) {
91 for (size_type i = j; i < adjA.rows(); i++) {
92 variRefA_[pos++]->adj_ += adjA(i, j);
93 }
94 }
95 } else if (TriView == Eigen::Upper) {
96 for (size_type j = 0; j < adjA.cols(); j++) {
97 for (size_type i = 0; i < j + 1; i++) {
98 variRefA_[pos++]->adj_ += adjA(i, j);
99 }
100 }
101 }
102 Map<matrix_vi>(variRefB_, M_, N_).adj() += adjB;
103 }
104};
105
106template <Eigen::UpLoType TriView, int R1, int C1, int R2, int C2>
108 public:
109 int M_; // A.rows() = A.cols() = B.rows()
110 int N_; // B.cols()
111 double *A_;
112 double *C_;
115
116 mdivide_left_tri_dv_vari(const Eigen::Matrix<double, R1, C1> &A,
117 const Eigen::Matrix<var, R2, C2> &B)
118 : vari(0.0),
119 M_(A.rows()),
120 N_(B.cols()),
121 A_(reinterpret_cast<double *>(
122 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * A.rows()
123 * A.cols()))),
124 C_(reinterpret_cast<double *>(
125 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * B.rows()
126 * B.cols()))),
127 variRefB_(reinterpret_cast<vari **>(
128 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
129 * B.cols()))),
130 variRefC_(reinterpret_cast<vari **>(
131 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
132 * B.cols()))) {
133 using Eigen::Map;
134
135 Map<matrix_d>(A_, M_, M_) = A;
136 Map<matrix_vi>(variRefB_, M_, N_) = B.vi();
137 Map<matrix_d> c_map(C_, M_, N_);
138 c_map = B.val();
139
140 c_map = Map<matrix_d>(A_, M_, M_)
141 .template triangularView<TriView>()
142 .solve(c_map);
143
144 Map<matrix_vi>(variRefC_, M_, N_)
145 = c_map.unaryExpr([](double x) { return new vari(x, false); });
146 }
147
148 virtual void chain() {
149 using Eigen::Map;
150
151 Map<matrix_vi>(variRefB_, M_, N_).adj()
152 += Map<matrix_d>(A_, M_, M_)
153 .template triangularView<TriView>()
154 .transpose()
155 .solve(Map<matrix_vi>(variRefC_, M_, N_).adj());
156 }
157};
158
159template <Eigen::UpLoType TriView, int R1, int C1, int R2, int C2>
161 public:
162 int M_; // A.rows() = A.cols() = B.rows()
163 int N_; // B.cols()
164 double *A_;
165 double *C_;
168
169 mdivide_left_tri_vd_vari(const Eigen::Matrix<var, R1, C1> &A,
170 const Eigen::Matrix<double, R2, C2> &B)
171 : vari(0.0),
172 M_(A.rows()),
173 N_(B.cols()),
174 A_(reinterpret_cast<double *>(
175 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * A.rows()
176 * A.cols()))),
177 C_(reinterpret_cast<double *>(
178 ChainableStack::instance_->memalloc_.alloc(sizeof(double) * B.rows()
179 * B.cols()))),
180 variRefA_(reinterpret_cast<vari **>(
181 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * A.rows()
182 * (A.rows() + 1) / 2))),
183 variRefC_(reinterpret_cast<vari **>(
184 ChainableStack::instance_->memalloc_.alloc(sizeof(vari *) * B.rows()
185 * B.cols()))) {
186 using Eigen::Map;
187 using Eigen::Matrix;
188
189 size_t pos = 0;
190 if (TriView == Eigen::Lower) {
191 for (size_type j = 0; j < M_; j++) {
192 for (size_type i = j; i < M_; i++) {
193 variRefA_[pos++] = A(i, j).vi_;
194 }
195 }
196 } else if (TriView == Eigen::Upper) {
197 for (size_type j = 0; j < M_; j++) {
198 for (size_type i = 0; i < j + 1; i++) {
199 variRefA_[pos++] = A(i, j).vi_;
200 }
201 }
202 }
203 Map<matrix_d> Ad(A_, M_, M_);
204 Map<matrix_d> Cd(C_, M_, N_);
205 Ad = A.val();
206
207 Cd = Ad.template triangularView<TriView>().solve(B);
208
209 Map<matrix_vi>(variRefC_, M_, N_)
210 = Cd.unaryExpr([](double x) { return new vari(x, false); });
211 }
212
213 virtual void chain() {
214 using Eigen::Map;
215 using Eigen::Matrix;
216 Matrix<double, R1, C1> adjA(M_, M_);
217 Matrix<double, R1, C2> adjC(M_, N_);
218
219 adjC = Map<matrix_vi>(variRefC_, M_, N_).adj();
220
221 adjA.noalias()
222 = -Map<Matrix<double, R1, C1>>(A_, M_, M_)
223 .template triangularView<TriView>()
224 .transpose()
225 .solve(adjC
226 * Map<Matrix<double, R1, C2>>(C_, M_, N_).transpose());
227
228 size_t pos = 0;
229 if (TriView == Eigen::Lower) {
230 for (size_type j = 0; j < adjA.cols(); j++) {
231 for (size_type i = j; i < adjA.rows(); i++) {
232 variRefA_[pos++]->adj_ += adjA(i, j);
233 }
234 }
235 } else if (TriView == Eigen::Upper) {
236 for (size_type j = 0; j < adjA.cols(); j++) {
237 for (size_type i = 0; i < j + 1; i++) {
238 variRefA_[pos++]->adj_ += adjA(i, j);
239 }
240 }
241 }
242 }
243};
244} // namespace internal
245
246template <Eigen::UpLoType TriView, typename T1, typename T2,
248inline Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime>
249mdivide_left_tri(const T1 &A, const T2 &b) {
250 check_square("mdivide_left_tri", "A", A);
251 check_multiplicable("mdivide_left_tri", "A", A, "b", b);
252 if (A.rows() == 0) {
253 return {0, b.cols()};
254 }
255
256 // NOTE: this is not a memory leak, this vari is used in the
257 // expression graph to evaluate the adjoint, but is not needed
258 // for the returned matrix. Memory will be cleaned up with the
259 // arena allocator.
260 auto *baseVari = new internal::mdivide_left_tri_vv_vari<
261 TriView, T1::RowsAtCompileTime, T1::ColsAtCompileTime,
262 T2::RowsAtCompileTime, T2::ColsAtCompileTime>(A, b);
263
264 Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime> res(
265 b.rows(), b.cols());
266 res.vi()
267 = Eigen::Map<matrix_vi>(&(baseVari->variRefC_[0]), b.rows(), b.cols());
268
269 return res;
270}
271template <Eigen::UpLoType TriView, typename T1, typename T2,
274inline Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime>
275mdivide_left_tri(const T1 &A, const T2 &b) {
276 check_square("mdivide_left_tri", "A", A);
277 check_multiplicable("mdivide_left_tri", "A", A, "b", b);
278 if (A.rows() == 0) {
279 return {0, b.cols()};
280 }
281
282 // NOTE: this is not a memory leak, this vari is used in the
283 // expression graph to evaluate the adjoint, but is not needed
284 // for the returned matrix. Memory will be cleaned up with the
285 // arena allocator.
286 auto *baseVari = new internal::mdivide_left_tri_dv_vari<
287 TriView, T1::RowsAtCompileTime, T1::ColsAtCompileTime,
288 T2::RowsAtCompileTime, T2::ColsAtCompileTime>(A, b);
289
290 Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime> res(
291 b.rows(), b.cols());
292 res.vi()
293 = Eigen::Map<matrix_vi>(&(baseVari->variRefC_[0]), b.rows(), b.cols());
294
295 return res;
296}
297template <Eigen::UpLoType TriView, typename T1, typename T2,
300inline Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime>
301mdivide_left_tri(const T1 &A, const T2 &b) {
302 check_square("mdivide_left_tri", "A", A);
303 check_multiplicable("mdivide_left_tri", "A", A, "b", b);
304 if (A.rows() == 0) {
305 return {0, b.cols()};
306 }
307
308 // NOTE: this is not a memory leak, this vari is used in the
309 // expression graph to evaluate the adjoint, but is not needed
310 // for the returned matrix. Memory will be cleaned up with the
311 // arena allocator.
312 auto *baseVari = new internal::mdivide_left_tri_vd_vari<
313 TriView, T1::RowsAtCompileTime, T1::ColsAtCompileTime,
314 T2::RowsAtCompileTime, T2::ColsAtCompileTime>(A, b);
315
316 Eigen::Matrix<var, T1::RowsAtCompileTime, T2::ColsAtCompileTime> res(
317 b.rows(), b.cols());
318 res.vi()
319 = Eigen::Map<matrix_vi>(&(baseVari->variRefC_[0]), b.rows(), b.cols());
320
321 return res;
322}
323
343template <Eigen::UpLoType TriView, typename T1, typename T2,
344 require_all_matrix_t<T1, T2> * = nullptr,
345 require_any_var_matrix_t<T1, T2> * = nullptr>
346inline auto mdivide_left_tri(const T1 &A, const T2 &B) {
347 using ret_val_type = plain_type_t<decltype(value_of(A) * value_of(B))>;
348 using ret_type = var_value<ret_val_type>;
349
350 if (A.size() == 0) {
351 return ret_type(ret_val_type(0, B.cols()));
352 }
353
354 check_square("mdivide_left_tri", "A", A);
355 check_multiplicable("mdivide_left_tri", "A", A, "B", B);
356
357 if (!is_constant<T1>::value && !is_constant<T2>::value) {
358 arena_t<promote_scalar_t<var, T1>> arena_A = A;
359 arena_t<promote_scalar_t<var, T2>> arena_B = B;
360 auto arena_A_val = to_arena(arena_A.val());
361
362 arena_t<ret_type> res
363 = arena_A_val.template triangularView<TriView>().solve(arena_B.val());
364
365 reverse_pass_callback([arena_A, arena_B, arena_A_val, res]() mutable {
366 promote_scalar_t<double, T2> adjB
367 = arena_A_val.template triangularView<TriView>().transpose().solve(
368 res.adj());
369
370 arena_B.adj() += adjB;
371 arena_A.adj() -= (adjB * res.val().transpose().eval())
372 .template triangularView<TriView>();
373 });
374
375 return ret_type(res);
376 } else if (!is_constant<T1>::value) {
377 arena_t<promote_scalar_t<var, T1>> arena_A = A;
378 auto arena_A_val = to_arena(arena_A.val());
379
380 arena_t<ret_type> res
381 = arena_A_val.template triangularView<TriView>().solve(value_of(B));
382
383 reverse_pass_callback([arena_A, arena_A_val, res]() mutable {
384 promote_scalar_t<double, T2> adjB
385 = arena_A_val.template triangularView<TriView>().transpose().solve(
386 res.adj());
387
388 arena_A.adj() -= (adjB * res.val().transpose().eval())
389 .template triangularView<TriView>();
390 });
391
392 return ret_type(res);
393 } else {
394 arena_t<promote_scalar_t<double, T1>> arena_A = value_of(A);
395 arena_t<promote_scalar_t<var, T2>> arena_B = B;
396
397 arena_t<ret_type> res
398 = arena_A.template triangularView<TriView>().solve(arena_B.val());
399
400 reverse_pass_callback([arena_A, arena_B, res]() mutable {
401 promote_scalar_t<double, T2> adjB
402 = arena_A.template triangularView<TriView>().transpose().solve(
403 res.adj());
404
405 arena_B.adj() += adjB;
406 });
407
408 return ret_type(res);
409 }
410}
411
412} // namespace math
413} // namespace stan
414#endif
mdivide_left_tri_dv_vari(const Eigen::Matrix< double, R1, C1 > &A, const Eigen::Matrix< var, R2, C2 > &B)
mdivide_left_tri_vd_vari(const Eigen::Matrix< var, R1, C1 > &A, const Eigen::Matrix< double, R2, C2 > &B)
mdivide_left_tri_vv_vari(const Eigen::Matrix< var, R1, C1 > &A, const Eigen::Matrix< var, R2, C2 > &B)
require_all_t< container_type_check_base< is_eigen, value_type_t, TypeCheck, Check >... > require_all_eigen_vt
Require all of the types satisfy is_eigen.
Definition is_eigen.hpp:188
require_t< container_type_check_base< is_eigen, value_type_t, TypeCheck, Check... > > require_eigen_vt
Require type satisfies is_eigen.
Definition is_eigen.hpp:152
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_square(const char *function, const char *name, const T_y &y)
Check if the specified matrix is square.
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.
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
vari_value< double > vari
Definition vari.hpp:197
arena_t< T > to_arena(const T &a)
Converts given argument into a type that either has any dynamic allocation on AD stack or schedules i...
Definition to_arena.hpp:25
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
Type for sizes and indexes in an Eigen matrix with double elements.
Definition typedefs.hpp:11
auto mdivide_left_tri(const T1 &A, const T2 &b)
Returns the solution of the system Ax=b when A is triangular.
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
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
This struct always provides access to the autodiff stack using the singleton pattern.