Stan Math Library  2.9.0
reverse mode automatic differentiation
lkj_corr_cholesky_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_CHOLESKY_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_CORR_CHOLESKY_LOG_HPP
3 
47 
48 namespace stan {
49  namespace math {
50 
51  // LKJ_Corr(L|eta) [ L Cholesky factor of correlation matrix
52  // eta > 0; eta == 1 <-> uniform]
53  template <bool propto,
54  typename T_covar, typename T_shape>
55  typename boost::math::tools::promote_args<T_covar, T_shape>::type
56  lkj_corr_cholesky_log(const Eigen::Matrix
57  <T_covar, Eigen::Dynamic, Eigen::Dynamic>& L,
58  const T_shape& eta) {
59  static const char* function("stan::math::lkj_corr_cholesky_log");
60 
61  using boost::math::tools::promote_args;
64  using stan::math::sum;
65 
66  typedef typename promote_args<T_covar, T_shape>::type lp_ret;
67  lp_ret lp(0.0);
68  check_positive(function, "Shape parameter", eta);
69  check_lower_triangular(function, "Random variable", L);
70 
71  const unsigned int K = L.rows();
72  if (K == 0)
73  return 0.0;
74 
76  lp += do_lkj_constant(eta, K);
78  const int Km1 = K - 1;
79  Eigen::Matrix<T_covar, Eigen::Dynamic, 1> log_diagonals =
80  L.diagonal().tail(Km1).array().log();
81  Eigen::Matrix<lp_ret, Eigen::Dynamic, 1> values(Km1);
82  for (int k = 0; k < Km1; k++)
83  values(k) = (Km1 - k - 1) * log_diagonals(k);
84  if ( (eta == 1.0) &&
86  lp += sum(values);
87  return(lp);
88  }
89  values += multiply(2.0 * eta - 2.0, log_diagonals);
90  lp += sum(values);
91  }
92 
93  return lp;
94  }
95 
96  template <typename T_covar, typename T_shape>
97  inline
98  typename boost::math::tools::promote_args<T_covar, T_shape>::type
99  lkj_corr_cholesky_log(const Eigen::Matrix
100  <T_covar, Eigen::Dynamic, Eigen::Dynamic>& L,
101  const T_shape& eta) {
102  return lkj_corr_cholesky_log<false>(L, eta);
103  }
104 
105  }
106 }
107 #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
Eigen::Matrix< fvar< T >, R1, C1 > multiply(const Eigen::Matrix< fvar< T >, R1, C1 > &m, const fvar< T > &c)
Definition: multiply.hpp:20
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
bool check_lower_triangular(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y)
Return true if the specified matrix is lower triangular.
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)
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)

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