Stan  1.0
probability, sampling & optimization
matrix.hpp
Go to the documentation of this file.
1 #ifndef __STAN__AGRAD__MATRIX_HPP__
2 #define __STAN__AGRAD__MATRIX_HPP__
3 
4 #include <stan/math/matrix.hpp>
6 
8 
9 namespace stan {
10  namespace agrad {
12  protected:
16  double dotval_;
17  size_t length_;
18  inline static double eval_gevv(const stan::agrad::var* alpha,
19  const stan::agrad::var* v1, int stride1,
20  const stan::agrad::var* v2, int stride2,
21  size_t length, double *dotprod) {
22  double result = 0;
23  for (size_t i = 0; i < length; i++)
24  result += v1[i*stride1].vi_->val_ * v2[i*stride2].vi_->val_;
25  *dotprod = result;
26  return alpha->vi_->val_ * result;
27  }
28  public:
30  const stan::agrad::var* v1, int stride1,
31  const stan::agrad::var* v2, int stride2, size_t length) :
32  vari(eval_gevv(alpha,v1,stride1,v2,stride2,length,&dotval_)), length_(length) {
33  alpha_ = alpha->vi_;
35  v2_ = v1_ + length_;
36  for (size_t i = 0; i < length_; i++)
37  v1_[i] = v1[i*stride1].vi_;
38  for (size_t i = 0; i < length_; i++)
39  v2_[i] = v2[i*stride2].vi_;
40  }
41  void chain() {
42  const double adj_alpha = adj_ * alpha_->val_;
43  for (size_t i = 0; i < length_; i++) {
44  v1_[i]->adj_ += adj_alpha * v2_[i]->val_;
45  v2_[i]->adj_ += adj_alpha * v1_[i]->val_;
46  }
47  alpha_->adj_ += adj_ * dotval_;
48  }
49  };
50  }
51 }
52 
53 namespace Eigen {
54 
59  template <> struct NumTraits<stan::agrad::var>
60  {
67 
74 
81 
88  inline static Real epsilon() {
90  }
91 
95  inline static Real dummy_precision() {
96  return 1e-12; // copied from NumTraits.h values for double
97  }
98 
105  inline static Real highest() {
107  }
108 
115  inline static Real lowest() {
117  }
118 
123  enum {
124  IsInteger = 0,
125  IsSigned = 1,
126  IsComplex = 0,
127  RequireInitialization = 0,
128  ReadCost = 1,
129  AddCost = 1,
130  MulCost = 1,
131  HasFloatingPoint = 1,
132  };
133  };
134 
135  namespace internal {
139  template<>
140  struct significant_decimals_default_impl<stan::agrad::var,false>
141  {
142  static inline int run()
143  {
144  using std::ceil;
145  return cast<double,int>(ceil(-log(NumTraits<stan::agrad::var>::epsilon().val())
146  /log(10.0)));
147  }
148  };
149 
154  template <>
155  struct scalar_product_traits<stan::agrad::var,double> {
157  };
158 
163  template <>
164  struct scalar_product_traits<double,stan::agrad::var> {
166  };
167 
171  template<typename Index, bool ConjugateLhs, bool ConjugateRhs>
172  struct general_matrix_vector_product<Index,stan::agrad::var,ColMajor,ConjugateLhs,stan::agrad::var,ConjugateRhs>
173  {
176  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
177  enum { LhsStorageOrder = ColMajor };
178 
179  EIGEN_DONT_INLINE static void run(
180  Index rows, Index cols,
181  const LhsScalar* lhs, Index lhsStride,
182  const RhsScalar* rhs, Index rhsIncr,
183  ResScalar* res, Index resIncr, const ResScalar &alpha)
184  {
185  for (Index i = 0; i < rows; i++) {
186  res[i*resIncr] += stan::agrad::var(new stan::agrad::gevv_vvv_vari(&alpha,((int)LhsStorageOrder == (int)ColMajor)?(&lhs[i]):(&lhs[i*lhsStride]),((int)LhsStorageOrder == (int)ColMajor)?(lhsStride):(1),rhs,rhsIncr,cols));
187  }
188  }
189  };
190  template<typename Index, bool ConjugateLhs, bool ConjugateRhs>
191  struct general_matrix_vector_product<Index,stan::agrad::var,RowMajor,ConjugateLhs,stan::agrad::var,ConjugateRhs>
192  {
195  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
196  enum { LhsStorageOrder = RowMajor };
197 
198  EIGEN_DONT_INLINE static void run(
199  Index rows, Index cols,
200  const LhsScalar* lhs, Index lhsStride,
201  const RhsScalar* rhs, Index rhsIncr,
202  ResScalar* res, Index resIncr, const RhsScalar &alpha)
203  {
204  for (Index i = 0; i < rows; i++) {
205  res[i*resIncr] += stan::agrad::var(new stan::agrad::gevv_vvv_vari(&alpha,((int)LhsStorageOrder == (int)ColMajor)?(&lhs[i]):(&lhs[i*lhsStride]),((int)LhsStorageOrder == (int)ColMajor)?(lhsStride):(1),rhs,rhsIncr,cols));
206  }
207  }
208  };
209  template<typename Index, int LhsStorageOrder, bool ConjugateLhs, int RhsStorageOrder, bool ConjugateRhs>
210  struct general_matrix_matrix_product<Index,stan::agrad::var,LhsStorageOrder,ConjugateLhs,stan::agrad::var,RhsStorageOrder,ConjugateRhs,ColMajor>
211  {
214  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
215  static void run(Index rows, Index cols, Index depth,
216  const LhsScalar* _lhs, Index lhsStride,
217  const RhsScalar* _rhs, Index rhsStride,
218  ResScalar* res, Index resStride,
219  const ResScalar &alpha,
220  level3_blocking<LhsScalar,RhsScalar>& blocking,
221  GemmParallelInfo<Index>* info = 0)
222  {
223  for (Index i = 0; i < cols; i++) {
224  general_matrix_vector_product<Index,LhsScalar,LhsStorageOrder,ConjugateLhs,RhsScalar,ConjugateRhs>::run(
225  rows,depth,_lhs,lhsStride,
226  &_rhs[((int)RhsStorageOrder == (int)ColMajor)?(i*rhsStride):(i)],((int)RhsStorageOrder == (int)ColMajor)?(1):(rhsStride),
227  &res[i*resStride],1,alpha);
228  }
229  }
230  };
231  }
232 }
233 
234 namespace stan {
235 
236  namespace agrad {
237 
242  typedef
243  Eigen::Matrix<var,Eigen::Dynamic,Eigen::Dynamic>
245 
250  typedef
251  Eigen::Matrix<var,Eigen::Dynamic,1>
253 
258  typedef
259  Eigen::Matrix<var,1,Eigen::Dynamic>
261 
266  inline void initialize_variable(var& variable, const var& value) {
267  variable = value;
268  }
269 
274  template <int R, int C>
275  inline void initialize_variable(Eigen::Matrix<var,R,C>& matrix, const var& value) {
276  for (int i = 0; i < matrix.size(); ++i)
277  matrix(i) = value;
278  }
279 
283  template <typename T>
284  inline void initialize_variable(std::vector<T>& variables, const var& value) {
285  for (size_t i = 0; i < variables.size(); ++i)
286  initialize_variable(variables[i],value);
287  }
288 
289 
298  inline var to_var(const double& x) {
299  return var(x);
300  }
309  inline var to_var(const var& x) {
310  return x;
311  }
321  matrix_v m_v(m.rows(), m.cols());
322  for (int j = 0; j < m.cols(); ++j)
323  for (int i = 0; i < m.rows(); ++i)
324  m_v(i,j) = m(i,j);
325  return m_v;
326  }
335  inline matrix_v to_var(const matrix_v& m) {
336  return m;
337  }
348  vector_v v_v(v.size());
349  for (int i = 0; i < v.size(); ++i)
350  v_v[i] = v[i];
351  return v_v;
352  }
362  inline vector_v to_var(const vector_v& v) {
363  return v;
364  }
375  row_vector_v rv_v(rv.size());
376  for (int i = 0; i < rv.size(); ++i)
377  rv_v[i] = rv[i];
378  return rv_v;
379  }
389  inline row_vector_v to_var(const row_vector_v& rv) {
390  return rv;
391  }
392 
393 
394 
395  // scalar returns
396 
397  // /**
398  // * Determinant of the matrix.
399  // *
400  // * Returns the determinant of the specified
401  // * square matrix.
402  // *
403  // * @param m Specified matrix.
404  // * @return Determinant of the matrix.
405  // * @throw std::domain_error if m is not a square matrix
406  // */
407  // var determinant(const Eigen::Matrix<var, Eigen::Dynamic, Eigen::Dynamic>& m);
408 
409  namespace {
410  class dot_self_vari : public vari {
411  protected:
412  vari** v_;
413  size_t size_;
414  public:
415  dot_self_vari(vari** v, size_t size)
416  : vari(var_dot_self(v,size)),
417  v_(v),
418  size_(size) {
419  }
420  template <int R, int C>
421  dot_self_vari(const Eigen::Matrix<var,R,C>& v) :
422  vari(var_dot_self(v)), size_(v.size()) {
423  v_ = (vari**) memalloc_.alloc(size_ * sizeof(vari*));
424  for (size_t i = 0; i < size_; ++i)
425  v_[i] = v(i).vi_;
426  }
427  inline static double square(double x) { return x * x; }
428  inline static double var_dot_self(vari** v, size_t size) {
429  double sum = 0.0;
430  for (size_t i = 0; i < size; ++i)
431  sum += square(v[i]->val_);
432  return sum;
433  }
434  template <int R, int C>
435  inline static double var_dot_self(const Eigen::Matrix<var,R,C> &v) {
436  double sum = 0.0;
437  for (int i = 0; i < v.size(); ++i)
438  sum += square(v(i).vi_->val_);
439  return sum;
440  }
441  void chain() {
442  for (size_t i = 0; i < size_; ++i)
443  v_[i]->adj_ += adj_ * 2.0 * v_[i]->val_;
444  }
445  };
446 
447  class sum_v_vari : public vari{
448  protected:
449  vari** v_;
450  size_t length_;
451  inline static double var_sum(const var *v, size_t length) {
452  double result = 0;
453  for (size_t i = 0; i < length; i++)
454  result += v[i].vi_->val_;
455  return result;
456  }
457  template<typename Derived>
458  inline static double var_sum(const Eigen::DenseBase<Derived> &v) {
459  double result = 0;
460  for (int i = 0; i < v.size(); i++)
461  result += v(i).vi_->val_;
462  return result;
463  }
464  public:
465  template<typename Derived>
466  sum_v_vari(const Eigen::DenseBase<Derived> &v) :
467  vari(var_sum(v)), length_(v.size()) {
468  v_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
469  for (size_t i = 0; i < length_; i++)
470  v_[i] = v(i).vi_;
471  }
472  template<int R1,int C1>
473  sum_v_vari(const Eigen::Matrix<var,R1,C1> &v1) :
474  vari(var_sum(v1)), length_(v1.size()) {
475  v_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
476  for (size_t i = 0; i < length_; i++)
477  v_[i] = v1(i).vi_;
478  }
479  sum_v_vari(const var *v, size_t len) :
480  vari(var_sum(v,len)), length_(len) {
481  v_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
482  for (size_t i = 0; i < length_; i++)
483  v_[i] = v[i].vi_;
484  }
485  void chain() {
486  for (size_t i = 0; i < length_; i++) {
487  v_[i]->adj_ += adj_;
488  }
489  }
490  };
491  class dot_product_vv_vari : public vari {
492  protected:
493  vari** v1_;
494  vari** v2_;
495  size_t length_;
496  inline static double var_dot(const var* v1, const var* v2,
497  size_t length) {
498  double result = 0;
499  for (size_t i = 0; i < length; i++)
500  result += v1[i].vi_->val_ * v2[i].vi_->val_;
501  return result;
502  }
503  template<typename Derived1,typename Derived2>
504  inline static double var_dot(const Eigen::DenseBase<Derived1> &v1,
505  const Eigen::DenseBase<Derived2> &v2) {
506  double result = 0;
507  for (int i = 0; i < v1.size(); i++)
508  result += v1[i].vi_->val_ * v2[i].vi_->val_;
509  return result;
510  }
511  inline static double var_dot(vari** v1, vari** v2, size_t length) {
512  double result = 0;
513  for (size_t i = 0; i < length; ++i)
514  result += v1[i]->val_ * v2[i]->val_;
515  return result;
516  }
517  public:
518  dot_product_vv_vari(vari** v1, vari** v2, size_t length)
519  : vari(var_dot(v1,v2,length)),
520  v1_(v1),
521  v2_(v2),
522  length_(length) {
523 
524  }
525  dot_product_vv_vari(const var* v1, const var* v2, size_t length,
526  dot_product_vv_vari* shared_v1 = NULL,
527  dot_product_vv_vari* shared_v2 = NULL) :
528  vari(var_dot(v1, v2, length)), length_(length) {
529  if (shared_v1 == NULL) {
530  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
531  for (size_t i = 0; i < length_; i++)
532  v1_[i] = v1[i].vi_;
533  }
534  else {
535  v1_ = shared_v1->v1_;
536  }
537  if (shared_v2 == NULL) {
538  v2_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
539  for (size_t i = 0; i < length_; i++)
540  v2_[i] = v2[i].vi_;
541  }
542  else {
543  v2_ = shared_v2->v2_;
544  }
545  }
546  template<typename Derived1,typename Derived2>
547  dot_product_vv_vari(const Eigen::DenseBase<Derived1> &v1,
548  const Eigen::DenseBase<Derived2> &v2,
549  dot_product_vv_vari* shared_v1 = NULL,
550  dot_product_vv_vari* shared_v2 = NULL) :
551  vari(var_dot(v1, v2)), length_(v1.size()) {
552  if (shared_v1 == NULL) {
553  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
554  for (size_t i = 0; i < length_; i++)
555  v1_[i] = v1[i].vi_;
556  }
557  else {
558  v1_ = shared_v1->v1_;
559  }
560  if (shared_v2 == NULL) {
561  v2_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
562  for (size_t i = 0; i < length_; i++)
563  v2_[i] = v2[i].vi_;
564  }
565  else {
566  v2_ = shared_v2->v2_;
567  }
568  }
569  template<int R1,int C1,int R2,int C2>
570  dot_product_vv_vari(const Eigen::Matrix<var,R1,C1> &v1,
571  const Eigen::Matrix<var,R2,C2> &v2,
572  dot_product_vv_vari* shared_v1 = NULL,
573  dot_product_vv_vari* shared_v2 = NULL) :
574  vari(var_dot(v1, v2)), length_(v1.size()) {
575  if (shared_v1 == NULL) {
576  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
577  for (size_t i = 0; i < length_; i++)
578  v1_[i] = v1[i].vi_;
579  }
580  else {
581  v1_ = shared_v1->v1_;
582  }
583  if (shared_v2 == NULL) {
584  v2_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
585  for (size_t i = 0; i < length_; i++)
586  v2_[i] = v2[i].vi_;
587  }
588  else {
589  v2_ = shared_v2->v2_;
590  }
591  }
592  void chain() {
593  for (size_t i = 0; i < length_; i++) {
594  v1_[i]->adj_ += adj_ * v2_[i]->val_;
595  v2_[i]->adj_ += adj_ * v1_[i]->val_;
596  }
597  }
598  };
599 
600  class dot_product_vd_vari : public vari {
601  protected:
602  vari** v1_;
603  double* v2_;
604  size_t length_;
605  inline static double var_dot(const var* v1, const double* v2,
606  size_t length) {
607  double result = 0;
608  for (size_t i = 0; i < length; i++)
609  result += v1[i].vi_->val_ * v2[i];
610  return result;
611  }
612  template<typename Derived1,typename Derived2>
613  inline static double var_dot(const Eigen::DenseBase<Derived1> &v1,
614  const Eigen::DenseBase<Derived2> &v2) {
615  double result = 0;
616  for (int i = 0; i < v1.size(); i++)
617  result += v1[i].vi_->val_ * v2[i];
618  return result;
619  }
620  public:
621  dot_product_vd_vari(const var* v1, const double* v2, size_t length,
622  dot_product_vd_vari *shared_v1 = NULL,
623  dot_product_vd_vari *shared_v2 = NULL) :
624  vari(var_dot(v1, v2, length)), length_(length) {
625  if (shared_v1 == NULL) {
626  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
627  for (size_t i = 0; i < length_; i++)
628  v1_[i] = v1[i].vi_;
629  } else {
630  v1_ = shared_v1->v1_;
631  }
632  if (shared_v2 == NULL) {
633  v2_ = (double*)memalloc_.alloc(length_*sizeof(double));
634  for (size_t i = 0; i < length_; i++)
635  v2_[i] = v2[i];
636  } else {
637  v2_ = shared_v2->v2_;
638  }
639  }
640  template<typename Derived1,typename Derived2>
641  dot_product_vd_vari(const Eigen::DenseBase<Derived1> &v1,
642  const Eigen::DenseBase<Derived2> &v2,
643  dot_product_vd_vari *shared_v1 = NULL,
644  dot_product_vd_vari *shared_v2 = NULL) :
645  vari(var_dot(v1, v2)), length_(v1.size()) {
646  if (shared_v1 == NULL) {
647  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
648  for (size_t i = 0; i < length_; i++)
649  v1_[i] = v1[i].vi_;
650  } else {
651  v1_ = shared_v1->v1_;
652  }
653  if (shared_v2 == NULL) {
654  v2_ = (double*)memalloc_.alloc(length_*sizeof(double));
655  for (size_t i = 0; i < length_; i++)
656  v2_[i] = v2[i];
657  } else {
658  v2_ = shared_v2->v2_;
659  }
660  }
661  template<int R1,int C1,int R2,int C2>
662  dot_product_vd_vari(const Eigen::Matrix<var,R1,C1> &v1,
663  const Eigen::Matrix<double,R2,C2> &v2,
664  dot_product_vd_vari *shared_v1 = NULL,
665  dot_product_vd_vari *shared_v2 = NULL) :
666  vari(var_dot(v1, v2)), length_(v1.size()) {
667  if (shared_v1 == NULL) {
668  v1_ = (vari**)memalloc_.alloc(length_*sizeof(vari*));
669  for (size_t i = 0; i < length_; i++)
670  v1_[i] = v1[i].vi_;
671  } else {
672  v1_ = shared_v1->v1_;
673  }
674  if (shared_v2 == NULL) {
675  v2_ = (double*)memalloc_.alloc(length_*sizeof(double));
676  for (size_t i = 0; i < length_; i++)
677  v2_[i] = v2[i];
678  } else {
679  v2_ = shared_v2->v2_;
680  }
681  }
682  void chain() {
683  for (size_t i = 0; i < length_; i++) {
684  v1_[i]->adj_ += adj_ * v2_[i];
685  }
686  }
687  };
688 
689  }
690 
701  template<int R, int C>
702  inline var dot_self(const Eigen::Matrix<var, R, C>& v) {
703  stan::math::validate_vector(v,"dot_self");
704  return var(new dot_self_vari(v));
705  }
706 
715  template<int R1,int C1,int R2, int C2>
716  inline var dot_product(const Eigen::Matrix<var, R1, C1>& v1,
717  const Eigen::Matrix<var, R2, C2>& v2) {
718  stan::math::validate_vector(v1,"dot_product");
719  stan::math::validate_vector(v2,"dot_product");
720  stan::math::validate_matching_sizes(v1,v2,"dot_product");
721  return var(new dot_product_vv_vari(v1,v2));
722  }
732  template<int R1,int C1,int R2, int C2>
733  inline var dot_product(const Eigen::Matrix<var, R1, C1>& v1,
734  const Eigen::Matrix<double, R2, C2>& v2) {
735  stan::math::validate_vector(v1,"dot_product");
736  stan::math::validate_vector(v2,"dot_product");
737  stan::math::validate_matching_sizes(v1,v2,"dot_product");
738  return var(new dot_product_vd_vari(v1,v2));
739  }
749  template<int R1,int C1,int R2, int C2>
750  inline var dot_product(const Eigen::Matrix<double, R1, C1>& v1,
751  const Eigen::Matrix<var, R2, C2>& v2) {
752  stan::math::validate_vector(v1,"dot_product");
753  stan::math::validate_vector(v2,"dot_product");
754  stan::math::validate_matching_sizes(v1,v2,"dot_product");
755  return var(new dot_product_vd_vari(v2,v1));
756  }
765  inline var dot_product(const var* v1, const var* v2, size_t length) {
766  return var(new dot_product_vv_vari(v1, v2, length));
767  }
776  inline var dot_product(const var* v1, const double* v2, size_t length) {
777  return var(new dot_product_vd_vari(v1, v2, length));
778  }
787  inline var dot_product(const double* v1, const var* v2, size_t length) {
788  return var(new dot_product_vd_vari(v2, v1, length));
789  }
798  inline var dot_product(const std::vector<var>& v1,
799  const std::vector<var>& v2) {
800  stan::math::validate_matching_sizes(v1,v2,"dot_product");
801  return var(new dot_product_vv_vari(&v1[0], &v2[0], v1.size()));
802  }
811  inline var dot_product(const std::vector<var>& v1,
812  const std::vector<double>& v2) {
813  stan::math::validate_matching_sizes(v1,v2,"dot_product");
814  return var(new dot_product_vd_vari(&v1[0], &v2[0], v1.size()));
815  }
824  inline var dot_product(const std::vector<double>& v1,
825  const std::vector<var>& v2) {
826  stan::math::validate_matching_sizes(v1,v2,"dot_product");
827  return var(new dot_product_vd_vari(&v2[0], &v1[0], v1.size()));
828  }
829 
836  template <int R, int C>
837  inline var sum(const Eigen::Matrix<var,R,C>& m) {
838  if (m.size() == 0)
839  return 0.0;
840  return var(new sum_v_vari(m));
841  }
842 
843 
851  inline double
852  divide(double x, double y) {
853  return x / y;
854  }
855  template <typename T1, typename T2>
856  inline var
857  divide(const T1& v, const T2& c) {
858  return to_var(v) / to_var(c);
859  }
867  template <typename T1, typename T2, int R, int C>
868  inline Eigen::Matrix<var,R,C>
869  divide(const Eigen::Matrix<T1, R,C>& v, const T2& c) {
870  return to_var(v) / to_var(c);
871  }
872 
873 
874 
881  template <typename T1, typename T2>
882  inline
883  typename boost::math::tools::promote_args<T1,T2>::type
884  multiply(const T1& v, const T2& c) {
885  return v * c;
886  }
887 
894  template<typename T1,typename T2,int R2,int C2>
895  inline Eigen::Matrix<var,R2,C2> multiply(const T1& c,
896  const Eigen::Matrix<T2, R2, C2>& m) {
897  // FIXME: pull out to eliminate overpromotion of one side
898  // move to matrix.hpp w. promotion?
899  return to_var(m) * to_var(c);
900  }
901 
908  template<typename T1,int R1,int C1,typename T2>
909  inline Eigen::Matrix<var,R1,C1> multiply(const Eigen::Matrix<T1, R1, C1>& m,
910  const T2& c) {
911  return to_var(m) * to_var(c);
912  }
913 
924  template<int R1,int C1,int R2,int C2>
925  inline Eigen::Matrix<var,R1,C2> multiply(const Eigen::Matrix<var,R1,C1>& m1,
926  const Eigen::Matrix<var,R2,C2>& m2) {
927  stan::math::validate_multiplicable(m1,m2,"multiply");
928  Eigen::Matrix<var,R1,C2> result(m1.rows(),m2.cols());
929  for (int i = 0; i < m1.rows(); i++) {
930  typename Eigen::Matrix<var,R1,C1>::ConstRowXpr crow(m1.row(i));
931  for (int j = 0; j < m2.cols(); j++) {
932  typename Eigen::Matrix<var,R2,C2>::ConstColXpr ccol(m2.col(j));
933  if (j == 0) {
934  if (i == 0) {
935  result(i,j) = var(new dot_product_vv_vari(crow,ccol));
936  }
937  else {
938  dot_product_vv_vari *v2 = static_cast<dot_product_vv_vari*>(result(0,j).vi_);
939  result(i,j) = var(new dot_product_vv_vari(crow,ccol,NULL,v2));
940  }
941  }
942  else {
943  if (i == 0) {
944  dot_product_vv_vari *v1 = static_cast<dot_product_vv_vari*>(result(i,0).vi_);
945  result(i,j) = var(new dot_product_vv_vari(crow,ccol,v1));
946  }
947  else /* if (i != 0 && j != 0) */ {
948  dot_product_vv_vari *v1 = static_cast<dot_product_vv_vari*>(result(i,0).vi_);
949  dot_product_vv_vari *v2 = static_cast<dot_product_vv_vari*>(result(0,j).vi_);
950  result(i,j) = var(new dot_product_vv_vari(crow,ccol,v1,v2));
951  }
952  }
953  }
954  }
955  return result;
956  }
957 
968  template<int R1,int C1,int R2,int C2>
969  inline Eigen::Matrix<var,R1,C2> multiply(const Eigen::Matrix<double,R1,C1>& m1,
970  const Eigen::Matrix<var,R2,C2>& m2) {
971  stan::math::validate_multiplicable(m1,m2,"multiply");
972  Eigen::Matrix<var,R1,C2> result(m1.rows(),m2.cols());
973  for (int i = 0; i < m1.rows(); i++) {
974  typename Eigen::Matrix<double,R1,C1>::ConstRowXpr crow(m1.row(i));
975  for (int j = 0; j < m2.cols(); j++) {
976  typename Eigen::Matrix<var,R2,C2>::ConstColXpr ccol(m2.col(j));
977 // result(i,j) = dot_product(crow,ccol);
978  if (j == 0) {
979  if (i == 0) {
980  result(i,j) = var(new dot_product_vd_vari(ccol,crow));
981  }
982  else {
983  dot_product_vd_vari *v2 = static_cast<dot_product_vd_vari*>(result(0,j).vi_);
984  result(i,j) = var(new dot_product_vd_vari(ccol,crow,v2,NULL));
985  }
986  }
987  else {
988  if (i == 0) {
989  dot_product_vd_vari *v1 = static_cast<dot_product_vd_vari*>(result(i,0).vi_);
990  result(i,j) = var(new dot_product_vd_vari(ccol,crow,NULL,v1));
991  }
992  else /* if (i != 0 && j != 0) */ {
993  dot_product_vd_vari *v1 = static_cast<dot_product_vd_vari*>(result(i,0).vi_);
994  dot_product_vd_vari *v2 = static_cast<dot_product_vd_vari*>(result(0,j).vi_);
995  result(i,j) = var(new dot_product_vd_vari(ccol,crow,v2,v1));
996  }
997  }
998  }
999  }
1000  return result;
1001  }
1002 
1013  template<int R1,int C1,int R2,int C2>
1014  inline Eigen::Matrix<var,R1,C2> multiply(const Eigen::Matrix<var,R1,C1>& m1,
1015  const Eigen::Matrix<double,R2,C2>& m2) {
1016  stan::math::validate_multiplicable(m1,m2,"multiply");
1017  Eigen::Matrix<var,R1,C2> result(m1.rows(),m2.cols());
1018  for (int i = 0; i < m1.rows(); i++) {
1019  typename Eigen::Matrix<var,R1,C1>::ConstRowXpr crow(m1.row(i));
1020  for (int j = 0; j < m2.cols(); j++) {
1021  typename Eigen::Matrix<double,R2,C2>::ConstColXpr ccol(m2.col(j));
1022 // result(i,j) = dot_product(crow,ccol);
1023  if (j == 0) {
1024  if (i == 0) {
1025  result(i,j) = var(new dot_product_vd_vari(crow,ccol));
1026  }
1027  else {
1028  dot_product_vd_vari *v2 = static_cast<dot_product_vd_vari*>(result(0,j).vi_);
1029  result(i,j) = var(new dot_product_vd_vari(crow,ccol,NULL,v2));
1030  }
1031  }
1032  else {
1033  if (i == 0) {
1034  dot_product_vd_vari *v1 = static_cast<dot_product_vd_vari*>(result(i,0).vi_);
1035  result(i,j) = var(new dot_product_vd_vari(crow,ccol,v1,NULL));
1036  }
1037  else /* if (i != 0 && j != 0) */ {
1038  dot_product_vd_vari *v1 = static_cast<dot_product_vd_vari*>(result(i,0).vi_);
1039  dot_product_vd_vari *v2 = static_cast<dot_product_vd_vari*>(result(0,j).vi_);
1040  result(i,j) = var(new dot_product_vd_vari(crow,ccol,v1,v2));
1041  }
1042  }
1043  }
1044  }
1045  return result;
1046  }
1047 
1057  template <int C1,int R2>
1058  inline var multiply(const Eigen::Matrix<var, 1, C1>& rv,
1059  const Eigen::Matrix<var, R2, 1>& v) {
1060  if (rv.size() != v.size())
1061  throw std::domain_error("row vector and vector must be same length in multiply");
1062  return dot_product(rv, v);
1063  }
1073  template <int C1,int R2>
1074  inline var multiply(const Eigen::Matrix<double, 1, C1>& rv,
1075  const Eigen::Matrix<var, R2, 1>& v) {
1076  stan::math::validate_multiplicable(rv,v,"multiply");
1077  return dot_product(rv, v);
1078  }
1088  template <int C1,int R2>
1089  inline var multiply(const Eigen::Matrix<var, 1, C1>& rv,
1090  const Eigen::Matrix<double, R2, 1>& v) {
1091  stan::math::validate_multiplicable(rv,v,"multiply");
1092  return dot_product(rv, v);
1093  }
1094 
1095  inline matrix_v
1097 // stan::math::validate_square(L,"multiply_lower_tri_self_transpose");
1098  int K = L.rows();
1099  int J = L.cols();
1100  matrix_v LLt(K,K);
1101  if (K == 0) return LLt;
1102  // if (K == 1) {
1103  // LLt(0,0) = L(0,0) * L(0,0);
1104  // return LLt;
1105  // }
1106  int Knz;
1107  if (K >= J)
1108  Knz = (K-J)*J + (J * (J + 1)) / 2;
1109  else // if (K < J)
1110  Knz = (K * (K + 1)) / 2;
1111  vari** vs = (vari**)memalloc_.alloc( Knz * sizeof(vari*) );
1112  int pos = 0;
1113  for (int m = 0; m < K; ++m)
1114  for (int n = 0; n < ((J < (m+1))?J:(m+1)); ++n) {
1115  vs[pos++] = L(m,n).vi_;
1116  }
1117  for (int m = 0, mpos=0; m < K; ++m, mpos += (J < m)?J:m) {
1118  LLt(m,m) = var(new dot_self_vari(vs + mpos, (J < (m+1))?J:(m+1)));
1119  for (int n = 0, npos = 0; n < m; ++n, npos += (J < n)?J:n) {
1120  LLt(m,n) = LLt(n,m) = var(new dot_product_vv_vari(vs + mpos, vs + npos, (J < (n+1))?J:(n+1)));
1121  }
1122  }
1123  return LLt;
1124  }
1125 
1132  inline matrix_v
1133  tcrossprod(const matrix_v& M) {
1134  if (M.rows() == 0)
1135  return matrix_v(0,0);
1136  if (M.rows() == 1)
1137  return M * M.transpose();
1138 
1139  // WAS JUST THIS
1140  // matrix_v result(M.rows(),M.rows());
1141  // return result.setZero().selfadjointView<Eigen::Upper>().rankUpdate(M);
1142 
1143  matrix_v MMt(M.rows(),M.rows());
1144 
1145  vari** vs
1146  = (vari**)memalloc_.alloc((M.rows() * M.cols() ) * sizeof(vari*));
1147  int pos = 0;
1148  for (int m = 0; m < M.rows(); ++m)
1149  for (int n = 0; n < M.cols(); ++n)
1150  vs[pos++] = M(m,n).vi_;
1151  for (int m = 0; m < M.rows(); ++m)
1152  MMt(m,m) = var(new dot_self_vari(vs + m * M.cols(),M.cols()));
1153  for (int m = 0; m < M.rows(); ++m) {
1154  for (int n = 0; n < m; ++n) {
1155  MMt(m,n) = var(new dot_product_vv_vari(vs + m * M.cols(),
1156  vs + n * M.cols(),
1157  M.cols()));
1158  MMt(n,m) = MMt(m,n);
1159  }
1160  }
1161  return MMt;
1162  }
1163 
1170  inline matrix_v
1171  crossprod(const matrix_v& M) {
1172  return tcrossprod(M.transpose());
1173  }
1174 
1175 
1176  // FIXME: double val?
1177  inline void assign_to_var(stan::agrad::var& var, const double& val) {
1178  var = val;
1179  }
1181  var = val;
1182  }
1183  // FIXME: int val?
1184  inline void assign_to_var(int& n_lhs, const int& n_rhs) {
1185  n_lhs = n_rhs; // FIXME: no call -- just filler to instantiate
1186  }
1187  // FIXME: double val?
1188  inline void assign_to_var(double& n_lhs, const double& n_rhs) {
1189  n_lhs = n_rhs; // FIXME: no call -- just filler to instantiate
1190  }
1191 
1192  template <typename LHS, typename RHS>
1193  inline void assign_to_var(std::vector<LHS>& x, const std::vector<RHS>& y) {
1194  for (size_t i = 0; i < x.size(); ++i)
1195  assign_to_var(x[i],y[i]);
1196  }
1197  template <typename LHS, typename RHS>
1198  inline void assign_to_var(Eigen::Matrix<LHS,Eigen::Dynamic,1>& x,
1199  const Eigen::Matrix<RHS,Eigen::Dynamic,1>& y) {
1200  typedef
1202  size_type;
1203  for (size_type i = 0; i < x.size(); ++i)
1204  assign_to_var(x(i),y(i));
1205  }
1206  template <typename LHS, typename RHS>
1207  inline void assign_to_var(Eigen::Matrix<LHS,1,Eigen::Dynamic>& x,
1208  const Eigen::Matrix<RHS,1,Eigen::Dynamic>& y) {
1209  typedef
1211  size_type;
1212  for (size_type i = 0; i < x.size(); ++i)
1213  assign_to_var(x(i),y(i));
1214  }
1215  template <typename LHS, typename RHS>
1216  inline void assign_to_var(Eigen::Matrix<LHS,Eigen::Dynamic,Eigen::Dynamic>& x,
1217  const Eigen::Matrix<RHS,Eigen::Dynamic,Eigen::Dynamic>& y) {
1218  typedef
1220  size_type1;
1221  typedef
1223  size_type2;
1224  for (size_type1 n = 0; n < x.cols(); ++n)
1225  for (size_type2 m = 0; m < x.rows(); ++m)
1226  assign_to_var(x(m,n),y(m,n));
1227  }
1228 
1229 
1230  template <typename LHS, typename RHS>
1234  };
1235 
1236  template <bool PromoteRHS, typename LHS, typename RHS>
1237  struct assigner {
1238  static inline void assign(LHS& var, const RHS& val) {
1239  throw std::domain_error("should not call base class of assigner");
1240  }
1241  };
1242 
1243  template <typename LHS, typename RHS>
1244  struct assigner<false,LHS,RHS> {
1245  static inline void assign(LHS& var, const RHS& val) {
1246  var = val; // no promotion of RHS
1247  }
1248  };
1249 
1250  template <typename LHS, typename RHS>
1251  struct assigner<true,LHS,RHS> {
1252  static inline void assign(LHS& var, const RHS& val) {
1253  assign_to_var(var,val); // promote RHS
1254  }
1255  };
1256 
1257 
1258  template <typename LHS, typename RHS>
1259  inline void assign(LHS& var, const RHS& val) {
1261  }
1262 
1263  void stan_print(std::ostream* o, const var& x) {
1264  *o << x.val();
1265  }
1266 
1267  }
1268 }
1269 
1270 
1271 #endif
1272 
internal::traits< Derived >::Index size_type
vari ** v2_
Definition: matrix.hpp:494
size_t size_
Definition: matrix.hpp:413
vari ** v1_
Definition: matrix.hpp:493
vari ** v_
Definition: matrix.hpp:412
size_t length_
Definition: matrix.hpp:450
void chain()
Apply the chain rule to this variable based on the variables on which it depends.
Definition: matrix.hpp:41
stan::agrad::vari * alpha_
Definition: matrix.hpp:13
static double eval_gevv(const stan::agrad::var *alpha, const stan::agrad::var *v1, int stride1, const stan::agrad::var *v2, int stride2, size_t length, double *dotprod)
Definition: matrix.hpp:18
gevv_vvv_vari(const stan::agrad::var *alpha, const stan::agrad::var *v1, int stride1, const stan::agrad::var *v2, int stride2, size_t length)
Definition: matrix.hpp:29
stan::agrad::vari ** v1_
Definition: matrix.hpp:14
stan::agrad::vari ** v2_
Definition: matrix.hpp:15
Independent (input) and dependent (output) variables for gradients.
Definition: agrad.hpp:214
double val() const
Return the value of this variable.
Definition: agrad.hpp:396
vari * vi_
Pointer to the implementation of this variable.
Definition: agrad.hpp:227
The variable implementation base class.
Definition: agrad.hpp:104
const double val_
The value of this variable.
Definition: agrad.hpp:113
double adj_
The adjoint of this variable, which is the partial derivative of this variable with respect to the ro...
Definition: agrad.hpp:119
void * alloc(size_t len)
Return a newly allocated block of memory of the appropriate size managed by the stack allocator.
(Expert) Numerical traits for algorithmic differentiation variables.
Definition: matrix.hpp:53
void stan_print(std::ostream *o, const var &x)
Definition: matrix.hpp:1263
matrix_v multiply_lower_tri_self_transpose(const matrix_v &L)
Definition: matrix.hpp:1096
var square(const var &x)
Return the square of the input variable.
var to_var(const double &x)
Converts argument to an automatic differentiation variable.
Definition: matrix.hpp:298
Eigen::Matrix< var, Eigen::Dynamic, Eigen::Dynamic > matrix_v
The type of a matrix holding stan::agrad::var values.
Definition: matrix.hpp:244
var dot_product(const Eigen::Matrix< var, R1, C1 > &v1, const Eigen::Matrix< var, R2, C2 > &v2)
Returns the dot product.
Definition: matrix.hpp:716
var ceil(const var &a)
Return the ceiling of the specified variable (cmath).
Definition: agrad.hpp:2073
matrix_v crossprod(const matrix_v &M)
Returns the result of pre-multiplying a matrix by its own transpose.
Definition: matrix.hpp:1171
matrix_v tcrossprod(const matrix_v &M)
Returns the result of post-multiplying a matrix by its own transpose.
Definition: matrix.hpp:1133
void assign_to_var(stan::agrad::var &var, const double &val)
Definition: matrix.hpp:1177
var dot_self(const Eigen::Matrix< var, R, C > &v)
Returns the dot product of a vector with itself.
Definition: matrix.hpp:702
void assign(LHS &var, const RHS &val)
Definition: matrix.hpp:1259
memory::stack_alloc memalloc_
Definition: agrad.cpp:11
boost::math::tools::promote_args< T1, T2 >::type multiply(const T1 &v, const T2 &c)
Return the product of two scalars.
Definition: matrix.hpp:884
Eigen::Matrix< var, 1, Eigen::Dynamic > row_vector_v
The type of a row vector holding stan::agrad::var values.
Definition: matrix.hpp:260
Eigen::Matrix< var, Eigen::Dynamic, 1 > vector_v
The type of a (column) vector holding stan::agrad::var values.
Definition: matrix.hpp:252
var sum(const Eigen::Matrix< var, R, C > &m)
Returns the sum of the coefficients of the specified matrix, column vector or row vector.
Definition: matrix.hpp:837
void initialize_variable(var &variable, const var &value)
Initialize variable to value.
Definition: matrix.hpp:266
double divide(double x, double y)
Return the division of the first scalar by the second scalar.
Definition: matrix.hpp:852
Eigen::Matrix< double, Eigen::Dynamic, 1 > vector_d
Type for (column) vector of double values.
Definition: matrix.hpp:158
double epsilon()
Return minimum positive number representable.
void validate_multiplicable(const Eigen::Matrix< T1, R1, C1 > &x1, const Eigen::Matrix< T2, R2, C2 > &x2, const char *msg)
Definition: matrix.hpp:784
double e()
Return the base of the natural logarithm.
size_t cols(const Eigen::Matrix< T, R, C > &m)
Definition: matrix.hpp:628
void validate_vector(const Eigen::Matrix< T, R, C > &x, const char *msg)
Definition: matrix.hpp:806
Eigen::Matrix< double, 1, Eigen::Dynamic > row_vector_d
Type for (row) vector of double values.
Definition: matrix.hpp:165
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > matrix_d
Type for matrix of double values.
Definition: matrix.hpp:151
void validate_matching_sizes(const std::vector< T1 > &x1, const std::vector< T2 > &x2, const char *msg)
Definition: matrix.hpp:757
size_t rows(const Eigen::Matrix< T, R, C > &m)
Definition: matrix.hpp:622
Eigen::Matrix< T, Rows, Cols > log(const Eigen::Matrix< T, Rows, Cols > &m)
Return the element-wise logarithm of the matrix or vector.
Definition: matrix.hpp:1198
T max(const Eigen::Matrix< T, R, C > &m)
Returns the maximum coefficient in the specified vector, row vector, or matrix.
Definition: matrix.hpp:1003
Probability, optimization and sampling library.
Definition: agrad.cpp:6
size_t length(const T &x)
Definition: traits.hpp:111
stan::agrad::var NonInteger
Non-integer valued variables.
Definition: matrix.hpp:73
static Real lowest()
Return standard library's lowest for double-precision floating point, &#45;std::numeric_limits<double>max...
Definition: matrix.hpp:115
stan::agrad::var Real
Real-valued variables.
Definition: matrix.hpp:66
static Real highest()
Return standard library's highest for double-precision floating point, std::numeric_limits<double>max...
Definition: matrix.hpp:105
static Real epsilon()
Return standard library's epsilon for double-precision floating point, std::numeric_limits<double>eps...
Definition: matrix.hpp:88
stan::agrad::var Nested
Nested variables.
Definition: matrix.hpp:80
static Real dummy_precision()
Return dummy precision.
Definition: matrix.hpp:95
static void run(Index rows, Index cols, Index depth, const LhsScalar *_lhs, Index lhsStride, const RhsScalar *_rhs, Index rhsStride, ResScalar *res, Index resStride, const ResScalar &alpha, level3_blocking< LhsScalar, RhsScalar > &blocking, GemmParallelInfo< Index > *info=0)
Definition: matrix.hpp:215
static EIGEN_DONT_INLINE void run(Index rows, Index cols, const LhsScalar *lhs, Index lhsStride, const RhsScalar *rhs, Index rhsIncr, ResScalar *res, Index resIncr, const ResScalar &alpha)
Definition: matrix.hpp:179
static EIGEN_DONT_INLINE void run(Index rows, Index cols, const LhsScalar *lhs, Index lhsStride, const RhsScalar *rhs, Index rhsIncr, ResScalar *res, Index resIncr, const RhsScalar &alpha)
Definition: matrix.hpp:198
static void assign(LHS &var, const RHS &val)
Definition: matrix.hpp:1245
static void assign(LHS &var, const RHS &val)
Definition: matrix.hpp:1252
static void assign(LHS &var, const RHS &val)
Definition: matrix.hpp:1238
Metaprogram to determine if a type has a base scalar type that can be assigned to type double.
Definition: traits.hpp:39

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