Stan Math Library  2.9.0
reverse mode automatic differentiation
multi_normal_rng.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_RNG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_RNG_HPP
3 
4 #include <boost/random/normal_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
14 
17 
18 namespace stan {
19 
20  namespace math {
21 
22  template <class RNG>
23  inline Eigen::VectorXd
25  const Eigen::Matrix<double, Eigen::Dynamic, 1>& mu,
26  const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& S,
27  RNG& rng
28  ) {
29  using boost::variate_generator;
30  using boost::normal_distribution;
31 
32  static const char* function("stan::math::multi_normal_rng");
33 
37 
38  check_positive(function, "Covariance matrix rows", S.rows());
39  check_symmetric(function, "Covariance matrix", S);
40  check_finite(function, "Location parameter", mu);
41 
42  variate_generator<RNG&, normal_distribution<> >
43  std_normal_rng(rng, normal_distribution<>(0, 1));
44 
45  Eigen::VectorXd z(S.cols());
46  for (int i = 0; i < S.cols(); i++)
47  z(i) = std_normal_rng();
48 
49  return mu + S.llt().matrixL() * z;
50  }
51  }
52 }
53 
54 #endif
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
bool check_finite(const char *function, const char *name, const T_y &y)
Return true if y is finite.
Eigen::VectorXd multi_normal_rng(const Eigen::Matrix< double, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
bool check_symmetric(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is symmetric.

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