1 #ifndef STAN_MATH_REV_MAT_FUN_MDIVIDE_LEFT_HPP 2 #define STAN_MATH_REV_MAT_FUN_MDIVIDE_LEFT_HPP 16 template <
int R1,
int C1,
int R2,
int C2>
17 class mdivide_left_vv_vari :
public vari {
27 mdivide_left_vv_vari(
const Eigen::Matrix<var, R1, C1> &A,
28 const Eigen::Matrix<var, R2, C2> &B)
32 A_(reinterpret_cast<double*>
34 .alloc(sizeof(double) * A.
rows() * A.
cols()))),
35 C_(reinterpret_cast<double*>
37 .alloc(sizeof(double) * B.
rows() * B.
cols()))),
38 variRefA_(reinterpret_cast<vari**>
40 .alloc(sizeof(vari*) * A.
rows() * A.
cols()))),
41 variRefB_(reinterpret_cast<vari**>
43 .alloc(sizeof(vari*) * B.
rows() * B.
cols()))),
44 variRefC_(reinterpret_cast<vari**>
46 .alloc(sizeof(vari*) * B.
rows() * B.
cols()))) {
53 variRefA_[pos] = A(i, j).vi_;
54 A_[pos++] = A(i, j).val();
61 variRefB_[pos] = B(i, j).vi_;
62 C_[pos++] = B(i, j).val();
66 Matrix<double, R1, C2> C(M_, N_);
67 C = Map<Matrix<double, R1, C2> >(
C_,
M_,
N_);
69 C = Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
70 .colPivHouseholderQr().solve(C);
76 variRefC_[pos] =
new vari(C_[pos],
false);
82 virtual void chain() {
85 Eigen::Matrix<double, R1, C1> adjA(M_, M_);
86 Eigen::Matrix<double, R2, C2> adjB(M_, N_);
87 Eigen::Matrix<double, R1, C2> adjC(M_, N_);
90 for (
size_type j = 0; j < adjC.cols(); j++)
91 for (
size_type i = 0; i < adjC.rows(); i++)
92 adjC(i, j) = variRefC_[pos++]->adj_;
94 adjB = Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
95 .
transpose().colPivHouseholderQr().solve(adjC);
96 adjA.noalias() = -adjB
100 for (
size_type j = 0; j < adjA.cols(); j++)
101 for (
size_type i = 0; i < adjA.rows(); i++)
102 variRefA_[pos++]->adj_ += adjA(i, j);
105 for (
size_type j = 0; j < adjB.cols(); j++)
106 for (
size_type i = 0; i < adjB.rows(); i++)
107 variRefB_[pos++]->adj_ += adjB(i, j);
111 template <
int R1,
int C1,
int R2,
int C2>
112 class mdivide_left_dv_vari :
public vari {
121 mdivide_left_dv_vari(
const Eigen::Matrix<double, R1, C1> &A,
122 const Eigen::Matrix<var, R2, C2> &B)
126 A_(reinterpret_cast<double*>
128 .alloc(
sizeof(
double) * A.rows() * A.cols()))),
129 C_(reinterpret_cast<double*>
131 .alloc(
sizeof(
double) * B.rows() * B.cols()))),
134 .alloc(
sizeof(vari*) * B.rows() * B.cols()))),
137 .alloc(
sizeof(vari*) * B.rows() * B.cols()))) {
151 variRefB_[pos] = B(i, j).vi_;
152 C_[pos++] = B(i, j).val();
156 Matrix<double, R1, C2> C(M_, N_);
157 C = Map<Matrix<double, R1, C2> >(
C_,
M_,
N_);
159 C = Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
160 .colPivHouseholderQr().solve(C);
166 variRefC_[pos] =
new vari(C_[pos],
false);
172 virtual void chain() {
175 Eigen::Matrix<double, R2, C2> adjB(M_, N_);
176 Eigen::Matrix<double, R1, C2> adjC(M_, N_);
179 for (
size_type j = 0; j < adjC.cols(); j++)
180 for (
size_type i = 0; i < adjC.rows(); i++)
181 adjC(i, j) = variRefC_[pos++]->adj_;
183 adjB = Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
184 .
transpose().colPivHouseholderQr().solve(adjC);
187 for (
size_type j = 0; j < adjB.cols(); j++)
188 for (
size_type i = 0; i < adjB.rows(); i++)
189 variRefB_[pos++]->adj_ += adjB(i, j);
193 template <
int R1,
int C1,
int R2,
int C2>
194 class mdivide_left_vd_vari :
public vari {
203 mdivide_left_vd_vari(
const Eigen::Matrix<var, R1, C1> &A,
204 const Eigen::Matrix<double, R2, C2> &B)
208 A_(reinterpret_cast<double*>
210 .alloc(
sizeof(
double) * A.rows() * A.cols()))),
211 C_(reinterpret_cast<double*>
213 .alloc(
sizeof(
double) * B.rows() * B.cols()))),
216 .alloc(
sizeof(vari*) * A.rows() * A.cols()))),
219 .alloc(
sizeof(vari*) * B.rows() * B.cols()))) {
226 variRefA_[pos] = A(i, j).vi_;
227 A_[pos++] = A(i, j).val();
231 Matrix<double, R1, C2> C(M_, N_);
232 C = Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
233 .colPivHouseholderQr().solve(B);
239 variRefC_[pos] =
new vari(C_[pos],
false);
245 virtual void chain() {
248 Eigen::Matrix<double, R1, C1> adjA(M_, M_);
249 Eigen::Matrix<double, R1, C2> adjC(M_, N_);
252 for (
size_type j = 0; j < adjC.cols(); j++)
253 for (
size_type i = 0; i < adjC.rows(); i++)
254 adjC(i, j) = variRefC_[pos++]->adj_;
257 adjA = -Map<Matrix<double, R1, C1> >(
A_,
M_,
M_)
259 .colPivHouseholderQr()
260 .solve(adjC*Map<Matrix<double, R1, C2> >(C_, M_, N_).
transpose());
263 for (
size_type j = 0; j < adjA.cols(); j++)
264 for (
size_type i = 0; i < adjA.rows(); i++)
265 variRefA_[pos++]->adj_ += adjA(i, j);
270 template <
int R1,
int C1,
int R2,
int C2>
272 Eigen::Matrix<var, R1, C2>
274 const Eigen::Matrix<var, R2, C2> &b) {
275 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
286 mdivide_left_vv_vari<R1, C1, R2, C2> *baseVari
287 =
new mdivide_left_vv_vari<R1, C1, R2, C2>(A, b);
290 for (
size_type j = 0; j < res.cols(); j++)
291 for (
size_type i = 0; i < res.rows(); i++)
292 res(i, j).vi_ = baseVari->variRefC_[pos++];
297 template <
int R1,
int C1,
int R2,
int C2>
299 Eigen::Matrix<var, R1, C2>
301 const Eigen::Matrix<double, R2, C2> &b) {
302 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
313 mdivide_left_vd_vari<R1, C1, R2, C2> *baseVari
314 =
new mdivide_left_vd_vari<R1, C1, R2, C2>(A, b);
317 for (
size_type j = 0; j < res.cols(); j++)
318 for (
size_type i = 0; i < res.rows(); i++)
319 res(i, j).vi_ = baseVari->variRefC_[pos++];
324 template <
int R1,
int C1,
int R2,
int C2>
326 Eigen::Matrix<var, R1, C2>
328 const Eigen::Matrix<var, R2, C2> &b) {
329 Eigen::Matrix<var, R1, C2> res(b.rows(), b.cols());
340 mdivide_left_dv_vari<R1, C1, R2, C2> *baseVari
341 =
new mdivide_left_dv_vari<R1, C1, R2, C2>(A, b);
344 for (
size_type j = 0; j < res.cols(); j++)
345 for (
size_type i = 0; i < res.rows(); i++)
346 res(i, j).vi_ = baseVari->variRefC_[pos++];
Eigen::Matrix< fvar< T >, R1, C2 > mdivide_left(const Eigen::Matrix< fvar< T >, R1, C1 > &A, const Eigen::Matrix< fvar< T >, R2, C2 > &b)
int rows(const Eigen::Matrix< T, R, C > &m)
Return the number of rows in the specified matrix, vector, or row vector.
static stack_alloc memalloc_
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
Type for sizes and indexes in an Eigen matrix with double e.
int cols(const Eigen::Matrix< T, R, C > &m)
Return the number of columns in the specified matrix, vector, or row vector.
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 check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Check if the specified matrix is square.
Eigen::Matrix< T, C, R > transpose(const Eigen::Matrix< T, R, C > &m)
AutodiffStackStorage< vari, chainable_alloc > ChainableStack