Stan Math Library  2.14.0
reverse mode automatic differentiation
inv_wishart_lpdf.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LPDF_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LPDF_HPP
3 
15 
16 namespace stan {
17  namespace math {
45  template <bool propto,
46  typename T_y, typename T_dof, typename T_scale>
47  typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
49  Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
50  const T_dof& nu,
51  const Eigen::Matrix
52  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
53  static const char* function("inv_wishart_lpdf");
54 
55  using boost::math::tools::promote_args;
56  using Eigen::Dynamic;
57  using Eigen::Matrix;
58 
60  = S.rows();
61  typename promote_args<T_y, T_dof, T_scale>::type lp(0.0);
62 
63  check_greater(function, "Degrees of freedom parameter", nu, k-1);
64  check_square(function, "random variable", W);
65  check_square(function, "scale parameter", S);
66  check_size_match(function,
67  "Rows of random variable", W.rows(),
68  "columns of scale parameter", S.rows());
69 
71  check_ldlt_factor(function, "LDLT_Factor of random variable", ldlt_W);
73  check_ldlt_factor(function, "LDLT_Factor of scale parameter", ldlt_S);
74 
76  lp -= lmgamma(k, 0.5 * nu);
78  lp += 0.5 * nu * log_determinant_ldlt(ldlt_S);
79  }
81  lp -= 0.5 * (nu + k + 1.0) * log_determinant_ldlt(ldlt_W);
82  }
84  // L = crossprod(mdivide_left_tri_low(L));
85  // Eigen::Matrix<T_y, Eigen::Dynamic, 1> W_inv_vec = Eigen::Map<
86  // const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic> >(
87  // &L(0), L.size(), 1);
88  // Eigen::Matrix<T_scale, Eigen::Dynamic, 1> S_vec = Eigen::Map<
89  // const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic> >(
90  // &S(0), S.size(), 1);
91  // lp -= 0.5 * dot_product(S_vec, W_inv_vec); // trace(S * W^-1)
92  Eigen::Matrix<typename promote_args<T_y, T_scale>::type,
93  Eigen::Dynamic, Eigen::Dynamic>
94  Winv_S(mdivide_left_ldlt
95  (ldlt_W,
96  static_cast<Eigen::Matrix
97  <T_scale, Eigen::Dynamic, Eigen::Dynamic> >
98  (S.template selfadjointView<Eigen::Lower>())));
99  lp -= 0.5*trace(Winv_S);
100  }
102  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
103  return lp;
104  }
105 
106  template <typename T_y, typename T_dof, typename T_scale>
107  inline
108  typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
110  Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
111  const T_dof& nu,
112  const Eigen::Matrix
113  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
114  return inv_wishart_lpdf<false>(W, nu, S);
115  }
116 
117  }
118 }
119 #endif
Eigen::Matrix< fvar< T2 >, R1, C2 > mdivide_left_ldlt(const 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.
void check_ldlt_factor(const char *function, const char *name, LDLT_factor< T, R, C > &A)
Check if the argument is a valid LDLT_factor.
boost::math::tools::promote_args< T_y, T_dof, T_scale >::type inv_wishart_lpdf(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...
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Check if the provided sizes match.
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:18
void check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Check if the specified matrix is square.
T trace(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m)
Returns the trace of the specified matrix.
Definition: trace.hpp:19
const double NEG_LOG_TWO_OVER_TWO
Definition: constants.hpp:188
void check_greater(const char *function, const char *name, const T_y &y, const T_low &low)
Check if y is strictly greater than low.
fvar< typename stan::return_type< T, int >::type > lmgamma(int x1, const fvar< T > &x2)
Definition: lmgamma.hpp:15
T log_determinant_ldlt(LDLT_factor< T, R, C > &A)

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