Stan Math Library  2.11.0
reverse mode automatic differentiation
inv_wishart_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LOG_HPP
3 
11 
16 
17 namespace stan {
18  namespace math {
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
55  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
56  static const char* function("stan::math::inv_wishart_log");
57 
58  using boost::math::tools::promote_args;
59  using Eigen::Dynamic;
60  using Eigen::Matrix;
65 
67  = S.rows();
68  typename promote_args<T_y, T_dof, T_scale>::type lp(0.0);
69 
70  check_greater(function, "Degrees of freedom parameter", nu, k-1);
71  check_square(function, "random variable", W);
72  check_square(function, "scale parameter", S);
73  check_size_match(function,
74  "Rows of random variable", W.rows(),
75  "columns of scale parameter", S.rows());
76 
77  // FIXME: domain checks
78 
79  using stan::math::lmgamma;
82  using stan::math::trace;
85 
87  check_ldlt_factor(function, "LDLT_Factor of random variable", ldlt_W);
89  check_ldlt_factor(function, "LDLT_Factor of scale parameter", ldlt_S);
90 
92  lp -= lmgamma(k, 0.5 * nu);
94  lp += 0.5 * nu * log_determinant_ldlt(ldlt_S);
95  }
97  lp -= 0.5 * (nu + k + 1.0) * log_determinant_ldlt(ldlt_W);
98  }
100  // L = crossprod(mdivide_left_tri_low(L));
101  // Eigen::Matrix<T_y, Eigen::Dynamic, 1> W_inv_vec = Eigen::Map<
102  // const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic> >(
103  // &L(0), L.size(), 1);
104  // Eigen::Matrix<T_scale, Eigen::Dynamic, 1> S_vec = Eigen::Map<
105  // const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic> >(
106  // &S(0), S.size(), 1);
107  // lp -= 0.5 * dot_product(S_vec, W_inv_vec); // trace(S * W^-1)
108  Eigen::Matrix<typename promote_args<T_y, T_scale>::type,
109  Eigen::Dynamic, Eigen::Dynamic>
110  Winv_S(mdivide_left_ldlt
111  (ldlt_W,
112  static_cast<Eigen::Matrix
113  <T_scale, Eigen::Dynamic, Eigen::Dynamic> >
114  (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
128  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
129  return inv_wishart_log<false>(W, nu, S);
130  }
131 
132  }
133 }
134 #endif
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:19
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.
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...
bool check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Return true if the provided sizes match.
T trace(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m)
Returns the trace of the specified matrix.
Definition: trace.hpp:20
const double NEG_LOG_TWO_OVER_TWO
Definition: constants.hpp:191
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
bool check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is square.
bool check_greater(const char *function, const char *name, const T_y &y, const T_low &low)
Return true if y is strictly greater than low.
bool check_ldlt_factor(const char *function, const char *name, stan::math::LDLT_factor< T, R, C > &A)
Return true if the argument is a valid stan::math::LDLT_factor.

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