1 #ifndef __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP__
2 #define __STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP__
45 template <
bool propto,
46 typename T_y,
typename T_dof,
typename T_scale,
48 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
51 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S,
53 static const char*
function =
"stan::prob::inv_wishart_log(%1%)";
57 using boost::math::tools::promote_args;
60 typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
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%)",
79 Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> L = LLT_W.matrixL();
89 lp += nu * S.llt().matrixLLT().diagonal().array().log().sum();
92 lp -= (nu + k + 1.0) * L.diagonal().array().log().sum();
96 Eigen::Matrix<T_y,Eigen::Dynamic,1> W_inv_vec = Eigen::Map<
97 const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic> >(
99 Eigen::Matrix<T_scale,Eigen::Dynamic,1> S_vec = Eigen::Map<
100 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic> >(
105 lp += nu * k * NEG_LOG_TWO_OVER_TWO;
109 template <
bool propto,
110 typename T_y,
typename T_dof,
typename T_scale>
112 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
115 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
120 template <
typename T_y,
typename T_dof,
typename T_scale,
123 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
126 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S,
128 return inv_wishart_log<false>(W,nu,S,Policy());
132 template <
typename T_y,
typename T_dof,
typename T_scale>
134 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
137 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
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.
matrix_v crossprod(const matrix_v &M)
Returns the result of pre-multiplying a matrix by its own transpose.
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)
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.
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.
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.
Probability, optimization and sampling library.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...