Stan Math Library  2.12.0
reverse mode automatic differentiation
inv_wishart_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_INV_WISHART_LOG_HPP
3 
15 
16 namespace stan {
17  namespace math {
45  template <bool propto,
46  typename T_y, typename T_dof, typename T_scale>
47  typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
48  inv_wishart_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
49  const T_dof& nu,
50  const Eigen::Matrix
51  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
52  static const char* function("inv_wishart_log");
53 
54  using boost::math::tools::promote_args;
55  using Eigen::Dynamic;
56  using Eigen::Matrix;
57 
59  = S.rows();
60  typename promote_args<T_y, T_dof, T_scale>::type lp(0.0);
61 
62  check_greater(function, "Degrees of freedom parameter", nu, k-1);
63  check_square(function, "random variable", W);
64  check_square(function, "scale parameter", S);
65  check_size_match(function,
66  "Rows of random variable", W.rows(),
67  "columns of scale parameter", S.rows());
68 
70  check_ldlt_factor(function, "LDLT_Factor of random variable", ldlt_W);
72  check_ldlt_factor(function, "LDLT_Factor of scale parameter", ldlt_S);
73 
75  lp -= lmgamma(k, 0.5 * nu);
77  lp += 0.5 * nu * log_determinant_ldlt(ldlt_S);
78  }
80  lp -= 0.5 * (nu + k + 1.0) * log_determinant_ldlt(ldlt_W);
81  }
83  // L = crossprod(mdivide_left_tri_low(L));
84  // Eigen::Matrix<T_y, Eigen::Dynamic, 1> W_inv_vec = Eigen::Map<
85  // const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic> >(
86  // &L(0), L.size(), 1);
87  // Eigen::Matrix<T_scale, Eigen::Dynamic, 1> S_vec = Eigen::Map<
88  // const Eigen::Matrix<T_scale, Eigen::Dynamic, Eigen::Dynamic> >(
89  // &S(0), S.size(), 1);
90  // lp -= 0.5 * dot_product(S_vec, W_inv_vec); // trace(S * W^-1)
91  Eigen::Matrix<typename promote_args<T_y, T_scale>::type,
92  Eigen::Dynamic, Eigen::Dynamic>
93  Winv_S(mdivide_left_ldlt
94  (ldlt_W,
95  static_cast<Eigen::Matrix
96  <T_scale, Eigen::Dynamic, Eigen::Dynamic> >
97  (S.template selfadjointView<Eigen::Lower>())));
98  lp -= 0.5*trace(Winv_S);
99  }
101  lp += nu * k * NEG_LOG_TWO_OVER_TWO;
102  return lp;
103  }
104 
105  template <typename T_y, typename T_dof, typename T_scale>
106  inline
107  typename boost::math::tools::promote_args<T_y, T_dof, T_scale>::type
108  inv_wishart_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& W,
109  const T_dof& nu,
110  const Eigen::Matrix
111  <T_scale, Eigen::Dynamic, Eigen::Dynamic>& S) {
112  return inv_wishart_log<false>(W, nu, S);
113  }
114 
115  }
116 }
117 #endif
Eigen::Matrix< fvar< T2 >, R1, C2 > mdivide_left_ldlt(const 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.
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
bool check_ldlt_factor(const char *function, const char *name, LDLT_factor< T, R, C > &A)
Return true if the argument is a valid LDLT_factor.
Primary template class for the metaprogram to compute the index type of a container.
Definition: index_type.hpp:18
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...
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:19
const double NEG_LOG_TWO_OVER_TWO
Definition: constants.hpp:188
fvar< typename stan::return_type< T, int >::type > lmgamma(int x1, const fvar< T > &x2)
Definition: lmgamma.hpp:15
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.
T log_determinant_ldlt(LDLT_factor< T, R, C > &A)
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.

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