Stan Math Library  2.20.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 
9 #include <boost/random/normal_distribution.hpp>
10 #include <boost/random/variate_generator.hpp>
11 
12 namespace stan {
13 namespace math {
14 
31 template <typename T_loc, class RNG>
33 multi_normal_rng(const T_loc& mu,
34  const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& S,
35  RNG& rng) {
36  using boost::normal_distribution;
37  using boost::variate_generator;
38 
39  static const char* function = "multi_normal_rng";
40 
41  check_positive(function, "Covariance matrix rows", S.rows());
42  check_symmetric(function, "Covariance matrix", S);
43  Eigen::LLT<Eigen::MatrixXd> llt_of_S = S.llt();
44  check_pos_definite("multi_normal_rng", "covariance matrix argument",
45  llt_of_S);
46 
47  vector_seq_view<T_loc> mu_vec(mu);
48  size_t size_mu = mu_vec[0].size();
49 
50  size_t N = length_mvt(mu);
51  int size_mu_old = size_mu;
52  for (size_t i = 1; i < N; i++) {
53  int size_mu_new = mu_vec[i].size();
54  check_size_match(function,
55  "Size of one of the vectors of "
56  "the location variable",
57  size_mu_new,
58  "Size of another vector of the "
59  "location variable",
60  size_mu_old);
61  size_mu_old = size_mu_new;
62  }
63 
64  for (size_t i = 0; i < N; i++) {
65  check_finite(function, "Location parameter", mu_vec[i]);
66  }
67 
69 
70  variate_generator<RNG&, normal_distribution<> > std_normal_rng(
71  rng, normal_distribution<>(0, 1));
72 
73  for (size_t n = 0; n < N; ++n) {
74  Eigen::VectorXd z(S.cols());
75  for (int i = 0; i < S.cols(); i++)
76  z(i) = std_normal_rng();
77 
78  output[n] = Eigen::VectorXd(mu_vec[n]) + llt_of_S.matrixL() * z;
79  }
80 
81  return output.data();
82 }
83 
84 } // namespace math
85 } // namespace stan
86 #endif
void check_finite(const char *function, const char *name, const T_y &y)
Check if y is finite.
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Check if the provided sizes match.
StdVectorBuilder allocates type T1 values to be used as intermediate values.
StdVectorBuilder< true, Eigen::VectorXd, T_loc >::type multi_normal_rng(const T_loc &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
Return a multivariate normal random variate with the given location and covariance using the specifie...
This class provides a low-cost wrapper for situations where you either need an Eigen Vector or RowVec...
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
void check_pos_definite(const char *function, const char *name, const Eigen::Matrix< T_y, -1, -1 > &y)
Check if the specified square, symmetric matrix is positive definite.
void check_symmetric(const char *function, const char *name, const matrix_cl &y)
Check if the matrix_cl is symmetric.
size_t length_mvt(const Eigen::Matrix< T, R, C > &)
Definition: length_mvt.hpp:12

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