Stan Math Library  2.14.0
reverse mode automatic differentiation
wishart_lpdf.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_WISHART_LPDF_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_WISHART_LPDF_HPP
3 
20 
21 namespace stan {
22  namespace math {
23 
51  template <bool propto,
52  typename T_y, typename T_dof, typename T_scale>
53  typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
54  wishart_lpdf(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
55  const T_dof& nu,
56  const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic>&
57  S) {
58  static const char* function("wishart_lpdf");
59 
60  using boost::math::tools::promote_args;
61  using Eigen::Dynamic;
62  using Eigen::Lower;
63  using Eigen::Matrix;
64 
66  = W.rows();
67  typename promote_args<T_y, T_dof, T_scale>::type lp(0.0);
68  check_greater(function, "Degrees of freedom parameter", nu, k-1);
69  check_square(function, "random variable", W);
70  check_square(function, "scale parameter", S);
71  check_size_match(function,
72  "Rows of random variable", W.rows(),
73  "columns of scale parameter", S.rows());
74 
76  check_ldlt_factor(function, "LDLT_Factor of random variable",
77  ldlt_W);
78 
80  check_ldlt_factor(function, "LDLT_Factor of scale parameter",
81  ldlt_S);
82 
84  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
85 
87  lp -= lmgamma(k, 0.5 * nu);
88 
90  lp -= 0.5 * nu * log_determinant_ldlt(ldlt_S);
91 
93  Matrix<typename promote_args<T_y, T_scale>::type, Dynamic, Dynamic>
94  Sinv_W(mdivide_left_ldlt
95  (ldlt_S,
96  static_cast<Matrix<T_y, Dynamic, Dynamic> >
97  (W.template selfadjointView<Lower>())));
98  lp -= 0.5 * trace(Sinv_W);
99  }
100 
101  if (include_summand<propto, T_y, T_dof>::value && nu != (k + 1))
102  lp += 0.5 * (nu - k - 1.0) * log_determinant_ldlt(ldlt_W);
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
109  wishart_lpdf(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
110  const T_dof& nu,
111  const Eigen::Matrix
112  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
113  return wishart_lpdf<false>(W, nu, S);
114  }
115 
116  }
117 }
118 #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.
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
boost::math::tools::promote_args< T_y, T_dof, T_scale >::type 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 Wishart density for the given W, degrees of freedom, and scale matrix.
T log_determinant_ldlt(LDLT_factor< T, R, C > &A)

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