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