Stan Math Library  2.8.0
reverse mode automatic differentiation
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Groups
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 
49 
50 namespace stan {
51  namespace math {
52 
53  // LKJ_Corr(L|eta) [ L Cholesky factor of correlation matrix
54  // eta > 0; eta == 1 <-> uniform]
55  template <bool propto,
56  typename T_covar, typename T_shape>
57  typename boost::math::tools::promote_args<T_covar, T_shape>::type
58  lkj_corr_cholesky_log(const Eigen::Matrix
59  <T_covar, Eigen::Dynamic, Eigen::Dynamic>& L,
60  const T_shape& eta) {
61  static const char* function("stan::math::lkj_corr_cholesky_log");
62 
63  using boost::math::tools::promote_args;
66  using stan::math::sum;
67 
68  typedef typename promote_args<T_covar, T_shape>::type lp_ret;
69  lp_ret lp(0.0);
70  check_positive(function, "Shape parameter", eta);
71  check_lower_triangular(function, "Random variable", L);
72 
73  const unsigned int K = L.rows();
74  if (K == 0)
75  return 0.0;
76 
78  lp += do_lkj_constant(eta, K);
80  const int Km1 = K - 1;
81  Eigen::Matrix<T_covar, Eigen::Dynamic, 1> log_diagonals =
82  L.diagonal().tail(Km1).array().log();
83  Eigen::Matrix<lp_ret, Eigen::Dynamic, 1> values(Km1);
84  for (int k = 0; k < Km1; k++)
85  values(k) = (Km1 - k - 1) * log_diagonals(k);
86  if ( (eta == 1.0) &&
88  lp += sum(values);
89  return(lp);
90  }
91  values += multiply(2.0 * eta - 2.0, log_diagonals);
92  lp += sum(values);
93  }
94 
95  return lp;
96  }
97 
98  template <typename T_covar, typename T_shape>
99  inline
100  typename boost::math::tools::promote_args<T_covar, T_shape>::type
101  lkj_corr_cholesky_log(const Eigen::Matrix
102  <T_covar, Eigen::Dynamic, Eigen::Dynamic>& L,
103  const T_shape& eta) {
104  return lkj_corr_cholesky_log<false>(L, eta);
105  }
106 
107  }
108 }
109 #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...
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.
bool check_lower_triangular(const char *function, const char *name, const Eigen::Matrix< T_y, Dynamic, Dynamic > &y)
Return true if the specified matrix is lower triangular.
T_shape do_lkj_constant(const T_shape &eta, const unsigned int &K)

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