Stan Math Library  2.8.0
reverse mode automatic differentiation
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Groups
lkj_corr_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_LOG_HPP
3 
50 
51 namespace stan {
52  namespace math {
53 
54  template <typename T_shape>
55  T_shape do_lkj_constant(const T_shape& eta, const unsigned int& K) {
56  using stan::math::sum;
57  using stan::math::lgamma;
58 
59  // Lewandowski, Kurowicka, and Joe (2009) theorem 5
60  T_shape constant;
61  const int Km1 = K - 1;
62  if (eta == 1.0) {
63  // C++ integer division is appropriate in this block
64  Eigen::VectorXd numerator(Km1 / 2);
65  for (int k = 1; k <= numerator.rows(); k++)
66  numerator(k-1) = lgamma(2 * k);
67  constant = sum(numerator);
68  if ( (K % 2) == 1 )
69  constant += 0.25 * (K * K - 1) * LOG_PI
70  - 0.25 * (Km1 * Km1) * LOG_TWO - Km1 * lgamma((K + 1) / 2);
71  else
72  constant += 0.25 * K * (K - 2) * LOG_PI
73  + 0.25 * (3 * K * K - 4 * K) * LOG_TWO
74  + K * lgamma(K / 2) - Km1 * lgamma(K);
75  } else {
76  constant = -Km1 * lgamma(eta + 0.5 * Km1);
77  for (int k = 1; k <= Km1; k++)
78  constant += 0.5 * k * LOG_PI + lgamma(eta + 0.5 * (Km1 - k));
79  }
80  return constant;
81  }
82 
83  // LKJ_Corr(y|eta) [ y correlation matrix (not covariance matrix)
84  // eta > 0; eta == 1 <-> uniform]
85  template <bool propto,
86  typename T_y, typename T_shape>
87  typename boost::math::tools::promote_args<T_y, T_shape>::type
88  lkj_corr_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
89  const T_shape& eta) {
90  static const char* function("stan::math::lkj_corr_log");
91 
94  using stan::math::sum;
95  using boost::math::tools::promote_args;
96 
97  typename promote_args<T_y, T_shape>::type lp(0.0);
98  check_positive(function, "Shape parameter", eta);
99  check_corr_matrix(function, "Correlation matrix", y);
100 
101  const unsigned int K = y.rows();
102  if (K == 0)
103  return 0.0;
104 
106  lp += do_lkj_constant(eta, K);
107 
108  if ( (eta == 1.0) &&
110  return lp;
111 
113  return lp;
114 
115  Eigen::Matrix<T_y, Eigen::Dynamic, 1> values =
116  y.ldlt().vectorD().array().log().matrix();
117  lp += (eta - 1.0) * sum(values);
118  return lp;
119  }
120 
121  template <typename T_y, typename T_shape>
122  inline
123  typename boost::math::tools::promote_args<T_y, T_shape>::type
124  lkj_corr_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
125  const T_shape& eta) {
126  return lkj_corr_log<false>(y, eta);
127  }
128 
129  }
130 }
131 #endif
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition: sum.hpp:20
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Definition: is_constant.hpp:22
Metaprogram structure to determine the base scalar type of a template argument.
Definition: scalar_type.hpp:37
boost::math::tools::promote_args< T_y, T_shape >::type lkj_corr_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const T_shape &eta)
fvar< T > lgamma(const fvar< T > &x)
Definition: lgamma.hpp:15
const double LOG_PI
Definition: constants.hpp:170
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
const double LOG_TWO
Definition: constants.hpp:177
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)
bool check_corr_matrix(const char *function, const char *name, const Eigen::Matrix< T_y, Dynamic, Dynamic > &y)
Return true if the specified matrix is a valid correlation matrix.

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