1 #ifndef STAN_MATH_PRIM_MAT_PROB_MULTI_STUDENT_T_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_MULTI_STUDENT_T_LOG_HPP
19 #include <boost/math/special_functions/gamma.hpp>
20 #include <boost/random/variate_generator.hpp>
33 template <
bool propto,
34 typename T_y,
typename T_dof,
typename T_loc,
typename T_scale>
35 typename return_type<T_y, T_dof, T_loc, T_scale>::type
39 const T_scale& Sigma) {
40 static const char*
function(
"multi_student_t");
61 int size_y = y_vec[0].size();
62 int size_mu = mu_vec[0].size();
64 int size_y_old = size_y;
67 int size_y_new = y_vec[i].size();
69 "Size of one of the vectors of the random variable",
71 "Size of another vector of the random variable",
73 size_y_old = size_y_new;
75 int size_mu_old = size_mu;
78 int size_mu_new = mu_vec[i].size();
80 "Size of one of the vectors "
81 "of the location variable",
83 "Size of another vector of "
84 "the location variable",
86 size_mu_old = size_mu_new;
95 "Size of random variable", size_y,
96 "size of location parameter", size_mu);
98 "Size of random variable", size_y,
99 "rows of scale parameter", Sigma.rows());
101 "Size of random variable", size_y,
102 "columns of scale parameter", Sigma.cols());
104 for (
size_t i = 0; i < size_vec; i++) {
105 check_finite(
function,
"Location parameter", mu_vec[i]);
111 Eigen::Dynamic, Eigen::Dynamic> ldlt_Sigma(Sigma);
118 lp +=
lgamma(0.5 * (nu + size_y)) * size_vec;
119 lp -=
lgamma(0.5 * nu) * size_vec;
120 lp -= (0.5 * size_y) *
log(nu) * size_vec;
124 lp -= (0.5 * size_y) *
LOG_PI * size_vec;
133 lp_type sum_lp_vec(0.0);
134 for (
size_t i = 0; i < size_vec; i++) {
135 Eigen::Matrix<typename return_type<T_y, T_loc>::type,
136 Eigen::Dynamic, 1> y_minus_mu(size_y);
137 for (
int j = 0; j < size_y; j++)
138 y_minus_mu(j) = y_vec[i](j)-mu_vec[i](j);
142 lp -= 0.5 * (nu + size_y) * sum_lp_vec;
147 template <
typename T_y,
typename T_dof,
typename T_loc,
typename T_scale>
151 const T_scale& Sigma) {
152 return multi_student_t_log<false>(y, nu, mu, Sigma);
boost::enable_if_c<!stan::is_var< T1 >::value &&!stan::is_var< T2 >::value, typename boost::math::tools::promote_args< T1, T2 >::type >::type trace_inv_quad_form_ldlt(const LDLT_factor< T1, R2, C2 > &A, const Eigen::Matrix< T2, R3, C3 > &B)
fvar< T > lgamma(const fvar< T > &x)
size_t max_size_mvt(const T1 &x1, const T2 &x2)
bool check_not_nan(const char *function, const char *name, const T_y &y)
Return true if y is not NaN.
fvar< T > log(const fvar< T > &x)
scalar_type_helper< is_vector< T >::value, T >::type type
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
bool check_ldlt_factor(const char *function, const char *name, LDLT_factor< T, R, C > &A)
Return true if the argument is a valid LDLT_factor.
boost::math::tools::promote_args< typename scalar_type< T1 >::type, typename scalar_type< T2 >::type, typename scalar_type< T3 >::type, typename scalar_type< T4 >::type, typename scalar_type< T5 >::type, typename scalar_type< T6 >::type >::type type
return_type< T_y, T_dof, T_loc, T_scale >::type multi_student_t_log(const T_y &y, const T_dof &nu, const T_loc &mu, const T_scale &Sigma)
Return the log of the multivariate Student t distribution at the specified arguments.
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
bool check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Return true if the provided sizes match.
int is_inf(const fvar< T > &x)
Returns 1 if the input's value is infinite and 0 otherwise.
return_type< T_y, T_loc, T_covar >::type multi_normal_log(const T_y &y, const T_loc &mu, const T_covar &Sigma)
bool check_finite(const char *function, const char *name, const T_y &y)
Return true if y is finite.
VectorViewMvt is a template expression that wraps either an Eigen::Matrix or a std::vector
fvar< T > log1p(const fvar< T > &x)
bool check_symmetric(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is symmetric.
size_t length_mvt(const Eigen::Matrix< T, R, C > &)
T log_determinant_ldlt(LDLT_factor< T, R, C > &A)