1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__WISHART_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__WISHART_HPP
4 #include <boost/concept_check.hpp>
59 template <
bool propto,
60 typename T_y,
typename T_dof,
typename T_scale>
61 typename boost::math::tools::promote_args<T_y,T_dof,T_scale>::type
62 wishart_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
64 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
65 static const char*
function =
"stan::prob::wishart_log(%1%)";
67 using boost::math::tools::promote_args;
80 typename index_type<Matrix<T_scale,Dynamic,Dynamic> >::type k
82 typename promote_args<T_y,T_dof,T_scale>::type lp(0.0);
84 "Degrees of freedom parameter", &lp);
86 W.rows(),
"Rows of random variable",
87 W.cols(),
"columns of random variable",
90 S.rows(),
"Rows of scale parameter",
91 S.cols(),
"columns of scale parameter",
94 W.rows(),
"Rows of random variable",
95 S.rows(),
"columns of scale parameter",
99 LDLT_factor<T_y,Eigen::Dynamic,Eigen::Dynamic> ldlt_W(W);
101 "LDLT_Factor of random variable",&lp))
104 LDLT_factor<T_scale,Eigen::Dynamic,Eigen::Dynamic> ldlt_S(S);
106 "LDLT_Factor of scale parameter",&lp))
112 lp += nu * k * NEG_LOG_TWO_OVER_TWO;
121 Matrix<typename promote_args<T_y,T_scale>::type,Dynamic,Dynamic>
123 static_cast<Matrix<T_y,Dynamic,Dynamic>
>(W.template selfadjointView<Lower>())));
124 lp -= 0.5 *
trace(Sinv_W);
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
135 wishart_log(
const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& W,
137 const Eigen::Matrix<T_scale,Eigen::Dynamic,Eigen::Dynamic>& S) {
138 return wishart_log<false>(W,nu,S);
142 inline Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>
144 const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
147 using Eigen::MatrixXd;
152 static const char*
function =
"stan::prob::wishart_rng(%1%)";
154 typename index_type<MatrixXd>::type k = S.rows();
158 S.rows(),
"Rows of scale parameter",
159 S.cols(),
"columns of scale parameter",
162 MatrixXd B = MatrixXd::Zero(k, k);
164 for (
int j = 0; j < k; ++j) {
165 for (
int i = 0; i < j; ++i)
bool check_greater(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result)
boost::math::tools::promote_args< T_y, T_dof, T_scale >::type 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 Wishart density for the given W, degrees of freedom, and scale matrix.
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< 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...
fvar< T > sqrt(const fvar< T > &x)
double normal_rng(const double mu, const double sigma, RNG &rng)
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
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)
matrix_d crossprod(const matrix_d &M)
Returns the result of pre-multiplying a matrix by its own transpose.
double chi_square_rng(const double nu, RNG &rng)