Stan Math Library  2.6.3
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros
wishart_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_WISHART_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_WISHART_LOG_HPP
3 
20 
21 namespace stan {
22 
23  namespace math {
24 
25  // Wishart(Sigma|n, Omega) [Sigma, Omega symmetric, non-neg, definite;
26  // Sigma.dims() = Omega.dims();
27  // n > Sigma.rows() - 1]
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,
59  const T_dof& nu,
60  const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic>&
61  S) {
62  static const char* function("stan::math::wishart_log");
63 
64  using boost::math::tools::promote_args;
65  using Eigen::Dynamic;
66  using Eigen::Lower;
67  using Eigen::Matrix;
76 
77 
79  = W.rows();
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);
82  check_square(function, "random variable", W);
83  check_square(function, "scale parameter", S);
84  check_size_match(function,
85  "Rows of random variable", W.rows(),
86  "columns of scale parameter", S.rows());
87  // FIXME: domain checks
88 
90  if (!check_ldlt_factor(function, "LDLT_Factor of random variable",
91  ldlt_W))
92  return lp;
93 
95  if (!check_ldlt_factor(function, "LDLT_Factor of scale parameter",
96  ldlt_S))
97  return lp;
98 
99  using stan::math::trace;
100  using stan::math::lmgamma;
102  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
103 
105  lp -= lmgamma(k, 0.5 * nu);
106 
108  lp -= 0.5 * nu * log_determinant_ldlt(ldlt_S);
109 
111  Matrix<typename promote_args<T_y, T_scale>::type, Dynamic, Dynamic>
112  Sinv_W(mdivide_left_ldlt
113  (ldlt_S,
114  static_cast<Matrix<T_y, Dynamic, Dynamic> >
115  (W.template selfadjointView<Lower>())));
116  lp -= 0.5 * trace(Sinv_W);
117  }
118 
119  if (include_summand<propto, T_y, T_dof>::value && nu != (k + 1))
120  lp += 0.5 * (nu - k - 1.0) * log_determinant_ldlt(ldlt_W);
121  return lp;
122  }
123 
124  template <typename T_y, typename T_dof, typename T_scale>
125  inline
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,
128  const T_dof& nu,
129  const Eigen::Matrix
130  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
131  return wishart_log<false>(W, nu, S);
132  }
133 
134  }
135 
136 }
137 #endif
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.
Definition: wishart_log.hpp:58
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:19
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.
Definition: trace.hpp:20
const double NEG_LOG_TWO_OVER_TWO
Definition: constants.hpp:191
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)
Definition: lmgamma.hpp:16
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.

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