Automatic Differentiation
 
Loading...
Searching...
No Matches
multi_student_t_rng.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_PROB_MULTI_STUDENT_T_RNG_HPP
2#define STAN_MATH_PRIM_PROB_MULTI_STUDENT_T_RNG_HPP
3
11#include <boost/random/normal_distribution.hpp>
12#include <boost/random/variate_generator.hpp>
13#include <cmath>
14
15namespace stan {
16namespace math {
17
40template <typename T_loc, class RNG>
43 double nu, const T_loc& mu,
44 const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& S, RNG& rng) {
45 using boost::normal_distribution;
46 using boost::variate_generator;
47 using boost::random::gamma_distribution;
48
49 static constexpr const char* function = "multi_student_t_rng";
50 check_not_nan(function, "Degrees of freedom parameter", nu);
51 check_positive(function, "Degrees of freedom parameter", nu);
52 check_positive(function, "Covariance matrix rows", S.rows());
53 vector_seq_view<T_loc> mu_vec(mu);
54 size_t size_mu = mu_vec[0].size();
55
56 size_t N = size_mvt(mu);
57 for (size_t i = 1; i < N; i++) {
58 check_size_match(function,
59 "Size of one of the vectors of "
60 "the location variable",
61 mu_vec[i].size(),
62 "Size of the first vector of the "
63 "location variable",
64 size_mu);
65 }
66
67 check_size_match(function, "Size of random variable", size_mu,
68 "rows of scale parameter", S.rows());
69
70 for (size_t i = 0; i < N; i++) {
71 check_finite(function, "Location parameter", mu_vec[i]);
72 }
73 const auto& S_ref = to_ref(S);
74 check_not_nan(function, "Covariance matrix", S_ref);
75 check_symmetric(function, "Covariance matrix", S_ref);
76 Eigen::LLT<Eigen::MatrixXd> llt_of_S = S_ref.llt();
77 check_pos_definite(function, "covariance matrix argument", llt_of_S);
78
80
81 variate_generator<RNG&, normal_distribution<> > std_normal_rng(
82 rng, normal_distribution<>(0, 1));
83
84 double w = inv_gamma_rng(nu / 2, nu / 2, rng);
85 for (size_t n = 0; n < N; ++n) {
86 Eigen::VectorXd z(S.cols());
87 for (int i = 0; i < S.cols(); i++) {
88 z(i) = std_normal_rng();
89 }
90 z *= std::sqrt(w);
91 output[n] = as_column_vector_or_scalar(mu_vec[n]) + llt_of_S.matrixL() * z;
92 }
93
94 return output.data();
95}
96
97} // namespace math
98} // namespace stan
99#endif
typename helper::type type
StdVectorBuilder allocates type T1 values to be used as intermediate values.
This class provides a low-cost wrapper for situations where you either need an Eigen Vector or RowVec...
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
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)
Return a multivariate student-t random variate with the given degrees of freedom location and Cholesk...
auto as_column_vector_or_scalar(T &&a)
as_column_vector_or_scalar of a kernel generator expression.
VectorBuilder< true, double, T_shape, T_scale >::type inv_gamma_rng(const T_shape &alpha, const T_scale &beta, RNG &rng)
Return a pseudorandom inverse gamma variate for the given shape and scale parameters using the specif...
double std_normal_rng(RNG &rng)
Return a standard Normal random variate using the specified random number generator.
int64_t size_mvt(const ScalarT &)
Provides the size of a multivariate argument.
Definition size_mvt.hpp:25
int64_t size(const T &m)
Returns the size (number of the elements) of a matrix_cl or var_value<matrix_cl<T>>.
Definition size.hpp:19
void check_pos_definite(const char *function, const char *name, const EigMat &y)
Check if the specified square, symmetric matrix is positive definite.
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:17
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
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.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...