Stan Math Library  2.20.0
reverse mode automatic differentiation
multi_student_t_lpdf.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_MULTI_STUDENT_T_LPDF_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_MULTI_STUDENT_T_LPDF_HPP
3 
17 #include <cmath>
18 #include <cstdlib>
19 
20 namespace stan {
21 namespace math {
22 
29 template <bool propto, typename T_y, typename T_dof, typename T_loc,
30  typename T_scale>
32  const T_y& y, const T_dof& nu, const T_loc& mu, const T_scale& Sigma) {
33  static const char* function = "multi_student_t";
34  using std::log;
35  typedef typename scalar_type<T_scale>::type T_scale_elem;
36  typedef typename return_type<T_y, T_dof, T_loc, T_scale>::type lp_type;
37 
38  check_not_nan(function, "Degrees of freedom parameter", nu);
39  check_positive(function, "Degrees of freedom parameter", nu);
40 
41  if (is_inf(nu))
42  return multi_normal_log(y, mu, Sigma);
43 
44  using Eigen::Matrix;
45  using std::vector;
46 
47  size_t number_of_y = length_mvt(y);
48  size_t number_of_mu = length_mvt(mu);
49  if (number_of_y == 0 || number_of_mu == 0)
50  return 0;
51  check_consistent_sizes_mvt(function, "y", y, "mu", mu);
52 
53  vector_seq_view<T_y> y_vec(y);
54  vector_seq_view<T_loc> mu_vec(mu);
55  size_t size_vec = max_size_mvt(y, mu);
56 
57  int size_y = y_vec[0].size();
58  int size_mu = mu_vec[0].size();
59  if (size_vec > 1) {
60  int size_y_old = size_y;
61  int size_y_new;
62  for (size_t i = 1, size_ = length_mvt(y); i < size_; i++) {
63  int size_y_new = y_vec[i].size();
65  function, "Size of one of the vectors of the random variable",
66  size_y_new, "Size of another vector of the random variable",
67  size_y_old);
68  size_y_old = size_y_new;
69  }
70  int size_mu_old = size_mu;
71  int size_mu_new;
72  for (size_t i = 1, size_ = length_mvt(mu); i < size_; i++) {
73  int size_mu_new = mu_vec[i].size();
74  check_size_match(function,
75  "Size of one of the vectors "
76  "of the location variable",
77  size_mu_new,
78  "Size of another vector of "
79  "the location variable",
80  size_mu_old);
81  size_mu_old = size_mu_new;
82  }
83  (void)size_y_old;
84  (void)size_y_new;
85  (void)size_mu_old;
86  (void)size_mu_new;
87  }
88 
89  check_size_match(function, "Size of random variable", size_y,
90  "size of location parameter", size_mu);
91  check_size_match(function, "Size of random variable", size_y,
92  "rows of scale parameter", Sigma.rows());
93  check_size_match(function, "Size of random variable", size_y,
94  "columns of scale parameter", Sigma.cols());
95 
96  for (size_t i = 0; i < size_vec; i++) {
97  check_finite(function, "Location parameter", mu_vec[i]);
98  check_not_nan(function, "Random variable", y_vec[i]);
99  }
100  check_symmetric(function, "Scale parameter", Sigma);
101 
103  check_ldlt_factor(function, "LDLT_Factor of scale parameter", ldlt_Sigma);
104 
105  if (size_y == 0)
106  return 0;
107 
108  lp_type lp(0);
109 
111  lp += lgamma(0.5 * (nu + size_y)) * size_vec;
112  lp -= lgamma(0.5 * nu) * size_vec;
113  lp -= (0.5 * size_y) * log(nu) * size_vec;
114  }
115 
117  lp -= (0.5 * size_y) * LOG_PI * size_vec;
118 
119  using Eigen::Array;
120 
122  lp -= 0.5 * log_determinant_ldlt(ldlt_Sigma) * size_vec;
123  }
124 
126  lp_type sum_lp_vec(0.0);
127  for (size_t i = 0; i < size_vec; i++) {
128  Eigen::Matrix<typename return_type<T_y, T_loc>::type, Eigen::Dynamic, 1>
129  y_minus_mu(size_y);
130  for (int j = 0; j < size_y; j++)
131  y_minus_mu(j) = y_vec[i](j) - mu_vec[i](j);
132  sum_lp_vec
133  += log1p(trace_inv_quad_form_ldlt(ldlt_Sigma, y_minus_mu) / nu);
134  }
135  lp -= 0.5 * (nu + size_y) * sum_lp_vec;
136  }
137  return lp;
138 }
139 
140 template <typename T_y, typename T_dof, typename T_loc, typename T_scale>
142 multi_student_t_lpdf(const T_y& y, const T_dof& nu, const T_loc& mu,
143  const T_scale& Sigma) {
144  return multi_student_t_lpdf<false>(y, nu, mu, Sigma);
145 }
146 
147 } // namespace math
148 } // namespace stan
149 #endif
void check_finite(const char *function, const char *name, const T_y &y)
Check if y is finite.
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition: lgamma.hpp:21
size_t max_size_mvt(const T1 &x1, const T2 &x2)
const double LOG_PI
Definition: constants.hpp:144
void check_ldlt_factor(const char *function, const char *name, LDLT_factor< T, R, C > &A)
Raise domain error if the specified LDLT factor is invalid.
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:12
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.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
return_type< T_y, T_dof, T_loc, T_scale >::type multi_student_t_lpdf(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.
LDLT_factor is a thin wrapper on Eigen::LDLT to allow for reusing factorizations and efficient autodi...
Definition: LDLT_factor.hpp:63
boost::math::tools::promote_args< double, typename scalar_type< T >::type, typename return_type< Types_pack... >::type >::type type
Definition: return_type.hpp:36
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
std::enable_if< !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)
int is_inf(const fvar< T > &x)
Returns 1 if the input&#39;s value is infinite and 0 otherwise.
Definition: is_inf.hpp:20
return_type< T_y, T_loc, T_covar >::type multi_normal_log(const T_y &y, const T_loc &mu, const T_covar &Sigma)
This class provides a low-cost wrapper for situations where you either need an Eigen Vector or RowVec...
fvar< T > log1p(const fvar< T > &x)
Definition: log1p.hpp:12
void check_consistent_sizes_mvt(const char *function, const char *name1, const T1 &x1, const char *name2, const T2 &x2)
Check if the dimension of x1 is consistent with x2.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
void check_symmetric(const char *function, const char *name, const matrix_cl &y)
Check if the matrix_cl is symmetric.
size_t length_mvt(const Eigen::Matrix< T, R, C > &)
Definition: length_mvt.hpp:12
T log_determinant_ldlt(LDLT_factor< T, R, C > &A)

     [ Stan Home Page ] © 2011–2018, Stan Development Team.