Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
multi_normal_cholesky.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__MULTI_NORMAL_CHOLESKY_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__MULTI_NORMAL_CHOLESKY_HPP
3 
4 #include <boost/random/normal_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
6 
7 #include <boost/random/normal_distribution.hpp>
8 #include <boost/random/variate_generator.hpp>
9 
10 #include <stan/agrad/rev.hpp>
18 #include <stan/math/matrix/log.hpp>
24 #include <stan/math/matrix/sum.hpp>
25 #include <stan/meta/traits.hpp>
26 #include <stan/prob/constants.hpp>
27 #include <stan/prob/traits.hpp>
28 
29 namespace stan {
30 
31  namespace prob {
32 
50  template <bool propto,
51  typename T_y, typename T_loc, typename T_covar>
52  typename boost::math::tools::promote_args<typename scalar_type<T_y>::type, typename scalar_type<T_loc>::type, T_covar>::type
54  const T_loc& mu,
55  const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L) {
56  static const char* function = "stan::prob::multi_normal_cholesky_log(%1%)";
57  typedef typename boost::math::tools::promote_args<typename scalar_type<T_y>::type, typename scalar_type<T_loc>::type, T_covar>::type lp_type;
58  lp_type lp(0.0);
59 
64  using stan::math::sum;
65 
69 
70  VectorViewMvt<const T_y> y_vec(y);
71  VectorViewMvt<const T_loc> mu_vec(mu);
72  //size of std::vector of Eigen vectors
73  size_t size_vec = max_size_mvt(y, mu);
74 
75  //Check if every vector of the array has the same size
76  int size_y = y_vec[0].size();
77  int size_mu = mu_vec[0].size();
78  if (size_vec > 1) {
79  int size_y_old = size_y;
80  int size_y_new;
81  for (size_t i = 1, size_ = length_mvt(y); i < size_; i++) {
82  int size_y_new = y_vec[i].size();
83  check_size_match(function,
84  size_y_new, "Size of one of the vectors of the random variable",
85  size_y_old, "Size of another vector of the random variable",
86  &lp);
87  size_y_old = size_y_new;
88  }
89  int size_mu_old = size_mu;
90  int size_mu_new;
91  for (size_t i = 1, size_ = length_mvt(mu); i < size_; i++) {
92  int size_mu_new = mu_vec[i].size();
93  check_size_match(function,
94  size_mu_new, "Size of one of the vectors of the location variable",
95  size_mu_old, "Size of another vector of the location variable",
96  &lp);
97  size_mu_old = size_mu_new;
98  }
99  (void) size_y_old;
100  (void) size_y_new;
101  (void) size_mu_old;
102  (void) size_mu_new;
103  }
104 
105  check_size_match(function,
106  size_y, "Size of random variable",
107  size_mu, "size of location parameter",
108  &lp);
109  check_size_match(function,
110  size_y, "Size of random variable",
111  L.rows(), "rows of covariance parameter",
112  &lp);
113  check_size_match(function,
114  size_y, "Size of random variable",
115  L.cols(), "columns of covariance parameter",
116  &lp);
117 
118  for (size_t i = 0; i < size_vec; i++) {
119  check_finite(function, mu_vec[i], "Location parameter", &lp);
120  check_not_nan(function, y_vec[i], "Random variable", &lp);
121  }
122 
123  if (size_y == 0)
124  return lp;
125 
127  lp += NEG_LOG_SQRT_TWO_PI * size_y * size_vec;
128 
130  lp -= L.diagonal().array().log().sum() * size_vec;
131 
133  lp_type sum_lp_vec(0.0);
134  for (size_t i = 0; i < size_vec; i++) {
135  Eigen::Matrix<typename
136  boost::math::tools::promote_args<typename scalar_type<T_y>::type, typename scalar_type<T_loc>::type>::type,
137  Eigen::Dynamic, 1> y_minus_mu(size_y);
138  for (int j = 0; j < size_y; j++)
139  y_minus_mu(j) = y_vec[i](j)-mu_vec[i](j);
140  Eigen::Matrix<typename
141  boost::math::tools::promote_args<T_covar,typename scalar_type<T_loc>::type,typename scalar_type<T_y>::type>::type,
142  Eigen::Dynamic, 1>
143  half(mdivide_left_tri_low(L,y_minus_mu));
144  // FIXME: this code does not compile. revert after fixing subtract()
145  // Eigen::Matrix<typename
146  // boost::math::tools::promote_args<T_covar,typename scalar_type<T_loc>::type,typename scalar_type<T_y>::type>::type>::type,
147  // Eigen::Dynamic, 1>
148  // half(mdivide_left_tri_low(L,subtract(y,mu)));
149  sum_lp_vec += dot_self(half);
150  }
151  lp -= 0.5*sum_lp_vec;
152  }
153  return lp;
154  }
155 
156  template <typename T_y, typename T_loc, typename T_covar>
157  inline
158  typename boost::math::tools::promote_args<typename scalar_type<T_y>::type, typename scalar_type<T_loc>::type, T_covar>::type
160  const T_loc& mu,
161  const Eigen::Matrix<T_covar,Eigen::Dynamic,Eigen::Dynamic>& L) {
162  return multi_normal_cholesky_log<false>(y,mu,L);
163  }
164 
165  template <class RNG>
166  inline Eigen::VectorXd
167  multi_normal_cholesky_rng(const Eigen::Matrix<double,Eigen::Dynamic,1>& mu,
168  const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& S,
169  RNG& rng) {
170  using boost::variate_generator;
171  using boost::normal_distribution;
172 
173  static const char* function = "stan::prob::multi_normal_cholesky_rng(%1%)";
174 
176 
177  check_finite(function, mu, "Location parameter", (double*)0);
178 
179  variate_generator<RNG&, normal_distribution<> >
180  std_normal_rng(rng, normal_distribution<>(0,1));
181 
182  Eigen::VectorXd z(S.cols());
183  for(int i = 0; i < S.cols(); i++)
184  z(i) = std_normal_rng();
185 
186  return mu + S * z;
187  }
188  }
189 }
190 
191 #endif
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R1, C2 > mdivide_left_tri_low(const Eigen::Matrix< T1, R1, C1 > &A, const Eigen::Matrix< T2, R2, C2 > &b)
fvar< T > dot_self(const Eigen::Matrix< fvar< T >, R, C > &v)
Definition: dot_self.hpp:16
double dot_self(const std::vector< double > &x)
Definition: dot_self.hpp:11
Eigen::VectorXd multi_normal_cholesky_rng(const Eigen::Matrix< double, Eigen::Dynamic, 1 > &mu, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &S, RNG &rng)
size_t max_size_mvt(const T1 &x1, const T2 &x2)
Definition: traits.hpp:575
boost::math::tools::promote_args< typename scalar_type< T_y >::type, typename scalar_type< T_loc >::type, T_covar >::type multi_normal_cholesky_log(const T_y &y, const T_loc &mu, const Eigen::Matrix< T_covar, Eigen::Dynamic, Eigen::Dynamic > &L)
The log of the multivariate normal density for the given y, mu, and a Cholesky factor L of the varian...
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.
scalar_type_helper< is_vector< T >::value, T >::type type
Definition: traits.hpp:139
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Definition: traits.hpp:35
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R, C > subtract(const Eigen::Matrix< T1, R, C > &m1, const Eigen::Matrix< T2, R, C > &m2)
Return the result of subtracting the second specified matrix from the first specified matrix...
Definition: subtract.hpp:27
size_t length_mvt(const T &)
Definition: traits.hpp:561
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.
size_t size_
Definition: dot_self.hpp:18
boost::enable_if_c< boost::is_arithmetic< T >::value, Eigen::Matrix< double, R, C > >::type multiply(const Eigen::Matrix< double, R, C > &m, T c)
Return specified matrix multiplied by specified scalar.
Definition: multiply.hpp:25
Eigen::Matrix< fvar< T >, R1, C1 > mdivide_left_tri_low(const Eigen::Matrix< fvar< T >, R1, C1 > &A, const Eigen::Matrix< fvar< T >, R2, C2 > &b)

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