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 
10 #include <boost/random/normal_distribution.hpp>
11 #include <boost/random/variate_generator.hpp>
12 
13 namespace stan {
14 namespace math {
15 
16 /**
17  * Return a multivariate normal random variate with the given location and
18  * covariance using the specified random number generator.
19  *
20  * mu can be either an Eigen::VectorXd, an Eigen::RowVectorXd, or a std::vector
21  * of either of those types.
22  *
23  * @tparam T_loc Type of location paramater
24  * @tparam RNG Type of pseudo-random number generator
25  * @param mu (Sequence of) location parameter(s)
26  * @param S Covariance matrix
27  * @param rng random number generator
28  * @throw std::domain_error if S is not positive definite, or
29  * std::invalid_argument if the length of (each) mu is not equal to the number
30  * of rows and columns in S
31  */
32 template <typename T_loc, class RNG>
34 multi_normal_rng(const T_loc& mu,
35  const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& S,
36  RNG& rng) {
37  using boost::normal_distribution;
38  using boost::variate_generator;
39 
40  static const char* function = "multi_normal_rng";
41 
42  check_positive(function, "Covariance matrix rows", S.rows());
43  check_symmetric(function, "Covariance matrix", S);
44  Eigen::LLT<Eigen::MatrixXd> llt_of_S = S.llt();
45  check_pos_definite("multi_normal_rng", "covariance matrix argument",
46  llt_of_S);
47 
48  vector_seq_view<T_loc> mu_vec(mu);
49  size_t size_mu = mu_vec[0].size();
50 
51  size_t N = length_mvt(mu);
52  int size_mu_old = size_mu;
53  for (size_t i = 1; i < N; i++) {
54  int size_mu_new = mu_vec[i].size();
55  check_size_match(function,
56  "Size of one of the vectors of "
57  "the location variable",
58  size_mu_new,
59  "Size of another vector of the "
60  "location variable",
61  size_mu_old);
62  size_mu_old = size_mu_new;
63  }
64 
65  for (size_t i = 0; i < N; i++) {
66  check_finite(function, "Location parameter", mu_vec[i]);
67  }
68 
70 
71  variate_generator<RNG&, normal_distribution<> > std_normal_rng(
72  rng, normal_distribution<>(0, 1));
73 
74  for (size_t n = 0; n < N; ++n) {
75  Eigen::VectorXd z(S.cols());
76  for (int i = 0; i < S.cols(); i++)
77  z(i) = std_normal_rng();
78 
79  output[n] = Eigen::VectorXd(mu_vec[n]) + llt_of_S.matrixL() * z;
80  }
81 
82  return output.data();
83 }
84 
85 } // namespace math
86 } // namespace stan
87 #endif
ofstream output
void check_finite(const char *function, const char *name, const T_y &y)
#define S(x, n)
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
void check_symmetric(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
void check_pos_definite(const char *function, const char *name, const Eigen::Matrix< T_y,-1,-1 > &y)
z
Definition: test.py:28
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)
void check_positive(const char *function, const char *name, const T_y &y)
size_t length_mvt(const Eigen::Matrix< T, R, C > &)
Definition: length_mvt.hpp:12