1 #ifndef STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_PREC_RNG_HPP 2 #define STAN_MATH_PRIM_MAT_PROB_MULTI_NORMAL_PREC_RNG_HPP 10 #include <boost/random/normal_distribution.hpp> 11 #include <boost/random/variate_generator.hpp> 32 template <
typename T_loc,
class RNG>
35 using boost::normal_distribution;
36 using boost::variate_generator;
38 static const char *
function =
"multi_normal_prec_rng";
44 Eigen::LLT<Eigen::MatrixXd> llt_of_S = S.llt();
50 size_t size_mu = mu_vec[0].size();
52 size_t N = mu_vec.size();
54 for (
size_t i = 1; i < N; i++) {
55 int size_mu_new = mu_vec[i].size();
57 "Size of one of the vectors of " 58 "the location variable",
60 "Size of another vector of the " 65 for (
size_t i = 0; i < N; i++) {
69 check_size_match(
function,
"Rows of location parameter", size_mu,
"Rows of S",
74 variate_generator<RNG &, normal_distribution<>> std_normal_rng(
75 rng, normal_distribution<>(0, 1));
77 for (
size_t n = 0; n < N; ++n) {
78 Eigen::VectorXd z(S.cols());
79 for (
int i = 0; i < S.cols(); i++)
80 z(i) = std_normal_rng();
82 output[n] = Eigen::VectorXd(mu_vec[n]) + llt_of_S.matrixU().solve(z);
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_prec_rng(const T_loc &mu, const Eigen::MatrixXd &S, RNG &rng)
Return a multivariate normal random variate with the given location and precision using the specified...
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.