1 #ifndef STAN_MATH_PRIM_MAT_PROB_WISHART_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_WISHART_LOG_HPP
55 template <
bool propto,
56 typename T_y,
typename T_dof,
typename T_scale>
57 typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
58 wishart_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
60 const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic>&
62 static const char*
function(
"stan::math::wishart_log");
64 using boost::math::tools::promote_args;
80 typename promote_args<T_y, T_dof, T_scale>::type lp(0.0);
81 check_greater(
function,
"Degrees of freedom parameter", nu, k-1);
85 "Rows of random variable", W.rows(),
86 "columns of scale parameter", S.rows());
111 Matrix<typename promote_args<T_y, T_scale>::type, Dynamic, Dynamic>
114 static_cast<Matrix<T_y, Dynamic, Dynamic>
>
115 (W.template selfadjointView<Lower>())));
116 lp -= 0.5 *
trace(Sinv_W);
124 template <
typename T_y,
typename T_dof,
typename T_scale>
126 typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
127 wishart_log(
const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
130 <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
131 return wishart_log<false>(W, nu, S);
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
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.
Primary template class for the metaprogram to compute the index type of a container.
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_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Return true if the provided sizes match.
T trace(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m)
Returns the trace of the specified matrix.
const double NEG_LOG_TWO_OVER_TWO
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)
bool check_square(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is square.
bool check_greater(const char *function, const char *name, const T_y &y, const T_low &low)
Return true if y is strictly greater than low.
bool check_ldlt_factor(const char *function, const char *name, stan::math::LDLT_factor< T, R, C > &A)
Return true if the argument is a valid stan::math::LDLT_factor.