Stan Math Library  2.8.0
reverse mode automatic differentiation
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Groups
lkj_cov_log.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_LKJ_COV_LOG_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_LKJ_COV_LOG_HPP
3 
8 
13 
14 namespace stan {
15 
16  namespace math {
17 
18  // LKJ_cov(y|mu, sigma, eta) [ y covariance matrix (not correlation matrix)
19  // mu vector, sigma > 0 vector, eta > 0 ]
20  template <bool propto,
21  typename T_y, typename T_loc, typename T_scale, typename T_shape>
22  typename
23  boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
24  lkj_cov_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
25  const Eigen::Matrix<T_loc, Eigen::Dynamic, 1>& mu,
26  const Eigen::Matrix<T_scale, Eigen::Dynamic, 1>& sigma,
27  const T_shape& eta) {
28  static const char* function("stan::math::lkj_cov_log");
29 
34  using boost::math::tools::promote_args;
35 
36  typename promote_args<T_y, T_loc, T_scale, T_shape>::type lp(0.0);
37  check_size_match(function,
38  "Rows of location parameter", mu.rows(),
39  "columns of scale parameter", sigma.rows());
40  check_square(function, "random variable", y);
41  check_size_match(function,
42  "Rows of random variable", y.rows(),
43  "rows of location parameter", mu.rows());
44  check_positive(function, "Shape parameter", eta);
45  check_finite(function, "Location parameter", mu);
46  check_finite(function, "Scale parameter", sigma);
47  // FIXME: build vectorized versions
48  for (int m = 0; m < y.rows(); ++m)
49  for (int n = 0; n < y.cols(); ++n)
50  check_finite(function, "Covariance matrix", y(m, n));
51 
52  const unsigned int K = y.rows();
53  const Eigen::Array<T_y, Eigen::Dynamic, 1> sds
54  = y.diagonal().array().sqrt();
55  for (unsigned int k = 0; k < K; k++) {
56  lp += lognormal_log<propto>(sds(k), mu(k), sigma(k));
57  }
58  if (stan::is_constant<typename stan::scalar_type<T_shape> >::value
59  && eta == 1.0) {
60  // no need to rescale y into a correlation matrix
61  lp += lkj_corr_log<propto, T_y, T_shape>(y, eta);
62  return lp;
63  }
64  Eigen::DiagonalMatrix<T_y, Eigen::Dynamic> D(K);
65  D.diagonal() = sds.inverse();
66  lp += lkj_corr_log<propto, T_y, T_shape>(D * y * D, eta);
67  return lp;
68  }
69 
70  template <typename T_y, typename T_loc, typename T_scale, typename T_shape>
71  inline
72  typename
73  boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
74  lkj_cov_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
75  const Eigen::Matrix<T_loc, Eigen::Dynamic, 1>& mu,
76  const Eigen::Matrix<T_scale, Eigen::Dynamic, 1>& sigma,
77  const T_shape& eta) {
78  return lkj_cov_log<false>(y, mu, sigma, eta);
79  }
80 
81  // LKJ_Cov(y|mu, sigma, eta) [ y covariance matrix (not correlation matrix)
82  // mu scalar, sigma > 0 scalar, eta > 0 ]
83  template <bool propto,
84  typename T_y, typename T_loc, typename T_scale, typename T_shape>
85  typename
86  boost::math::tools::promote_args<T_y, T_loc, T_scale, T_shape>::type
87  lkj_cov_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
88  const T_loc& mu,
89  const T_scale& sigma,
90  const T_shape& eta) {
91  static const char* function("stan::math::lkj_cov_log");
92 
95  using boost::math::tools::promote_args;
96 
97  typename promote_args<T_y, T_loc, T_scale, T_shape>::type lp(0.0);
98  check_positive(function, "Shape parameter", eta);
99  check_finite(function, "Location parameter", mu);
100  check_finite(function, "Scale parameter", sigma);
101 
102  const unsigned int K = y.rows();
103  const Eigen::Array<T_y, Eigen::Dynamic, 1> sds
104  = y.diagonal().array().sqrt();
105  for (unsigned int k = 0; k < K; k++) {
106  lp += lognormal_log<propto>(sds(k), mu, sigma);
107  }
108  if (stan::is_constant<typename stan::scalar_type<T_shape> >::value
109  && eta == 1.0) {
110  // no need to rescale y into a correlation matrix
111  lp += lkj_corr_log<propto>(y, eta);
112  return lp;
113  }
114  Eigen::DiagonalMatrix<T_y, Eigen::Dynamic> D(K);
115  D.diagonal() = sds.inverse();
116  lp += lkj_corr_log<propto, T_y, T_shape>(D * y * D, eta);
117  return lp;
118  }
119 
120  template <typename T_y, typename T_loc, typename T_scale, typename T_shape>
121  inline
122  typename boost::math::tools::promote_args
123  <T_y, T_loc, T_scale, T_shape>::type
124  lkj_cov_log(const Eigen::Matrix<T_y, Eigen::Dynamic, Eigen::Dynamic>& y,
125  const T_loc& mu,
126  const T_scale& sigma,
127  const T_shape& eta) {
128  return lkj_cov_log<false>(y, mu, sigma, eta);
129  }
130 
131 
132  }
133 }
134 #endif
boost::math::tools::promote_args< T_y, T_loc, T_scale, T_shape >::type lkj_cov_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const Eigen::Matrix< T_loc, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< T_scale, Eigen::Dynamic, 1 > &sigma, const T_shape &eta)
Definition: lkj_cov_log.hpp:24
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
bool check_positive(const char *function, const char *name, const T_y &y)
Return true if y is positive.
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.
bool check_finite(const char *function, const char *name, const T_y &y)
Return true if y is finite.
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.

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