Stan  1.0
probability, sampling & optimization
inv_wishart.hpp
Go to the documentation of this file.
1 #ifndef __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP__
2 #define __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP__
3 
8 #include <stan/prob/traits.hpp>
9 #include <stan/meta/traits.hpp>
10 #include <stan/agrad/agrad.hpp>
11 #include <stan/agrad/matrix.hpp>
12 
13 namespace stan {
14  namespace prob {
15  // InvWishart(Sigma|n,Omega) [W, S symmetric, non-neg, definite;
16  // W.dims() = S.dims();
17  // n > S.rows() - 1]
45  template <bool propto,
46  typename T_y, typename T_dof, typename T_scale,
47  class Policy>
48  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
49  inv_wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
50  const T_dof& nu,
51  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S,
52  const Policy&) {
53  static const char* function = "stan::prob::inv_wishart_log(%1%)";
54 
57  using boost::math::tools::promote_args;
58 
60  typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
61  if(!check_greater_or_equal(function, nu, k-1, "Degrees of freedom parameter",
62  &lp, Policy()))
63  return lp;
64  if (!check_size_match(function, W.rows(), W.cols(), &lp, Policy()))
65  return lp;
66  if (!check_size_match(function, S.rows(), S.cols(), &lp, Policy()))
67  return lp;
68  if (!check_size_match(function, W.rows(), S.rows(), &lp, Policy()))
69  return lp;
70  // FIXME: domain checks
71 
72  Eigen::LLT< Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> > LLT_W = W.llt();
73  if (LLT_W.info() != Eigen::Success) {
74  lp = stan::math::policies::raise_domain_error<T_y>(function,
75  "W is not positive definite (%1%)",
76  0,Policy());
77  return lp;
78  }
79  Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> L = LLT_W.matrixL();
80 
84  using stan::math::lmgamma;
85 
87  lp -= lmgamma(k, 0.5 * nu);
89  lp += nu * S.llt().matrixLLT().diagonal().array().log().sum();
90  }
92  lp -= (nu + k + 1.0) * L.diagonal().array().log().sum();
93  }
96  Eigen::Matrix<T_y,Eigen::Dynamic,1> W_inv_vec = Eigen::Map<
97  const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> >(
98  &L(0), L.size(), 1);
99  Eigen::Matrix<T_scale,Eigen::Dynamic,1> S_vec = Eigen::Map<
100  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic> >(
101  &S(0), S.size(), 1);
102  lp -= 0.5 * dot_product(S_vec, W_inv_vec); // trace(S * W^-1)
103  }
105  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
106  return lp;
107  }
108 
109  template <bool propto,
110  typename T_y, typename T_dof, typename T_scale>
111  inline
112  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
113  inv_wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
114  const T_dof& nu,
115  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
116  return inv_wishart_log<propto>(W,nu,S,stan::math::default_policy());
117  }
118 
119 
120  template <typename T_y, typename T_dof, typename T_scale,
121  class Policy>
122  inline
123  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
124  inv_wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
125  const T_dof& nu,
126  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S,
127  const Policy&) {
128  return inv_wishart_log<false>(W,nu,S,Policy());
129  }
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  inv_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 inv_wishart_log<false>(W,nu,S,stan::math::default_policy());
139  }
140 
141 
142  }
143 }
144 #endif
internal::traits< Derived >::Index size_type
var dot_product(const Eigen::Matrix< var, R1, C1 > &v1, const Eigen::Matrix< var, R2, C2 > &v2)
Returns the dot product.
Definition: matrix.hpp:716
matrix_v crossprod(const matrix_v &M)
Returns the result of pre-multiplying a matrix by its own transpose.
Definition: matrix.hpp:1171
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R1, C2 > mdivide_left_tri_low(const Eigen::Matrix< T1, R1, C1 > &A, const Eigen::Matrix< T2, R2, C2 > &b)
Definition: matrix.hpp:1708
boost::math::tools::promote_args< T >::type lmgamma(unsigned int k, T x)
Return the natural logarithm of the multivariate gamma function with the speciifed dimensions and arg...
bool check_greater_or_equal(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result, const Policy &)
matrix_d crossprod(const matrix_d &M)
Returns the result of pre-multiplying a matrix by its own transpose.
Definition: matrix.hpp:1565
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: matrix.hpp:875
boost::math::policies::policy default_policy
Default error-handling policy from Boost.
bool check_size_match(const char *function, T_size1 i, T_size2 j, T_result *result, const Policy &)
boost::math::tools::promote_args< T_y, T_dof, T_scale >::type inv_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, const Policy &)
The log of the Inverse-Wishart density for the given W, degrees of freedom, and scale matrix.
Definition: inv_wishart.hpp:49
Probability, optimization and sampling library.
Definition: agrad.cpp:6
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Definition: traits.hpp:33

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