1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTI_STUDENT_T_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTI_STUDENT_T_HPP
6 #include <boost/math/special_functions/gamma.hpp>
7 #include <boost/math/special_functions/fpclassify.hpp>
20 #include <boost/random/variate_generator.hpp>
32 template <
bool propto,
33 typename T_y,
typename T_dof,
typename T_loc,
typename T_scale>
34 typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,T_dof,
typename scalar_type<T_loc>::type,T_scale>::type
39 Eigen::Matrix<T_scale,
40 Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
41 static const char*
function =
"stan::prob::multi_student_t(%1%)";
48 using boost::math::tools::promote_args;
54 typedef typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,T_dof,
typename scalar_type<T_loc>::type,T_scale>::type lp_type;
59 "Degrees of freedom parameter", &lp);
61 "Degrees of freedom parameter", &lp);
78 int size_y = y_vec[0].size();
79 int size_mu = mu_vec[0].size();
81 int size_y_old = size_y;
84 int size_y_new = y_vec[i].size();
86 size_y_new,
"Size of one of the vectors of the random variable",
87 size_y_old,
"Size of another vector of the random variable",
89 size_y_old = size_y_new;
91 int size_mu_old = size_mu;
94 int size_mu_new = mu_vec[i].size();
96 size_mu_new,
"Size of one of the vectors of the location variable",
97 size_mu_old,
"Size of another vector of the location variable",
99 size_mu_old = size_mu_new;
109 size_y,
"Size of random variable",
110 size_mu,
"size of location parameter",
113 size_y,
"Size of random variable",
114 Sigma.rows(),
"rows of scale parameter",
117 size_y,
"Size of random variable",
118 Sigma.cols(),
"columns of scale parameter",
121 for (
size_t i = 0; i < size_vec; i++) {
122 check_finite(
function, mu_vec[i],
"Location parameter", &lp);
128 LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_Sigma(Sigma);
136 lp +=
lgamma(0.5 * (nu + size_y)) * size_vec;
137 lp -=
lgamma(0.5 * nu) * size_vec;
138 lp -= (0.5 * size_y) *
log(nu) * size_vec;
142 lp -= (0.5 * size_y) * LOG_PI * size_vec;
155 lp_type sum_lp_vec(0.0);
156 for (
size_t i = 0; i < size_vec; i++) {
158 boost::math::tools::promote_args<typename scalar_type<T_y>::type,
typename scalar_type<T_loc>::type>::type,
159 Dynamic, 1> y_minus_mu(size_y);
160 for (
int j = 0; j < size_y; j++)
161 y_minus_mu(j) = y_vec[i](j)-mu_vec[i](j);
164 lp -= 0.5 * (nu + size_y) * sum_lp_vec;
169 template <
typename T_y,
typename T_dof,
typename T_loc,
typename T_scale>
171 typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,T_dof,
typename scalar_type<T_loc>::type,T_scale>::type
176 Eigen::Matrix<T_scale,
177 Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
178 return multi_student_t_log<false>(y,nu,mu,Sigma);
183 inline Eigen::VectorXd
185 const Eigen::Matrix<double,Eigen::Dynamic,1>& mu,
186 const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& s,
189 static const char*
function =
"stan::prob::multi_student_t_rng(%1%)";
196 check_finite(
function, mu,
"Location parameter", (
double*)0);
199 "Degrees of freedom parameter", (
double*)0);
201 "Degrees of freedom parameter", (
double*)0);
203 Eigen::VectorXd z(s.cols());
size_t max_size_mvt(const T1 &x1, const T2 &x2)
boost::math::tools::promote_args< typename scalar_type< T_y >::type, T_dof, typename scalar_type< T_loc >::type, T_scale >::type multi_student_t_log(const T_y &y, const T_dof &nu, const T_loc &mu, const Eigen::Matrix< T_scale, Eigen::Dynamic, Eigen::Dynamic > &Sigma)
Return the log of the multivariate Student t distribution at the specified arguments.
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
fvar< T > lgamma(const fvar< T > &x)
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.
bool isinf(const stan::agrad::var &v)
Checks if the given number is infinite.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R, C > subtract(const Eigen::Matrix< T1, R, C > &m1, const Eigen::Matrix< T2, R, C > &m2)
Return the result of subtracting the second specified matrix from the first specified matrix...
fvar< T > sqrt(const fvar< T > &x)
double inv_gamma_rng(const double alpha, const double beta, RNG &rng)
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_student_t_rng(const double nu, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &s, RNG &rng)
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.
T log_determinant_ldlt(stan::math::LDLT_factor< T, R, C > &A)
boost::enable_if_c< boost::is_arithmetic< T >::value, Eigen::Matrix< double, R, C > >::type multiply(const Eigen::Matrix< double, R, C > &m, T c)
Return specified matrix multiplied by specified scalar.
fvar< T > log(const fvar< T > &x)
double dot_product(const Eigen::Matrix< double, R1, C1 > &v1, const Eigen::Matrix< double, R2, C2 > &v2)
Returns the dot product of the specified vectors.