Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mdivide_left.hpp
Go to the documentation of this file.
1 #ifndef STAN__AGRAD__REV__MATRIX__COLUMNS_MDIVIDE_LEFT_HPP
2 #define STAN__AGRAD__REV__MATRIX__COLUMNS_MDIVIDE_LEFT_HPP
3 
4 #include <vector>
8 #include <stan/agrad/rev/var.hpp>
11 
12 namespace stan {
13  namespace agrad {
14 
15  namespace {
16  template <int R1,int C1,int R2,int C2>
17  class mdivide_left_vv_vari : public vari {
18  public:
19  int M_; // A.rows() = A.cols() = B.rows()
20  int N_; // B.cols()
21  double* A_;
22  double* C_;
23  vari** _variRefA;
24  vari** _variRefB;
25  vari** _variRefC;
26 
27  mdivide_left_vv_vari(const Eigen::Matrix<var,R1,C1> &A,
28  const Eigen::Matrix<var,R2,C2> &B)
29  : vari(0.0),
30  M_(A.rows()),
31  N_(B.cols()),
32  A_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
33  * A.rows() * A.cols())),
34  C_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
35  * B.rows() * B.cols())),
36  _variRefA((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
37  * A.rows() * A.cols())),
38  _variRefB((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
39  * B.rows() * B.cols())),
40  _variRefC((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
41  * B.rows() * B.cols()))
42  {
43  using Eigen::Matrix;
44  using Eigen::Map;
45 
46  size_t pos = 0;
47  for (size_type j = 0; j < M_; j++) {
48  for (size_type i = 0; i < M_; i++) {
49  _variRefA[pos] = A(i,j).vi_;
50  A_[pos++] = A(i,j).val();
51  }
52  }
53 
54  pos = 0;
55  for (size_type j = 0; j < N_; j++) {
56  for (size_type i = 0; i < M_; i++) {
57  _variRefB[pos] = B(i,j).vi_;
58  C_[pos++] = B(i,j).val();
59  }
60  }
61 
62  Matrix<double,R1,C2> C(M_,N_);
63  C = Map<Matrix<double,R1,C2> >(C_,M_,N_);
64 
65  C = Map<Matrix<double,R1,C1> >(A_,M_,M_)
66  .colPivHouseholderQr().solve(C);
67 
68  pos = 0;
69  for (size_type j = 0; j < N_; j++) {
70  for (size_type i = 0; i < M_; i++) {
71  C_[pos] = C(i,j);
72  _variRefC[pos] = new vari(C_[pos],false);
73  pos++;
74  }
75  }
76  }
77 
78  virtual void chain() {
79  using Eigen::Matrix;
80  using Eigen::Map;
81  Eigen::Matrix<double,R1,C1> adjA(M_,M_);
82  Eigen::Matrix<double,R2,C2> adjB(M_,N_);
83  Eigen::Matrix<double,R1,C2> adjC(M_,N_);
84 
85  size_t pos = 0;
86  for (size_type j = 0; j < adjC.cols(); j++)
87  for (size_type i = 0; i < adjC.rows(); i++)
88  adjC(i,j) = _variRefC[pos++]->adj_;
89 
90 
91  adjB = Map<Matrix<double,R1,C1> >(A_,M_,M_)
92  .transpose().colPivHouseholderQr().solve(adjC);
93  adjA.noalias() = -adjB
94  * Map<Matrix<double,R1,C2> >(C_,M_,N_).transpose();
95 
96  pos = 0;
97  for (size_type j = 0; j < adjA.cols(); j++)
98  for (size_type i = 0; i < adjA.rows(); i++)
99  _variRefA[pos++]->adj_ += adjA(i,j);
100 
101  pos = 0;
102  for (size_type j = 0; j < adjB.cols(); j++)
103  for (size_type i = 0; i < adjB.rows(); i++)
104  _variRefB[pos++]->adj_ += adjB(i,j);
105  }
106  };
107 
108  template <int R1,int C1,int R2,int C2>
109  class mdivide_left_dv_vari : public vari {
110  public:
111  int M_; // A.rows() = A.cols() = B.rows()
112  int N_; // B.cols()
113  double* A_;
114  double* C_;
115  vari** _variRefB;
116  vari** _variRefC;
117 
118  mdivide_left_dv_vari(const Eigen::Matrix<double,R1,C1> &A,
119  const Eigen::Matrix<var,R2,C2> &B)
120  : vari(0.0),
121  M_(A.rows()),
122  N_(B.cols()),
123  A_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
124  * A.rows() * A.cols())),
125  C_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
126  * B.rows() * B.cols())),
127  _variRefB((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
128  * B.rows() * B.cols())),
129  _variRefC((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
130  * B.rows() * B.cols()))
131  {
132  using Eigen::Matrix;
133  using Eigen::Map;
134 
135  size_t pos = 0;
136  for (size_type j = 0; j < M_; j++) {
137  for (size_type i = 0; i < M_; i++) {
138  A_[pos++] = A(i,j);
139  }
140  }
141 
142  pos = 0;
143  for (size_type j = 0; j < N_; j++) {
144  for (size_type i = 0; i < M_; i++) {
145  _variRefB[pos] = B(i,j).vi_;
146  C_[pos++] = B(i,j).val();
147  }
148  }
149 
150  Matrix<double,R1,C2> C(M_,N_);
151  C = Map<Matrix<double,R1,C2> >(C_,M_,N_);
152 
153  C = Map<Matrix<double,R1,C1> >(A_,M_,M_)
154  .colPivHouseholderQr().solve(C);
155 
156  pos = 0;
157  for (size_type j = 0; j < N_; j++) {
158  for (size_type i = 0; i < M_; i++) {
159  C_[pos] = C(i,j);
160  _variRefC[pos] = new vari(C_[pos],false);
161  pos++;
162  }
163  }
164  }
165 
166  virtual void chain() {
167  using Eigen::Matrix;
168  using Eigen::Map;
169  Eigen::Matrix<double,R2,C2> adjB(M_,N_);
170  Eigen::Matrix<double,R1,C2> adjC(M_,N_);
171 
172  size_t pos = 0;
173  for (size_type j = 0; j < adjC.cols(); j++)
174  for (size_type i = 0; i < adjC.rows(); i++)
175  adjC(i,j) = _variRefC[pos++]->adj_;
176 
177  adjB = Map<Matrix<double,R1,C1> >(A_,M_,M_)
178  .transpose().colPivHouseholderQr().solve(adjC);
179 
180  pos = 0;
181  for (size_type j = 0; j < adjB.cols(); j++)
182  for (size_type i = 0; i < adjB.rows(); i++)
183  _variRefB[pos++]->adj_ += adjB(i,j);
184  }
185  };
186 
187  template <int R1,int C1,int R2,int C2>
188  class mdivide_left_vd_vari : public vari {
189  public:
190  int M_; // A.rows() = A.cols() = B.rows()
191  int N_; // B.cols()
192  double* A_;
193  double* C_;
194  vari** _variRefA;
195  vari** _variRefC;
196 
197  mdivide_left_vd_vari(const Eigen::Matrix<var,R1,C1> &A,
198  const Eigen::Matrix<double,R2,C2> &B)
199  : vari(0.0),
200  M_(A.rows()),
201  N_(B.cols()),
202  A_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
203  * A.rows() * A.cols())),
204  C_((double*)stan::agrad::memalloc_.alloc(sizeof(double)
205  * B.rows() * B.cols())),
206  _variRefA((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
207  * A.rows() * A.cols())),
208  _variRefC((vari**)stan::agrad::memalloc_.alloc(sizeof(vari*)
209  * B.rows() * B.cols()))
210  {
211  using Eigen::Matrix;
212  using Eigen::Map;
213 
214  size_t pos = 0;
215  for (size_type j = 0; j < M_; j++) {
216  for (size_type i = 0; i < M_; i++) {
217  _variRefA[pos] = A(i,j).vi_;
218  A_[pos++] = A(i,j).val();
219  }
220  }
221 
222  Matrix<double,R1,C2> C(M_,N_);
223  C = Map<Matrix<double,R1,C1> >(A_,M_,M_)
224  .colPivHouseholderQr().solve(B);
225 
226  pos = 0;
227  for (size_type j = 0; j < N_; j++) {
228  for (size_type i = 0; i < M_; i++) {
229  C_[pos] = C(i,j);
230  _variRefC[pos] = new vari(C_[pos],false);
231  pos++;
232  }
233  }
234  }
235 
236  virtual void chain() {
237  using Eigen::Matrix;
238  using Eigen::Map;
239  Eigen::Matrix<double,R1,C1> adjA(M_,M_);
240  Eigen::Matrix<double,R1,C2> adjC(M_,N_);
241 
242  size_t pos = 0;
243  for (size_type j = 0; j < adjC.cols(); j++)
244  for (size_type i = 0; i < adjC.rows(); i++)
245  adjC(i,j) = _variRefC[pos++]->adj_;
246 
247  // FIXME: add .noalias() to LHS
248  adjA = -Map<Matrix<double,R1,C1> >(A_,M_,M_)
249  .transpose()
250  .colPivHouseholderQr()
251  .solve(adjC*Map<Matrix<double,R1,C2> >(C_,M_,N_).transpose());
252 
253  pos = 0;
254  for (size_type j = 0; j < adjA.cols(); j++)
255  for (size_type i = 0; i < adjA.rows(); i++)
256  _variRefA[pos++]->adj_ += adjA(i,j);
257  }
258  };
259  }
260 
261  template <int R1,int C1,int R2,int C2>
262  inline
263  Eigen::Matrix<var,R1,C2>
264  mdivide_left(const Eigen::Matrix<var,R1,C1> &A,
265  const Eigen::Matrix<var,R2,C2> &b) {
266  Eigen::Matrix<var,R1,C2> res(b.rows(),b.cols());
267 
268  stan::math::check_square("mdivide_left(%1%)",A,"A",(double*)0);
269  stan::math::check_multiplicable("mdivide_left(%1%)",A,"A",
270  b,"b",(double*)0);
271 
272  // NOTE: this is not a memory leak, this vari is used in the
273  // expression graph to evaluate the adjoint, but is not needed
274  // for the returned matrix. Memory will be cleaned up with the arena allocator.
275  mdivide_left_vv_vari<R1,C1,R2,C2> *baseVari = new mdivide_left_vv_vari<R1,C1,R2,C2>(A,b);
276 
277  size_t pos = 0;
278  for (size_type j = 0; j < res.cols(); j++)
279  for (size_type i = 0; i < res.rows(); i++)
280  res(i,j).vi_ = baseVari->_variRefC[pos++];
281 
282  return res;
283  }
284 
285  template <int R1,int C1,int R2,int C2>
286  inline
287  Eigen::Matrix<var,R1,C2>
288  mdivide_left(const Eigen::Matrix<var,R1,C1> &A,
289  const Eigen::Matrix<double,R2,C2> &b) {
290  Eigen::Matrix<var,R1,C2> res(b.rows(),b.cols());
291 
292  stan::math::check_square("mdivide_left(%1%)",A,"A",(double*)0);
293  stan::math::check_multiplicable("mdivide_left(%1%)",A,"A",
294  b,"b",(double*)0);
295 
296  // NOTE: this is not a memory leak, this vari is used in the
297  // expression graph to evaluate the adjoint, but is not needed
298  // for the returned matrix. Memory will be cleaned up with the arena allocator.
299  mdivide_left_vd_vari<R1,C1,R2,C2> *baseVari = new mdivide_left_vd_vari<R1,C1,R2,C2>(A,b);
300 
301  size_t pos = 0;
302  for (size_type j = 0; j < res.cols(); j++)
303  for (size_type i = 0; i < res.rows(); i++)
304  res(i,j).vi_ = baseVari->_variRefC[pos++];
305 
306  return res;
307  }
308 
309  template <int R1,int C1,int R2,int C2>
310  inline
311  Eigen::Matrix<var,R1,C2>
312  mdivide_left(const Eigen::Matrix<double,R1,C1> &A,
313  const Eigen::Matrix<var,R2,C2> &b) {
314  Eigen::Matrix<var,R1,C2> res(b.rows(),b.cols());
315 
316  stan::math::check_square("mdivide_left(%1%)",A,"A",(double*)0);
317  stan::math::check_multiplicable("mdivide_left(%1%)",A,"A",
318  b,"b",(double*)0);
319 
320  // NOTE: this is not a memory leak, this vari is used in the
321  // expression graph to evaluate the adjoint, but is not needed
322  // for the returned matrix. Memory will be cleaned up with the arena allocator.
323  mdivide_left_dv_vari<R1,C1,R2,C2> *baseVari = new mdivide_left_dv_vari<R1,C1,R2,C2>(A,b);
324 
325  size_t pos = 0;
326  for (size_type j = 0; j < res.cols(); j++)
327  for (size_type i = 0; i < res.rows(); i++)
328  res(i,j).vi_ = baseVari->_variRefC[pos++];
329 
330  return res;
331  }
332 
333  }
334 }
335 #endif
int N_
double * A_
Eigen::Matrix< fvar< T >, R1, C2 > mdivide_left(const Eigen::Matrix< fvar< T >, R1, C1 > &A, const Eigen::Matrix< fvar< T >, R2, C2 > &b)
memory::stack_alloc memalloc_
Definition: var_stack.cpp:16
vari ** _variRefA
bool check_multiplicable(const char *function, const T1 &y1, const char *name1, const T2 &y2, const char *name2, T_result *result)
Eigen::Matrix< T, C, R > transpose(const Eigen::Matrix< T, R, C > &m)
Definition: transpose.hpp:12
vari ** _variRefB
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic >::Index size_type
Definition: typedefs.hpp:14
size_t rows(const Eigen::Matrix< T, R, C > &m)
Definition: rows.hpp:12
double * C_
size_t cols(const Eigen::Matrix< T, R, C > &m)
Definition: cols.hpp:12
vari ** _variRefC
bool check_square(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is square.
int M_

     [ Stan Home Page ] © 2011–2014, Stan Development Team.