Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
4 #include <stan/agrad/rev.hpp>
12 #include <stan/meta/traits.hpp>
13 #include <stan/prob/constants.hpp>
14 #include <stan/prob/traits.hpp>
16 
17 namespace stan {
18  namespace prob {
19  // InvWishart(Sigma|n,Omega) [W, S symmetric, non-neg, definite;
20  // W.dims() = S.dims();
21  // n > S.rows() - 1]
49  template <bool propto,
50  typename T_y, typename T_dof, typename T_scale>
51  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
52  inv_wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
53  const T_dof& nu,
54  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
55  static const char* function = "stan::prob::inv_wishart_log(%1%)";
56 
57  using boost::math::tools::promote_args;
58  using Eigen::Dynamic;
59  using Eigen::Matrix;
63 
64  typename index_type<Matrix<T_scale,Dynamic,Dynamic> >::type k
65  = S.rows();
66  typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
67 
68  check_greater(function, nu, k-1, "Degrees of freedom parameter",
69  &lp);
70  check_size_match(function,
71  W.rows(), "Rows of random variable",
72  W.cols(), "columns of random variable",
73  &lp);
74  check_size_match(function,
75  S.rows(), "Rows of scale parameter",
76  S.cols(), "columns of scale parameter",
77  &lp);
78  check_size_match(function,
79  W.rows(), "Rows of random variable",
80  S.rows(), "columns of scale parameter",
81  &lp);
82 
83  // FIXME: domain checks
84 
85  using stan::math::lmgamma;
88  using stan::math::trace;
91 
92  LDLT_factor<T_y,Eigen::Dynamic,Eigen::Dynamic> ldlt_W(W);
93  check_ldlt_factor(function,ldlt_W,"LDLT_Factor of random variable",&lp);
94  LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_S(S);
95  check_ldlt_factor(function,ldlt_S,"LDLT_Factor of scale parameter",&lp);
96 
98  lp -= lmgamma(k, 0.5 * nu);
100  lp += 0.5 * nu * log_determinant_ldlt(ldlt_S);
101  }
103  lp -= 0.5 * (nu + k + 1.0) * log_determinant_ldlt(ldlt_W);
104  }
106 // L = crossprod(mdivide_left_tri_low(L));
107 // Eigen::Matrix<T_y,Eigen::Dynamic,1> W_inv_vec = Eigen::Map<
108 // const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> >(
109 // &L(0), L.size(), 1);
110 // Eigen::Matrix<T_scale,Eigen::Dynamic,1> S_vec = Eigen::Map<
111 // const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic> >(
112 // &S(0), S.size(), 1);
113 // lp -= 0.5 * dot_product(S_vec, W_inv_vec); // trace(S * W^-1)
114  Eigen::Matrix<typename promote_args<T_y,T_scale>::type,Eigen::Dynamic,Eigen::Dynamic> Winv_S(mdivide_left_ldlt(ldlt_W, static_cast<Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic> >(S.template selfadjointView<Eigen::Lower>())));
115  lp -= 0.5*trace(Winv_S);
116  }
118  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
119  return lp;
120  }
121 
122  template <typename T_y, typename T_dof, typename T_scale>
123  inline
124  typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
125  inv_wishart_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
126  const T_dof& nu,
127  const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
128  return inv_wishart_log<false>(W,nu,S);
129  }
130 
131  template <class RNG>
132  inline Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>
133  inv_wishart_rng(const double nu,
134  const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
135  RNG& rng) {
136 
137  static const char* function = "stan::prob::inv_wishart_rng(%1%)";
138 
141  using Eigen::MatrixXd;
143 
144  typename index_type<MatrixXd>::type k = S.rows();
145 
146  check_greater(function, nu, k-1, "Degrees of freedom parameter",
147  (double*)0);
148  check_size_match(function,
149  S.rows(), "Rows of scale parameter",
150  S.cols(), "columns of scale parameter", (double*)0);
151 
152  MatrixXd S_inv = MatrixXd::Identity(k, k);
153  S_inv = S.ldlt().solve(S_inv);
154 
155  return wishart_rng(nu, S_inv, rng).inverse();
156  }
157  }
158 }
159 #endif
bool check_greater(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result)
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< double, Eigen::Dynamic, Eigen::Dynamic > inv_wishart_rng(const double nu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
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
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)
The log of the Inverse-Wishart density for the given W, degrees of freedom, and scale matrix...
Definition: inv_wishart.hpp:52
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

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