1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__MULTI_NORMAL_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__MULTI_NORMAL_HPP
4 #include <boost/random/normal_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
24 template <
bool propto,
25 typename T_y,
typename T_loc,
typename T_covar>
26 typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,
typename scalar_type<T_loc>::type, T_covar>::type
29 const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
30 static const char*
function =
"stan::prob::multi_normal_log(%1%)";
31 typedef typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,
typename scalar_type<T_loc>::type, T_covar>::type lp_type;
42 Sigma.rows(),
"Rows of covariance parameter",
43 Sigma.cols(),
"columns of covariance parameter",
45 check_positive(
function, Sigma.rows(),
"Covariance matrix rows", &lp);
57 int size_y = y_vec[0].size();
58 int size_mu = mu_vec[0].size();
60 int size_y_old = size_y;
63 int size_y_new = y_vec[i].size();
65 size_y_new,
"Size of one of the vectors of the random variable",
66 size_y_old,
"Size of another vector of the random variable",
68 size_y_old = size_y_new;
70 int size_mu_old = size_mu;
73 int size_mu_new = mu_vec[i].size();
75 size_mu_new,
"Size of one of the vectors of the location variable",
76 size_mu_old,
"Size of another vector of the location variable",
78 size_mu_old = size_mu_new;
87 size_y,
"Size of random variable",
88 size_mu,
"size of location parameter",
91 size_y,
"Size of random variable",
92 Sigma.rows(),
"rows of covariance parameter",
95 size_y,
"Size of random variable",
96 Sigma.cols(),
"columns of covariance parameter",
99 for (
size_t i = 0; i < size_vec; i++) {
100 check_finite(
function, mu_vec[i],
"Location parameter", &lp);
108 lp += NEG_LOG_SQRT_TWO_PI * size_y * size_vec;
114 lp_type sum_lp_vec(0.0);
115 for (
size_t i = 0; i < size_vec; i++) {
116 Eigen::Matrix<
typename
117 boost::math::tools::promote_args<typename scalar_type<T_y>::type,
typename scalar_type<T_loc>::type>::type,
118 Eigen::Dynamic, 1> y_minus_mu(size_y);
119 for (
int j = 0; j < size_y; j++)
120 y_minus_mu(j) = y_vec[i](j)-mu_vec[i](j);
123 lp -= 0.5*sum_lp_vec;
128 template <
typename T_y,
typename T_loc,
typename T_covar>
130 typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,
typename scalar_type<T_loc>::type, T_covar>::type
133 const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
134 return multi_normal_log<false>(y,mu,Sigma);
138 inline Eigen::VectorXd
140 const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
142 using boost::variate_generator;
143 using boost::normal_distribution;
145 static const char*
function =
"stan::prob::multi_normal_rng(%1%)";
151 check_positive(
function, S.rows(),
"Covariance matrix rows", (
double*)0);
153 check_finite(
function, mu,
"Location parameter", (
double*)0);
155 variate_generator<RNG&, normal_distribution<> >
156 std_normal_rng(rng, normal_distribution<>(0,1));
158 Eigen::VectorXd z(S.cols());
159 for(
int i = 0; i < S.cols(); i++)
160 z(i) = std_normal_rng();
162 return mu + S.llt().matrixL() * z;
size_t max_size_mvt(const T1 &x1, const T2 &x2)
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
scalar_type_helper< is_vector< T >::value, T >::type type
bool check_symmetric(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is symmetric.
bool check_ldlt_factor(const char *function, stan::math::LDLT_factor< T, R, C > &A, const char *name, T_result *result)
Return true if the underlying matrix is positive definite.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
size_t length_mvt(const T &)
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
Eigen::VectorXd multi_normal_rng(const Eigen::Matrix< double, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
bool check_not_nan(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is nan.
var log_determinant_ldlt(stan::math::LDLT_factor< var, R, C > &A)
boost::math::tools::promote_args< typename scalar_type< T_y >::type, typename scalar_type< T_loc >::type, T_covar >::type multi_normal_log(const T_y &y, const T_loc &mu, const Eigen::Matrix< T_covar, Eigen::Dynamic, Eigen::Dynamic > &Sigma)
boost::enable_if_c< boost::is_same< T2, var >::value||boost::is_same< T3, var >::value, var >::type trace_inv_quad_form_ldlt(const stan::math::LDLT_factor< T2, R2, C2 > &A, const Eigen::Matrix< T3, R3, C3 > &B)
Compute the trace of an inverse quadratic form.