1 #ifndef __STAN__AGRAD__MATRIX_HPP__
2 #define __STAN__AGRAD__MATRIX_HPP__
21 size_t length,
double *dotprod) {
23 for (
size_t i = 0; i <
length; i++)
24 result += v1[i*stride1].vi_->
val_ * v2[i*stride2].
vi_->
val_;
26 return alpha->vi_->
val_ * result;
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_;
43 for (
size_t i = 0; i <
length_; i++) {
59 template <>
struct NumTraits<
stan::agrad::var>
127 RequireInitialization = 0,
131 HasFloatingPoint = 1,
140 struct significant_decimals_default_impl<
stan::agrad::var,false>
155 struct scalar_product_traits<
stan::agrad::var,double> {
164 struct scalar_product_traits<double,
stan::agrad::var> {
171 template<
typename Index,
bool ConjugateLhs,
bool ConjugateRhs>
172 struct general_matrix_vector_product<Index,
stan::agrad::var,ColMajor,ConjugateLhs,
stan::agrad::var,ConjugateRhs>
176 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType
ResScalar;
177 enum { LhsStorageOrder = ColMajor };
179 EIGEN_DONT_INLINE
static void run(
185 for (Index i = 0; i <
rows; i++) {
190 template<
typename Index,
bool ConjugateLhs,
bool ConjugateRhs>
191 struct general_matrix_vector_product<Index,
stan::agrad::var,RowMajor,ConjugateLhs,
stan::agrad::var,ConjugateRhs>
195 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType
ResScalar;
196 enum { LhsStorageOrder = RowMajor };
198 EIGEN_DONT_INLINE
static void run(
204 for (Index i = 0; i <
rows; i++) {
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>
214 typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType
ResScalar;
220 level3_blocking<LhsScalar,RhsScalar>& blocking,
221 GemmParallelInfo<Index>* info = 0)
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);
243 Eigen::Matrix<var,Eigen::Dynamic,Eigen::Dynamic>
251 Eigen::Matrix<var,Eigen::Dynamic,1>
259 Eigen::Matrix<var,1,Eigen::Dynamic>
274 template <
int R,
int C>
276 for (
int i = 0; i < matrix.size(); ++i)
283 template <
typename T>
285 for (
size_t i = 0; i < variables.size(); ++i)
322 for (
int j = 0; j < m.cols(); ++j)
323 for (
int i = 0; i < m.rows(); ++i)
349 for (
int i = 0; i < v.size(); ++i)
376 for (
int i = 0; i < rv.size(); ++i)
410 class dot_self_vari :
public vari {
415 dot_self_vari(vari** v,
size_t size)
416 : vari(var_dot_self(v,size)),
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()) {
424 for (
size_t i = 0; i <
size_; ++i)
427 inline static double square(
double x) {
return x * x; }
428 inline static double var_dot_self(vari** v,
size_t size) {
430 for (
size_t i = 0; i < size; ++i)
434 template <
int R,
int C>
435 inline static double var_dot_self(
const Eigen::Matrix<var,R,C> &v) {
437 for (
int i = 0; i < v.size(); ++i)
442 for (
size_t i = 0; i <
size_; ++i)
443 v_[i]->adj_ += adj_ * 2.0 *
v_[i]->val_;
447 class sum_v_vari :
public vari{
451 inline static double var_sum(
const var *v,
size_t length) {
453 for (
size_t i = 0; i <
length; i++)
454 result += v[i].vi_->val_;
457 template<
typename Derived>
458 inline static double var_sum(
const Eigen::DenseBase<Derived> &v) {
460 for (
int i = 0; i < v.size(); i++)
461 result += v(i).vi_->val_;
465 template<
typename Derived>
466 sum_v_vari(
const Eigen::DenseBase<Derived> &v) :
467 vari(var_sum(v)),
length_(v.size()) {
469 for (
size_t i = 0; i <
length_; i++)
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()) {
476 for (
size_t i = 0; i <
length_; i++)
479 sum_v_vari(
const var *v,
size_t len) :
480 vari(var_sum(v,len)),
length_(len) {
482 for (
size_t i = 0; i <
length_; i++)
486 for (
size_t i = 0; i <
length_; i++) {
491 class dot_product_vv_vari :
public vari {
496 inline static double var_dot(
const var* v1,
const var* v2,
499 for (
size_t i = 0; i <
length; i++)
500 result += v1[i].vi_->val_ * v2[i].vi_->val_;
503 template<
typename Derived1,
typename Derived2>
504 inline static double var_dot(
const Eigen::DenseBase<Derived1> &v1,
505 const Eigen::DenseBase<Derived2> &v2) {
507 for (
int i = 0; i < v1.size(); i++)
508 result += v1[i].vi_->val_ * v2[i].vi_->val_;
511 inline static double var_dot(vari** v1, vari** v2,
size_t length) {
513 for (
size_t i = 0; i <
length; ++i)
514 result += v1[i]->val_ * v2[i]->val_;
518 dot_product_vv_vari(vari** v1, vari** v2,
size_t length)
519 : vari(var_dot(v1,v2,
length)),
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) :
529 if (shared_v1 == NULL) {
531 for (
size_t i = 0; i <
length_; i++)
535 v1_ = shared_v1->v1_;
537 if (shared_v2 == NULL) {
539 for (
size_t i = 0; i <
length_; i++)
543 v2_ = shared_v2->v2_;
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) {
554 for (
size_t i = 0; i <
length_; i++)
558 v1_ = shared_v1->v1_;
560 if (shared_v2 == NULL) {
562 for (
size_t i = 0; i <
length_; i++)
566 v2_ = shared_v2->v2_;
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) {
577 for (
size_t i = 0; i <
length_; i++)
581 v1_ = shared_v1->v1_;
583 if (shared_v2 == NULL) {
585 for (
size_t i = 0; i <
length_; i++)
589 v2_ = shared_v2->v2_;
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_;
600 class dot_product_vd_vari :
public vari {
605 inline static double var_dot(
const var* v1,
const double* v2,
608 for (
size_t i = 0; i <
length; i++)
609 result += v1[i].vi_->val_ * v2[i];
612 template<
typename Derived1,
typename Derived2>
613 inline static double var_dot(
const Eigen::DenseBase<Derived1> &v1,
614 const Eigen::DenseBase<Derived2> &v2) {
616 for (
int i = 0; i < v1.size(); i++)
617 result += v1[i].vi_->val_ * v2[i];
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) :
625 if (shared_v1 == NULL) {
627 for (
size_t i = 0; i <
length_; i++)
630 v1_ = shared_v1->v1_;
632 if (shared_v2 == NULL) {
634 for (
size_t i = 0; i <
length_; i++)
637 v2_ = shared_v2->v2_;
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) {
648 for (
size_t i = 0; i <
length_; i++)
651 v1_ = shared_v1->v1_;
653 if (shared_v2 == NULL) {
655 for (
size_t i = 0; i <
length_; i++)
658 v2_ = shared_v2->v2_;
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) {
669 for (
size_t i = 0; i <
length_; i++)
672 v1_ = shared_v1->v1_;
674 if (shared_v2 == NULL) {
676 for (
size_t i = 0; i <
length_; i++)
679 v2_ = shared_v2->v2_;
683 for (
size_t i = 0; i <
length_; i++) {
684 v1_[i]->adj_ += adj_ *
v2_[i];
701 template<
int R,
int C>
704 return var(
new dot_self_vari(v));
715 template<
int R1,
int C1,
int R2,
int C2>
717 const Eigen::Matrix<var, R2, C2>& v2) {
721 return var(
new dot_product_vv_vari(v1,v2));
732 template<
int R1,
int C1,
int R2,
int C2>
734 const Eigen::Matrix<double, R2, C2>& v2) {
738 return var(
new dot_product_vd_vari(v1,v2));
749 template<
int R1,
int C1,
int R2,
int C2>
751 const Eigen::Matrix<var, R2, C2>& v2) {
755 return var(
new dot_product_vd_vari(v2,v1));
766 return var(
new dot_product_vv_vari(v1, v2,
length));
777 return var(
new dot_product_vd_vari(v1, v2,
length));
788 return var(
new dot_product_vd_vari(v2, v1,
length));
799 const std::vector<var>& v2) {
801 return var(
new dot_product_vv_vari(&v1[0], &v2[0], v1.size()));
812 const std::vector<double>& v2) {
814 return var(
new dot_product_vd_vari(&v1[0], &v2[0], v1.size()));
825 const std::vector<var>& v2) {
827 return var(
new dot_product_vd_vari(&v2[0], &v1[0], v1.size()));
836 template <
int R,
int C>
837 inline var sum(
const Eigen::Matrix<var,R,C>& m) {
840 return var(
new sum_v_vari(m));
855 template <
typename T1,
typename T2>
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) {
881 template <
typename T1,
typename T2>
883 typename boost::math::tools::promote_args<T1,T2>::type
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) {
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,
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) {
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));
935 result(i,j) =
var(
new dot_product_vv_vari(crow,ccol));
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));
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));
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));
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) {
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));
980 result(i,j) =
var(
new dot_product_vd_vari(ccol,crow));
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));
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));
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));
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) {
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));
1025 result(i,j) =
var(
new dot_product_vd_vari(crow,ccol));
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));
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));
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));
1057 template <
int C1,
int R2>
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");
1073 template <
int C1,
int R2>
1075 const Eigen::Matrix<var, R2, 1>& v) {
1088 template <
int C1,
int R2>
1090 const Eigen::Matrix<double, R2, 1>& v) {
1101 if (K == 0)
return LLt;
1108 Knz = (K-J)*J + (J * (J + 1)) / 2;
1110 Knz = (K * (K + 1)) / 2;
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_;
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)));
1137 return M * M.transpose();
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(),
1158 MMt(n,m) = MMt(m,n);
1192 template <
typename LHS,
typename RHS>
1194 for (
size_t i = 0; i < x.size(); ++i)
1197 template <
typename LHS,
typename RHS>
1199 const Eigen::Matrix<RHS,Eigen::Dynamic,1>& y) {
1203 for (
size_type i = 0; i < x.size(); ++i)
1206 template <
typename LHS,
typename RHS>
1208 const Eigen::Matrix<RHS,1,Eigen::Dynamic>& y) {
1212 for (
size_type i = 0; i < x.size(); ++i)
1215 template <
typename LHS,
typename RHS>
1217 const Eigen::Matrix<RHS,Eigen::Dynamic,Eigen::Dynamic>& y) {
1224 for (size_type1 n = 0; n < x.cols(); ++n)
1225 for (size_type2 m = 0; m < x.rows(); ++m)
1230 template <
typename LHS,
typename RHS>
1236 template <
bool PromoteRHS,
typename LHS,
typename RHS>
1239 throw std::domain_error(
"should not call base class of assigner");
1243 template <
typename LHS,
typename RHS>
1250 template <
typename LHS,
typename RHS>
1258 template <
typename LHS,
typename RHS>
internal::traits< Derived >::Index size_type
void chain()
Apply the chain rule to this variable based on the variables on which it depends.
stan::agrad::vari * alpha_
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)
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)
Independent (input) and dependent (output) variables for gradients.
double val() const
Return the value of this variable.
vari * vi_
Pointer to the implementation of this variable.
The variable implementation base class.
const double val_
The value of this variable.
double adj_
The adjoint of this variable, which is the partial derivative of this variable with respect to the ro...
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.
void stan_print(std::ostream *o, const var &x)
matrix_v multiply_lower_tri_self_transpose(const matrix_v &L)
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.
Eigen::Matrix< var, Eigen::Dynamic, Eigen::Dynamic > matrix_v
The type of a matrix holding stan::agrad::var values.
var dot_product(const Eigen::Matrix< var, R1, C1 > &v1, const Eigen::Matrix< var, R2, C2 > &v2)
Returns the dot product.
var ceil(const var &a)
Return the ceiling of the specified variable (cmath).
matrix_v crossprod(const matrix_v &M)
Returns the result of pre-multiplying a matrix by its own transpose.
matrix_v tcrossprod(const matrix_v &M)
Returns the result of post-multiplying a matrix by its own transpose.
void assign_to_var(stan::agrad::var &var, const double &val)
var dot_self(const Eigen::Matrix< var, R, C > &v)
Returns the dot product of a vector with itself.
void assign(LHS &var, const RHS &val)
memory::stack_alloc memalloc_
boost::math::tools::promote_args< T1, T2 >::type multiply(const T1 &v, const T2 &c)
Return the product of two scalars.
Eigen::Matrix< var, 1, Eigen::Dynamic > row_vector_v
The type of a row vector holding stan::agrad::var values.
Eigen::Matrix< var, Eigen::Dynamic, 1 > vector_v
The type of a (column) vector holding stan::agrad::var values.
var sum(const Eigen::Matrix< var, R, C > &m)
Returns the sum of the coefficients of the specified matrix, column vector or row vector.
void initialize_variable(var &variable, const var &value)
Initialize variable to value.
double divide(double x, double y)
Return the division of the first scalar by the second scalar.
Eigen::Matrix< double, Eigen::Dynamic, 1 > vector_d
Type for (column) vector of double values.
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)
double e()
Return the base of the natural logarithm.
size_t cols(const Eigen::Matrix< T, R, C > &m)
void validate_vector(const Eigen::Matrix< T, R, C > &x, const char *msg)
Eigen::Matrix< double, 1, Eigen::Dynamic > row_vector_d
Type for (row) vector of double values.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > matrix_d
Type for matrix of double values.
void validate_matching_sizes(const std::vector< T1 > &x1, const std::vector< T2 > &x2, const char *msg)
size_t rows(const Eigen::Matrix< T, R, C > &m)
Eigen::Matrix< T, Rows, Cols > log(const Eigen::Matrix< T, Rows, Cols > &m)
Return the element-wise logarithm of the matrix or vector.
T max(const Eigen::Matrix< T, R, C > &m)
Returns the maximum coefficient in the specified vector, row vector, or matrix.
Probability, optimization and sampling library.
size_t length(const T &x)
stan::agrad::var NonInteger
Non-integer valued variables.
static Real lowest()
Return standard library's lowest for double-precision floating point, -std::numeric_limits<double>max...
stan::agrad::var Real
Real-valued variables.
static Real highest()
Return standard library's highest for double-precision floating point, std::numeric_limits<double>max...
static Real epsilon()
Return standard library's epsilon for double-precision floating point, std::numeric_limits<double>eps...
stan::agrad::var Nested
Nested variables.
static Real dummy_precision()
Return dummy precision.
stan::agrad::var LhsScalar
scalar_product_traits< LhsScalar, RhsScalar >::ReturnType ResScalar
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)
stan::agrad::var RhsScalar
stan::agrad::var RhsScalar
stan::agrad::var LhsScalar
scalar_product_traits< LhsScalar, RhsScalar >::ReturnType ResScalar
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)
scalar_product_traits< LhsScalar, RhsScalar >::ReturnType ResScalar
stan::agrad::var LhsScalar
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)
stan::agrad::var RhsScalar
stan::agrad::var ReturnType
stan::agrad::var ReturnType
static void assign(LHS &var, const RHS &val)
static void assign(LHS &var, const RHS &val)
static void assign(LHS &var, const RHS &val)
Metaprogram to determine if a type has a base scalar type that can be assigned to type double.