Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
lkj_corr.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_CORR_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_CORR_HPP
3 
7 #include <stan/prob/traits.hpp>
10 
11 namespace stan {
12  namespace prob {
13 
14  template <typename T_shape>
15  T_shape do_lkj_constant(const T_shape& eta, const unsigned int& K) {
16 
17  using stan::math::sum;
18  using stan::math::lgamma;
19 
20  // Lewandowski, Kurowicka, and Joe (2009) theorem 5
21  T_shape constant;
22  const int Km1 = K - 1;
23  if (eta == 1.0) {
24  // C++ integer division is appropriate in this block
25  Eigen::VectorXd numerator( Km1 / 2 );
26  for(size_t k = 1; k <= numerator.rows(); k++)
27  numerator(k-1) = lgamma(2 * k);
28  constant = sum(numerator);
29  if ( (K % 2) == 1 ) constant += 0.25 * (K * K - 1) * LOG_PI -
30  0.25 * (Km1 * Km1) * LOG_TWO - Km1 * lgamma( (K + 1) / 2);
31  else constant += 0.25 * K * (K - 2) * LOG_PI +
32  0.25 * (3 * K * K - 4 * K) * LOG_TWO +
33  K * lgamma(K / 2) - Km1 * lgamma(K);
34  }
35  else {
36  constant = -Km1 * lgamma(eta + 0.5 * Km1);
37  for (size_t k = 1; k <= Km1; k++)
38  constant += 0.5 * k * LOG_PI + lgamma(eta + 0.5 * (Km1 - k));
39  }
40  return constant;
41  }
42 
43  // LKJ_Corr(L|eta) [ L Cholesky factor of correlation matrix
44  // eta > 0; eta == 1 <-> uniform]
45  template <bool propto,
46  typename T_covar, typename T_shape>
47  typename boost::math::tools::promote_args<T_covar, T_shape>::type
49  const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L,
50  const T_shape& eta) {
51 
52  static const char* function
53  = "stan::prob::lkj_corr_cholesky_log(%1%)";
54 
55  using boost::math::tools::promote_args;
58  using stan::math::sum;
59 
60  typename promote_args<T_covar,T_shape>::type lp(0.0);
61  check_positive(function, eta, "Shape parameter", &lp);
62  check_lower_triangular(function, L, "Random variable", &lp);
63 
64  const unsigned int K = L.rows();
65  if (K == 0)
66  return 0.0;
67 
69  lp += do_lkj_constant(eta, K);
71  const int Km1 = K - 1;
72  Eigen::Matrix<T_covar,Eigen::Dynamic,1> log_diagonals =
73  L.diagonal().tail(Km1).array().log();
74  Eigen::Matrix<T_covar,Eigen::Dynamic,1> values(Km1);
75  for (size_t k = 0; k < Km1; k++)
76  values(k) = (Km1 - k - 1) * log_diagonals(k);
77  if ( (eta == 1.0) &&
79  lp += sum(values);
80  return(lp);
81  }
82  values += (2.0 * eta - 2.0) * log_diagonals;
83  lp += sum(values);
84  }
85 
86  return lp;
87  }
88 
89  template <typename T_covar, typename T_shape>
90  inline
91  typename boost::math::tools::promote_args<T_covar, T_shape>::type
93  const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L,
94  const T_shape& eta) {
95  return lkj_corr_cholesky_log<false>(L,eta);
96  }
97 
98  // LKJ_Corr(y|eta) [ y correlation matrix (not covariance matrix)
99  // eta > 0; eta == 1 <-> uniform]
100  template <bool propto,
101  typename T_y, typename T_shape>
102  typename boost::math::tools::promote_args<T_y, T_shape>::type
103  lkj_corr_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
104  const T_shape& eta) {
105  static const char* function = "stan::prob::lkj_corr_log(%1%)";
106 
111  using stan::math::sum;
112  using boost::math::tools::promote_args;
113 
114  typename promote_args<T_y,T_shape>::type lp(0.0);
115  check_positive(function, eta, "Shape parameter", &lp);
116  check_size_match(function,
117  y.rows(), "Rows of correlation matrix",
118  y.cols(), "columns of correlation matrix",
119  &lp);
120  check_not_nan(function, y, "Correlation matrix", &lp);
121  check_corr_matrix(function, y, "Correlation matrix", &lp);
122 
123  const unsigned int K = y.rows();
124  if (K == 0)
125  return 0.0;
126 
128  lp += do_lkj_constant(eta, K);
129 
130  if ( (eta == 1.0) &&
132  return lp;
133 
135  return lp;
136 
137  Eigen::Matrix<T_y,Eigen::Dynamic,1> values =
138  y.ldlt().vectorD().array().log().matrix();
139  lp += (eta - 1.0) * sum(values);
140  return lp;
141  }
142 
143  template <typename T_y, typename T_shape>
144  inline
145  typename boost::math::tools::promote_args<T_y, T_shape>::type
146  lkj_corr_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
147  const T_shape& eta) {
148  return lkj_corr_log<false>(y,eta);
149  }
150 
151  template <class RNG>
152  inline Eigen::MatrixXd
153  lkj_corr_cholesky_rng(const size_t K,
154  const double eta,
155  RNG& rng) {
156  static const char* function
157  = "stan::prob::lkj_corr_cholesky_rng(%1%)";
158 
160 
161  check_positive(function, eta, "Shape parameter", (double*)0);
162 
163  Eigen::ArrayXd CPCs( (K * (K - 1)) / 2 );
164  double alpha = eta + 0.5 * (K - 1);
165  unsigned int count = 0;
166  for (size_t i = 0; i < (K - 1); i++) {
167  alpha -= 0.5;
168  for (size_t j = i + 1; j < K; j++) {
169  CPCs(count) = 2.0 * stan::prob::beta_rng(alpha,alpha,rng) - 1.0;
170  count++;
171  }
172  }
173  return stan::prob::read_corr_L(CPCs, K);
174  }
175 
176  template <class RNG>
177  inline Eigen::MatrixXd
178  lkj_corr_rng(const size_t K,
179  const double eta,
180  RNG& rng) {
181 
182  static const char* function
183  = "stan::prob::lkj_corr_rng(%1%)";
184 
186 
187  check_positive(function, eta, "Shape parameter", (double*)0);
188 
191  lkj_corr_cholesky_rng(K, eta, rng) );
192  }
193 
194  }
195 }
196 #endif
Eigen::MatrixXd lkj_corr_cholesky_rng(const size_t K, const double eta, RNG &rng)
Definition: lkj_corr.hpp:153
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Definition: traits.hpp:43
Metaprogram structure to determine the base scalar type of a template argument.
Definition: traits.hpp:138
boost::math::tools::promote_args< T_covar, T_shape >::type lkj_corr_cholesky_log(const Eigen::Matrix< T_covar, Eigen::Dynamic, Eigen::Dynamic > &L, const T_shape &eta)
Definition: lkj_corr.hpp:48
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)
Definition: lkj_corr.hpp:15
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)
Definition: lkj_corr.hpp:103
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
Definition: sum.hpp:14
fvar< T > lgamma(const fvar< T > &x)
Definition: lgamma.hpp:15
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Definition: traits.hpp:35
Eigen::MatrixXd lkj_corr_rng(const size_t K, const double eta, RNG &rng)
Definition: lkj_corr.hpp:178
double lgamma(double x)
Definition: lgamma.hpp:31
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > read_corr_L(const Eigen::Array< T, Eigen::Dynamic, 1 > &CPCs, const size_t K)
Return the Cholesky factor of the correlation matrix of the specified dimensionality corresponding to...
Definition: transform.hpp:139
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)
bool check_corr_matrix(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is a valid correlation matrix.
double beta_rng(const double alpha, const double beta, RNG &rng)
Definition: beta.hpp:593
double sum(std::vector< double > &x)
Definition: sum.hpp:10
bool check_not_nan(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is nan.
Eigen::Matrix< fvar< T >, R, R > multiply_lower_tri_self_transpose(const Eigen::Matrix< fvar< T >, R, C > &m)
matrix_d multiply_lower_tri_self_transpose(const matrix_d &L)
Returns the result of multiplying the lower triangular portion of the input matrix by its own transpo...
bool check_lower_triangular(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is lower triangular.

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