Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wishart.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__WISHART_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__WISHART_HPP
3 
4 #include <boost/concept_check.hpp>
5 
20 #include <stan/prob/constants.hpp>
21 #include <stan/prob/traits.hpp>
24 
25 namespace stan {
26 
27  namespace prob {
28 
29  // Wishart(Sigma|n,Omega) [Sigma, Omega symmetric, non-neg, definite;
30  // Sigma.dims() = Omega.dims();
31  // n > Sigma.rows() - 1]
59  template <bool propto,
60  typename T_y, typename T_dof, typename T_scale>
61  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
62  wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
63  const T_dof& nu,
64  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
65  static const char* function = "stan::prob::wishart_log(%1%)";
66 
67  using boost::math::tools::promote_args;
68  using Eigen::Dynamic;
69  using Eigen::Lower;
70  using Eigen::Matrix;
78 
79 
80  typename index_type<Matrix<T_scale,Dynamic,Dynamic> >::type k
81  = W.rows();
82  typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
83  check_greater(function, nu, k-1,
84  "Degrees of freedom parameter", &lp);
85  check_size_match(function,
86  W.rows(), "Rows of random variable",
87  W.cols(), "columns of random variable",
88  &lp);
89  check_size_match(function,
90  S.rows(), "Rows of scale parameter",
91  S.cols(), "columns of scale parameter",
92  &lp);
93  check_size_match(function,
94  W.rows(), "Rows of random variable",
95  S.rows(), "columns of scale parameter",
96  &lp);
97  // FIXME: domain checks
98 
99  LDLT_factor<T_y,Eigen::Dynamic,Eigen::Dynamic> ldlt_W(W);
100  if (!check_ldlt_factor(function,ldlt_W,
101  "LDLT_Factor of random variable",&lp))
102  return lp;
103 
104  LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_S(S);
105  if (!check_ldlt_factor(function,ldlt_S,
106  "LDLT_Factor of scale parameter",&lp))
107  return lp;
108 
109  using stan::math::trace;
110  using stan::math::lmgamma;
112  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
113 
115  lp -= lmgamma(k, 0.5 * nu);
116 
118  lp -= 0.5 * nu * log_determinant_ldlt(ldlt_S);
119 
121  Matrix<typename promote_args<T_y,T_scale>::type,Dynamic,Dynamic>
122  Sinv_W(mdivide_left_ldlt(ldlt_S,
123  static_cast<Matrix<T_y,Dynamic,Dynamic> >(W.template selfadjointView<Lower>())));
124  lp -= 0.5 * trace(Sinv_W);
125  }
126 
127  if (include_summand<propto,T_y,T_dof>::value && nu != (k + 1))
128  lp += 0.5 * (nu - k - 1.0) * log_determinant_ldlt(ldlt_W);
129  return lp;
130  }
131 
132  template <typename T_y, typename T_dof, typename T_scale>
133  inline
134  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
135  wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
136  const T_dof& nu,
137  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
138  return wishart_log<false>(W,nu,S);
139  }
140 
141  template <class RNG>
142  inline Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>
143  wishart_rng(const double nu,
144  const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
145  RNG& rng) {
146 
147  using Eigen::MatrixXd;
151 
152  static const char* function = "stan::prob::wishart_rng(%1%)";
153 
154  typename index_type<MatrixXd>::type k = S.rows();
155 
156  check_positive(function,nu,"degrees of freedom",(double*)0);
157  check_size_match(function,
158  S.rows(), "Rows of scale parameter",
159  S.cols(), "columns of scale parameter",
160  (double*)0);
161 
162  MatrixXd B = MatrixXd::Zero(k, k);
163 
164  for (int j = 0; j < k; ++j) {
165  for (int i = 0; i < j; ++i)
166  B(i, j) = normal_rng(0, 1, rng);
167  B(j,j) = std::sqrt(chi_square_rng(nu - j, rng));
168  }
169 
170  return stan::math::crossprod(B * S.llt().matrixU());
171  }
172 
173 
174  }
175 
176 }
177 #endif
bool check_greater(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result)
boost::math::tools::promote_args< T_y, T_dof, T_scale >::type wishart_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &W, const T_dof &nu, const Eigen::Matrix< T_scale, Eigen::Dynamic, Eigen::Dynamic > &S)
The log of the Wishart density for the given W, degrees of freedom, and scale matrix.
Definition: wishart.hpp:62
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R1, C2 > mdivide_left_ldlt(const stan::math::LDLT_factor< T1, R1, C1 > &A, const Eigen::Matrix< T2, R2, C2 > &b)
Returns the solution of the system Ax=b given an LDLT_factor of A.
Eigen::Matrix< fvar< T2 >, R1, C2 > mdivide_left_ldlt(const stan::math::LDLT_factor< double, R1, C1 > &A, const Eigen::Matrix< fvar< T2 >, R2, C2 > &b)
Returns the solution of the system Ax=b given an LDLT_factor of A.
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.
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:21
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Definition: traits.hpp:35
boost::math::tools::promote_args< T >::type lmgamma(const int k, T x)
Return the natural logarithm of the multivariate gamma function with the speciifed dimensions and arg...
Definition: lmgamma.hpp:57
fvar< T > sqrt(const fvar< T > &x)
Definition: sqrt.hpp:15
double normal_rng(const double mu, const double sigma, RNG &rng)
Definition: normal.hpp:377
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > wishart_rng(const double nu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
Definition: wishart.hpp:143
var log_determinant_ldlt(stan::math::LDLT_factor< var, R, C > &A)
T trace(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m)
Returns the trace of the specified matrix.
Definition: trace.hpp:20
T log_determinant_ldlt(stan::math::LDLT_factor< T, R, C > &A)
fvar< typename stan::return_type< T, int >::type > lmgamma(int x1, const fvar< T > &x2)
Definition: lmgamma.hpp:16
matrix_d crossprod(const matrix_d &M)
Returns the result of pre-multiplying a matrix by its own transpose.
Definition: crossprod.hpp:17
double chi_square_rng(const double nu, RNG &rng)
Definition: chi_square.hpp:442

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