Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
lkj_cov.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__LKJ_COV_HPP
3 
5 #include <stan/math/matrix.hpp>
8 #include <stan/meta/traits.hpp>
9 #include <stan/prob/traits.hpp>
10 
13 
14 namespace stan {
15 
16  namespace prob {
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 boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
23  lkj_cov_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
24  const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
25  const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
26  const T_shape& eta) {
27  static const char* function = "stan::prob::lkj_cov_log(%1%)";
28 
32  using boost::math::tools::promote_args;
33 
34  typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
35  check_size_match(function,
36  mu.rows(), "Rows of location parameter",
37  sigma.rows(), "columns of scale parameter",
38  &lp);
39  check_size_match(function,
40  y.rows(), "Rows of random variable",
41  y.cols(), "columns of random variable",
42  &lp);
43  check_size_match(function,
44  y.rows(), "Rows of random variable",
45  mu.rows(), "rows of location parameter",
46  &lp);
47  check_positive(function, eta, "Shape parameter", &lp);
48  check_finite(function, mu, "Location parameter", &lp);
49  check_finite(function, sigma, "Scale parameter", &lp);
50  // FIXME: build vectorized versions
51  for (int m = 0; m < y.rows(); ++m)
52  for (int n = 0; n < y.cols(); ++n)
53  check_finite(function, y(m,n), "Covariance matrix", &lp);
54 
55  const unsigned int K = y.rows();
56  const Eigen::Array<T_y,Eigen::Dynamic,1> sds
57  = y.diagonal().array().sqrt();
58  for (unsigned int k = 0; k < K; k++) {
59  lp += lognormal_log<propto>(sds(k), mu(k), sigma(k));
60  }
61  if (stan::is_constant<typename stan::scalar_type<T_shape> >::value
62  && eta == 1.0) {
63  // no need to rescale y into a correlation matrix
64  lp += lkj_corr_log<propto,T_y,T_shape>(y, eta);
65  return lp;
66  }
67  Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
68  D.diagonal() = sds.inverse();
69  lp += lkj_corr_log<propto,T_y,T_shape>(D * y * D, eta);
70  return lp;
71  }
72 
73  template <typename T_y, typename T_loc, typename T_scale, typename T_shape>
74  inline
75  typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
76  lkj_cov_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
77  const Eigen::Matrix<T_loc,Eigen::Dynamic,1>& mu,
78  const Eigen::Matrix<T_scale,Eigen::Dynamic,1>& sigma,
79  const T_shape& eta) {
80  return lkj_cov_log<false>(y,mu,sigma,eta);
81  }
82 
83  // LKJ_Cov(y|mu,sigma,eta) [ y covariance matrix (not correlation matrix)
84  // mu scalar, sigma > 0 scalar, eta > 0 ]
85  template <bool propto,
86  typename T_y, typename T_loc, typename T_scale, typename T_shape>
87  typename boost::math::tools::promote_args<T_y,T_loc,T_scale,T_shape>::type
88  lkj_cov_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
89  const T_loc& mu,
90  const T_scale& sigma,
91  const T_shape& eta) {
92  static const char* function = "stan::prob::lkj_cov_log(%1%)";
93 
96  using boost::math::tools::promote_args;
97 
98  typename promote_args<T_y,T_loc,T_scale,T_shape>::type lp(0.0);
99  check_positive(function, eta, "Shape parameter", &lp);
100  check_finite(function, mu, "Location parameter", &lp);
101  check_finite(function, sigma, "Scale parameter", &lp);
102 
103  const unsigned int K = y.rows();
104  const Eigen::Array<T_y,Eigen::Dynamic,1> sds
105  = y.diagonal().array().sqrt();
106  for (unsigned int k = 0; k < K; k++) {
107  lp += lognormal_log<propto>(sds(k), mu, sigma);
108  }
109  if (stan::is_constant<typename stan::scalar_type<T_shape> >::value
110  && eta == 1.0) {
111  // no need to rescale y into a correlation matrix
112  lp += lkj_corr_log<propto>(y,eta);
113  return lp;
114  }
115  Eigen::DiagonalMatrix<T_y,Eigen::Dynamic> D(K);
116  D.diagonal() = sds.inverse();
117  lp += lkj_corr_log<propto,T_y,T_shape>(D * y * D, eta);
118  return lp;
119  }
120 
121  template <typename T_y, typename T_loc, typename T_scale, typename T_shape>
122  inline
123  typename boost::math::tools::promote_args<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
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_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.hpp:23
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result)

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