1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__INV_WISHART_HPP
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
54 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
55 static const char*
function =
"stan::prob::inv_wishart_log(%1%)";
57 using boost::math::tools::promote_args;
64 typename index_type<Matrix<T_scale,Dynamic,Dynamic> >::type k
66 typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
68 check_greater(
function, nu, k-1,
"Degrees of freedom parameter",
71 W.rows(),
"Rows of random variable",
72 W.cols(),
"columns of random variable",
75 S.rows(),
"Rows of scale parameter",
76 S.cols(),
"columns of scale parameter",
79 W.rows(),
"Rows of random variable",
80 S.rows(),
"columns of scale parameter",
92 LDLT_factor<T_y,Eigen::Dynamic,Eigen::Dynamic> ldlt_W(W);
94 LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_S(S);
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);
118 lp += nu * k * NEG_LOG_TWO_OVER_TWO;
122 template <
typename T_y,
typename T_dof,
typename T_scale>
124 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
127 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
128 return inv_wishart_log<false>(W,nu,S);
132 inline Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>
134 const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
137 static const char*
function =
"stan::prob::inv_wishart_rng(%1%)";
141 using Eigen::MatrixXd;
144 typename index_type<MatrixXd>::type k = S.rows();
146 check_greater(
function, nu, k-1,
"Degrees of freedom parameter",
149 S.rows(),
"Rows of scale parameter",
150 S.cols(),
"columns of scale parameter", (
double*)0);
152 MatrixXd S_inv = MatrixXd::Identity(k, k);
153 S_inv = S.ldlt().solve(S_inv);
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.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
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...
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...
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > wishart_rng(const double nu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
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.
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)