Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
multi_student_t.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTI_STUDENT_T_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTI_STUDENT_T_HPP
3 
4 #include <cstdlib>
5 
6 #include <boost/math/special_functions/gamma.hpp>
7 #include <boost/math/special_functions/fpclassify.hpp>
8 
9 #include <stan/math/matrix.hpp>
15 #include <stan/prob/constants.hpp>
16 #include <stan/prob/traits.hpp>
20 #include <boost/random/variate_generator.hpp>
21 
22 namespace stan {
23 
24  namespace prob {
25 
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
35  multi_student_t_log(const T_y& y,
36  const T_dof& nu,
37  const T_loc& mu,
38  const
39  Eigen::Matrix<T_scale,
40  Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
41  static const char* function = "stan::prob::multi_student_t(%1%)";
42 
48  using boost::math::tools::promote_args;
49  using boost::math::lgamma;
53 
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;
55  lp_type lp(0.0);
56 
57  // allows infinities
58  check_not_nan(function, nu,
59  "Degrees of freedom parameter", &lp);
60  check_positive(function, nu,
61  "Degrees of freedom parameter", &lp);
62 
63  using boost::math::isinf;
64 
65  if (isinf(nu)) // already checked nu > 0
66  return multi_normal_log(y,mu,Sigma);
67 
68  using Eigen::Matrix;
69  using Eigen::Dynamic;
70  using std::vector;
71  VectorViewMvt<const T_y> y_vec(y);
72  VectorViewMvt<const T_loc> mu_vec(mu);
73  //size of std::vector of Eigen vectors
74  size_t size_vec = max_size_mvt(y, mu);
75 
76 
77  //Check if every vector of the array has the same size
78  int size_y = y_vec[0].size();
79  int size_mu = mu_vec[0].size();
80  if (size_vec > 1) {
81  int size_y_old = size_y;
82  int size_y_new;
83  for (size_t i = 1, size_ = length_mvt(y); i < size_; i++) {
84  int size_y_new = y_vec[i].size();
85  check_size_match(function,
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",
88  &lp);
89  size_y_old = size_y_new;
90  }
91  int size_mu_old = size_mu;
92  int size_mu_new;
93  for (size_t i = 1, size_ = length_mvt(mu); i < size_; i++) {
94  int size_mu_new = mu_vec[i].size();
95  check_size_match(function,
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",
98  &lp);
99  size_mu_old = size_mu_new;
100  }
101  (void) size_y_old;
102  (void) size_y_new;
103  (void) size_mu_old;
104  (void) size_mu_new;
105  }
106 
107 
108  check_size_match(function,
109  size_y, "Size of random variable",
110  size_mu, "size of location parameter",
111  &lp);
112  check_size_match(function,
113  size_y, "Size of random variable",
114  Sigma.rows(), "rows of scale parameter",
115  &lp);
116  check_size_match(function,
117  size_y, "Size of random variable",
118  Sigma.cols(), "columns of scale parameter",
119  &lp);
120 
121  for (size_t i = 0; i < size_vec; i++) {
122  check_finite(function, mu_vec[i], "Location parameter", &lp);
123  check_not_nan(function, y_vec[i], "Random variable", &lp);
124  }
125  check_symmetric(function, Sigma, "Scale parameter", &lp);
126 
127 
128  LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_Sigma(Sigma);
129  check_ldlt_factor(function,ldlt_Sigma,"LDLT_Factor of scale parameter",
130  &lp);
131 
132  if (size_y == 0) //y_vec[0].size() == 0
133  return lp;
134 
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;
139  }
140 
142  lp -= (0.5 * size_y) * LOG_PI * size_vec;
143 
144  using stan::math::multiply;
146  using stan::math::subtract;
147  using Eigen::Array;
148 
149 
151  lp -= 0.5 * log_determinant_ldlt(ldlt_Sigma) * size_vec;
152  }
153 
155  lp_type sum_lp_vec(0.0);
156  for (size_t i = 0; i < size_vec; i++) {
157  Matrix<typename
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);
162  sum_lp_vec += log(1.0 + trace_inv_quad_form_ldlt(ldlt_Sigma,y_minus_mu) / nu);
163  }
164  lp -= 0.5 * (nu + size_y) * sum_lp_vec;
165  }
166  return lp;
167  }
168 
169  template <typename T_y, typename T_dof, typename T_loc, typename T_scale>
170  inline
171  typename boost::math::tools::promote_args<typename scalar_type<T_y>::type,T_dof,typename scalar_type<T_loc>::type,T_scale>::type
172  multi_student_t_log(const T_y& y,
173  const T_dof& nu,
174  const T_loc& mu,
175  const
176  Eigen::Matrix<T_scale,
177  Eigen::Dynamic,Eigen::Dynamic>& Sigma) {
178  return multi_student_t_log<false>(y,nu,mu,Sigma);
179  }
180 
181 
182  template <class RNG>
183  inline Eigen::VectorXd
184  multi_student_t_rng(const double nu,
185  const Eigen::Matrix<double,Eigen::Dynamic,1>& mu,
186  const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& s,
187  RNG& rng) {
188 
189  static const char* function = "stan::prob::multi_student_t_rng(%1%)";
190 
195 
196  check_finite(function, mu, "Location parameter", (double*)0);
197  check_symmetric(function, s, "Scale parameter", (double*)0);
198  check_not_nan(function, nu,
199  "Degrees of freedom parameter", (double*)0);
200  check_positive(function, nu,
201  "Degrees of freedom parameter", (double*)0);
202 
203  Eigen::VectorXd z(s.cols());
204  z.setZero();
205 
206  double w = stan::prob::inv_gamma_rng(nu / 2, nu / 2, rng);
207  return mu + std::sqrt(w) * stan::prob::multi_normal_rng(z, s, rng);
208  }
209  }
210 }
211 #endif
size_t max_size_mvt(const T1 &x1, const T2 &x2)
Definition: traits.hpp:575
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
Definition: traits.hpp:139
fvar< T > lgamma(const fvar< T > &x)
Definition: lgamma.hpp:15
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...
Definition: traits.hpp:35
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...
Definition: subtract.hpp:27
fvar< T > sqrt(const fvar< T > &x)
Definition: sqrt.hpp:15
double inv_gamma_rng(const double alpha, const double beta, RNG &rng)
Definition: inv_gamma.hpp:508
size_t length_mvt(const T &)
Definition: traits.hpp:561
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)
size_t size_
Definition: dot_self.hpp:18
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.
Definition: multiply.hpp:25
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:15
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.
Definition: dot_product.hpp:22

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