multi_normal_cholesky_rng.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_CHOLESKY_RNG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_CHOLESKY_RNG_HPP
3 
8 #include <boost/random/normal_distribution.hpp>
9 #include <boost/random/variate_generator.hpp>
10 
11 namespace stan {
12 namespace math {
13 
14 /**
15  * Return a multivariate normal random variate with the given location and
16  * Cholesky factorization of the covariance using the specified random number
17  * generator.
18  *
19  * mu can be either an Eigen::VectorXd, an Eigen::RowVectorXd, or a std::vector
20  * of either of those types.
21  *
22  * @tparam T_loc Type of location paramater
23  * @tparam RNG Type of pseudo-random number generator
24  * @param mu (Sequence of) location parameter(s)
25  * @param L Lower Cholesky factor of target covariance matrix
26  * @param rng random number generator
27  * @throw std::invalid_argument if the length of (each) mu is not equal to the
28  * number of rows and columns in L
29  */
30 template <typename T_loc, class RNG>
33  const T_loc& mu,
34  const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& L, RNG& rng) {
35  using boost::normal_distribution;
36  using boost::variate_generator;
37 
38  static const char* function = "multi_normal_cholesky_rng";
39  vector_seq_view<T_loc> mu_vec(mu);
40  size_t size_mu = mu_vec[0].size();
41 
42  size_t N = length_mvt(mu);
43  int size_mu_old = size_mu;
44  for (size_t i = 1; i < N; i++) {
45  int size_mu_new = mu_vec[i].size();
46  check_size_match(function,
47  "Size of one of the vectors of "
48  "the location variable",
49  size_mu_new,
50  "Size of another vector of the "
51  "location variable",
52  size_mu_old);
53  size_mu_old = size_mu_new;
54  }
55 
56  for (size_t i = 0; i < N; i++) {
57  check_finite(function, "Location parameter", mu_vec[i]);
58  }
59 
61 
62  variate_generator<RNG&, normal_distribution<> > std_normal_rng(
63  rng, normal_distribution<>(0, 1));
64 
65  for (size_t n = 0; n < N; ++n) {
66  Eigen::VectorXd z(L.cols());
67  for (int i = 0; i < L.cols(); i++)
68  z(i) = std_normal_rng();
69 
70  output[n] = Eigen::VectorXd(mu_vec[n]) + L * z;
71  }
72 
73  return output.data();
74 }
75 } // namespace math
76 } // namespace stan
77 #endif
ofstream output
void check_finite(const char *function, const char *name, const T_y &y)
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
StdVectorBuilder< true, Eigen::VectorXd, T_loc >::type multi_normal_cholesky_rng(const T_loc &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &L, RNG &rng)
static constexpr double L
z
Definition: test.py:28
size_t length_mvt(const Eigen::Matrix< T, R, C > &)
Definition: length_mvt.hpp:12