Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gaussian_dlm_obs.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__GAUSSIAN_DLM_OBS_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__GAUSSIAN_DLM_OBS_HPP
3 
4 #include <boost/random/normal_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
6 
10 #include <stan/prob/constants.hpp>
11 #include <stan/prob/traits.hpp>
12 #include <stan/agrad/rev.hpp>
13 #include <stan/meta/traits.hpp>
16 
17 #include <stan/math/matrix/add.hpp>
20 #include <stan/math/matrix/log.hpp>
28 
29 /*
30  TODO: time-varying system matrices
31  TODO: use sequential processing even for non-diagonal obs
32  covariance.
33  TODO: add constant terms in observation.
34  */
35 
36 namespace stan {
37  namespace prob {
71  template <bool propto,
72  typename T_y,
73  typename T_F, typename T_G,
74  typename T_V, typename T_W,
75  typename T_m0, typename T_C0
76  >
77  typename return_type<
78  T_y,
79  typename return_type<T_F,T_G,T_V,T_W,T_m0,T_C0>::type >::type
80  gaussian_dlm_obs_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
81  const Eigen::Matrix<T_F,Eigen::Dynamic,Eigen::Dynamic>& F,
82  const Eigen::Matrix<T_G,Eigen::Dynamic,Eigen::Dynamic>& G,
83  const Eigen::Matrix<T_V,Eigen::Dynamic,Eigen::Dynamic>& V,
84  const Eigen::Matrix<T_W,Eigen::Dynamic,Eigen::Dynamic>& W,
85  const Eigen::Matrix<T_m0,Eigen::Dynamic,1>& m0,
86  const Eigen::Matrix<T_C0,Eigen::Dynamic,Eigen::Dynamic>& C0) {
87  static const char* function = "stan::prob::gaussian_dlm_obs_log(%1%)";
88  typedef typename return_type<
89  T_y,
91  T_lp lp(0.0);
92 
93  using stan::math::add;
101  using stan::math::multiply;
103  using stan::math::subtract;
105  using stan::math::transpose;
106 
107  int r = y.rows(); // number of variables
108  int T = y.cols(); // number of observations
109  int n = G.rows(); // number of states
110 
111  // check y
112  check_finite(function, y, "y", &lp);
113  check_not_nan(function, y, "y", &lp);
114  // check F
115  check_size_match(function,
116  F.cols(), "columns of F",
117  y.rows(), "rows of y",
118  &lp);
119  check_size_match(function,
120  F.rows(), "rows of F",
121  G.rows(), "rows of G",
122  &lp);
123  check_finite(function, F, "F", &lp);
124  check_not_nan(function, F, "F", &lp);
125  // check G
126  check_size_match(function,
127  G.rows(), "rows of G",
128  G.cols(), "columns of G",
129  &lp);
130  check_finite(function, G, "G", &lp);
131  check_not_nan(function, G, "G", &lp);
132  // check V
133  check_size_match(function,
134  V.rows(), "rows of V",
135  y.rows(), "rows of y",
136  &lp);
137  // TODO: incorporate support for infinite V
138  check_finite(function, V, "V", &lp);
139  check_not_nan(function, V, "V", &lp);
140  check_spsd_matrix(function, V, "V", &lp);
141  // check W
142  check_size_match(function,
143  W.rows(), "rows of W",
144  G.rows(), "rows of G",
145  &lp);
146  // TODO: incorporate support for infinite W
147  check_finite(function, W, "W", &lp);
148  check_not_nan(function, W, "W", &lp);
149  check_spsd_matrix(function, W, "W", &lp);
150  // check m0
151  check_size_match(function,
152  m0.size(), "size of m0",
153  G.rows(), "rows of G",
154  &lp);
155  check_finite(function, m0, "m0", &lp);
156  check_not_nan(function, m0, "m0", &lp);
157  // check C0
158  check_size_match(function,
159  C0.rows(), "rows of C0",
160  G.rows(), "rows of G",
161  &lp);
162  check_cov_matrix(function, C0, "C0", &lp);
163  check_finite(function, C0, "C0", &lp);
164  check_not_nan(function, C0, "C0", &lp);
165 
166  if (y.cols() == 0 || y.rows() == 0)
167  return lp;
168 
170  lp -= 0.5 * LOG_TWO_PI * r * T;
171  }
172 
174  Eigen::Matrix<T_lp,Eigen::Dynamic, 1> m(n);
175  Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> C(n, n);
176 
177  // TODO: how to recast matrices
178  for (int i = 0; i < m0.size(); i ++ ) {
179  m(i) = m0(i);
180  }
181  for (int i = 0; i < C0.rows(); i ++ ) {
182  for (int j = 0; j < C0.cols(); j ++ ) {
183  C(i, j) = C0(i, j);
184  }
185  }
186 
187  Eigen::Matrix<typename return_type<T_y>::type,
188  Eigen::Dynamic, 1> yi(r);
189  Eigen::Matrix<T_lp,Eigen::Dynamic, 1> a(n);
190  Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> R(n, n);
191  Eigen::Matrix<T_lp,Eigen::Dynamic, 1> f(r);
192  Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> Q(r, r);
193  Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> Q_inv(r, r);
194  Eigen::Matrix<T_lp,Eigen::Dynamic, 1> e(r);
195  Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> A(n, r);
196 
197  for (int i = 0; i < y.cols(); i++) {
198  yi = y.col(i);
199  // // Predict state
200  // a_t = G_t m_{t-1}
201  a = multiply(G, m);
202  // R_t = G_t C_{t-1} G_t' + W_t
203  R = add(quad_form_sym(C, transpose(G)), W);
204  // // predict observation
205  // f_t = F_t' a_t
206  f = multiply(transpose(F), a);
207  // Q_t = F'_t R_t F_t + V_t
208  Q = add(quad_form_sym(R, F), V);
209  Q_inv = inverse_spd(Q);
210  // // filtered state
211  // e_t = y_t - f_t
212  e = subtract(yi, f);
213  // A_t = R_t F_t Q^{-1}_t
214  A = multiply(multiply(R, F), Q_inv);
215  // m_t = a_t + A_t e_t
216  m = add(a, multiply(A, e));
217  // C = R_t - A_t Q_t A_t'
218  C = subtract(R, quad_form_sym(Q, transpose(A)));
219  lp -= 0.5 * (log_determinant_spd(Q) + trace_quad_form(Q_inv, e));
220  }
221  }
222  return lp;
223  }
224 
225  template <typename T_y,
226  typename T_F, typename T_G,
227  typename T_V, typename T_W,
228  typename T_m0, typename T_C0
229  >
230  inline
231  typename return_type<
232  T_y,
234  gaussian_dlm_obs_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
235  const Eigen::Matrix<T_F,Eigen::Dynamic,Eigen::Dynamic>& F,
236  const Eigen::Matrix<T_G,Eigen::Dynamic,Eigen::Dynamic>& G,
237  const Eigen::Matrix<T_V,Eigen::Dynamic,Eigen::Dynamic>& V,
238  const Eigen::Matrix<T_W,Eigen::Dynamic,Eigen::Dynamic>& W,
239  const Eigen::Matrix<T_m0,Eigen::Dynamic,1>& m0,
240  const Eigen::Matrix<T_C0,Eigen::Dynamic,Eigen::Dynamic>& C0) {
241  return gaussian_dlm_obs_log<false>(y, F, G, V, W, m0, C0);
242  }
243 
279  template <bool propto,
280  typename T_y,
281  typename T_F, typename T_G,
282  typename T_V, typename T_W,
283  typename T_m0, typename T_C0
284  >
285  typename return_type<
286  T_y,
288  gaussian_dlm_obs_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
289  const Eigen::Matrix<T_F,Eigen::Dynamic,Eigen::Dynamic>& F,
290  const Eigen::Matrix<T_G,Eigen::Dynamic,Eigen::Dynamic>& G,
291  const Eigen::Matrix<T_V,Eigen::Dynamic,1>& V,
292  const Eigen::Matrix<T_W,Eigen::Dynamic,Eigen::Dynamic>& W,
293  const Eigen::Matrix<T_m0,Eigen::Dynamic,1>& m0,
294  const Eigen::Matrix<T_C0,Eigen::Dynamic,Eigen::Dynamic>& C0) {
295  static const char* function = "stan::prob::gaussian_dlm_obs_log(%1%)";
296  typedef typename return_type<
297  T_y,
298  typename return_type<T_F,T_G,T_V,T_W,T_m0,T_C0>::type >::type T_lp;
299  T_lp lp(0.0);
300 
301  using stan::math::add;
309  using stan::math::multiply;
311  using stan::math::subtract;
314  using stan::math::transpose;
315 
316  int r = y.rows(); // number of variables
317  int T = y.cols(); // number of observations
318  int n = G.rows(); // number of states
319 
320  // check y
321  check_finite(function, y, "y", &lp);
322  check_not_nan(function, y, "y", &lp);
323  // check F
324  check_size_match(function,
325  F.cols(), "columns of F",
326  y.rows(), "rows of y",
327  &lp);
328  check_size_match(function,
329  F.rows(), "rows of F",
330  G.rows(), "rows of G",
331  &lp);
332  check_finite(function, F, "F", &lp);
333  check_not_nan(function, F, "F", &lp);
334  // check G
335  check_size_match(function,
336  G.rows(), "rows of G",
337  G.cols(), "columns of G",
338  &lp);
339  check_finite(function, G, "G", &lp);
340  check_not_nan(function, G, "G", &lp);
341  // check V
342  check_nonnegative(function, V, "V", &lp);
343  check_size_match(function,
344  V.size(), "size of V",
345  y.rows(), "rows of y",
346  &lp);
347  // TODO: support infinite V
348  check_finite(function, V, "V", &lp);
349  check_not_nan(function, V, "V", &lp);
350  // check W
351  check_spsd_matrix(function, W, "W", &lp);
352  check_size_match(function,
353  W.rows(), "rows of W",
354  G.rows(), "rows of G",
355  &lp);
356  // TODO: support infinite W
357  check_finite(function, W, "W", &lp);
358  check_not_nan(function, W, "W", &lp);
359  // check m0
360  check_size_match(function,
361  m0.size(), "size of m0",
362  G.rows(), "rows of G",
363  &lp);
364  check_finite(function, m0, "m0", &lp);
365  check_not_nan(function, m0, "m0", &lp);
366  // check C0
367  check_cov_matrix(function, C0, "C0", &lp);
368  check_size_match(function,
369  C0.rows(), "rows of C0",
370  G.rows(), "rows of G",
371  &lp);
372  check_finite(function, C0, "C0", &lp);
373  check_not_nan(function, C0, "C0", &lp);
374 
375  if (y.cols() == 0 || y.rows() == 0)
376  return lp;
377 
379  lp += 0.5 * NEG_LOG_TWO_PI * r * T;
380  }
381 
383  T_lp f;
384  T_lp Q;
385  T_lp Q_inv;
386  T_lp e;
387  Eigen::Matrix<T_lp, Eigen::Dynamic, 1> A(n);
388  Eigen::Matrix<T_lp, Eigen::Dynamic, 1> Fj(n);
389  Eigen::Matrix<T_lp, Eigen::Dynamic, 1> m(n);
390  Eigen::Matrix<T_lp, Eigen::Dynamic, Eigen::Dynamic> C(n, n);
391 
392  // TODO: how to recast matrices
393  for (int i = 0; i < m0.size(); i ++ ) {
394  m(i) = m0(i);
395  }
396  for (int i = 0; i < C0.rows(); i ++ ) {
397  for (int j = 0; j < C0.cols(); j ++ ) {
398  C(i, j) = C0(i, j);
399  }
400  }
401 
402  for (int i = 0; i < y.cols(); i++) {
403  // Predict state
404  // reuse m and C instead of using a and R
405  m = multiply(G, m);
406  C = add(quad_form_sym(C, transpose(G)), W);
407  for (int j = 0; j < y.rows(); ++j) {
408  // predict observation
409  T_lp yij(y(j, i));
410  // dim Fj = (n, 1)
411  for (int k = 0; k < F.rows(); ++k) {
412  Fj(k) = F(k, j);
413  }
414  // // f_{t,i} = F_{t,i}' m_{t,i-1}
415  f = dot_product(Fj, m);
416  Q = trace_quad_form(C, Fj) + V(j);
417  Q_inv = 1.0 / Q;
418  // // filtered observation
419  // // e_{t,i} = y_{t,i} - f_{t,i}
420  e = yij - f;
421  // // A_{t,i} = C_{t,i-1} F_{t,i} Q_{t,i}^{-1}
422  A = multiply(multiply(C, Fj), Q_inv);
423  // // m_{t,i} = m_{t,i-1} + A_{t,i} e_{t,i}
424  m += multiply(A, e);
425  // // c_{t,i} = C_{t,i-1} - Q_{t,i} A_{t,i} A_{t,i}'
426  // // // tcrossprod throws an error (ambiguous)
427  // C = subtract(C, multiply(Q, tcrossprod(A)));
428  C -= multiply(Q, multiply(A, transpose(A)));
429  C = 0.5 * add(C, transpose(C));
430  lp -= 0.5 * (log(Q) + pow(e, 2) * Q_inv);
431  }
432  }
433  }
434  return lp;
435  }
436 
437  template <typename T_y,
438  typename T_F, typename T_G,
439  typename T_V, typename T_W,
440  typename T_m0, typename T_C0
441  >
442  inline
443  typename return_type<
444  T_y,
446  >::type
447  gaussian_dlm_obs_log(const Eigen::Matrix<T_y,Eigen::Dynamic,Eigen::Dynamic>& y,
448  const Eigen::Matrix<T_F,Eigen::Dynamic,Eigen::Dynamic>& F,
449  const Eigen::Matrix<T_G,Eigen::Dynamic,Eigen::Dynamic>& G,
450  const Eigen::Matrix<T_V,Eigen::Dynamic,1>& V,
451  const Eigen::Matrix<T_W,Eigen::Dynamic,Eigen::Dynamic>& W,
452  const Eigen::Matrix<T_m0,Eigen::Dynamic,1>& m0,
453  const Eigen::Matrix<T_C0,Eigen::Dynamic,Eigen::Dynamic>& C0) {
454  return gaussian_dlm_obs_log<false>(y, F, G, V, W, m0, C0);
455  }
456  }
457 
458 }
459 
460 #endif
bool check_cov_matrix(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is a valid covariance matrix.
bool check_spsd_matrix(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const char *name, T_result *result)
Return true if the specified matrix is a square, symmetric, and positive semi-definite.
Eigen::Matrix< typename boost::math::tools::promote_args< T1, T2 >::type, R, C > add(const Eigen::Matrix< T1, R, C > &m1, const Eigen::Matrix< T2, R, C > &m2)
Return the sum of the specified matrices.
Definition: add.hpp:27
return_type< T_y, typename return_type< T_F, T_G, T_V, T_W, T_m0, T_C0 >::type >::type gaussian_dlm_obs_log(const Eigen::Matrix< T_y, Eigen::Dynamic, Eigen::Dynamic > &y, const Eigen::Matrix< T_F, Eigen::Dynamic, Eigen::Dynamic > &F, const Eigen::Matrix< T_G, Eigen::Dynamic, Eigen::Dynamic > &G, const Eigen::Matrix< T_V, Eigen::Dynamic, Eigen::Dynamic > &V, const Eigen::Matrix< T_W, Eigen::Dynamic, Eigen::Dynamic > &W, const Eigen::Matrix< T_m0, Eigen::Dynamic, 1 > &m0, const Eigen::Matrix< T_C0, Eigen::Dynamic, Eigen::Dynamic > &C0)
The log of a Gaussian dynamic linear model (GDLM).
fvar< T > pow(const fvar< T > &x1, const fvar< T > &x2)
Definition: pow.hpp:17
bool check_size_match(const char *function, T_size1 i, const char *name_i, T_size2 j, const char *name_j, T_result *result)
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > inverse_spd(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m)
Returns the inverse of the specified symmetric, pos/neg-definite matrix.
Definition: inverse_spd.hpp:19
Eigen::Matrix< T, CB, CB > quad_form_sym(const Eigen::Matrix< T, RA, CA > &A, const Eigen::Matrix< T, RB, CB > &B)
Definition: quad_form.hpp:47
Metaprogram to calculate the base scalar return type resulting from promoting all the scalar types of...
Definition: traits.hpp:398
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
var log_determinant_spd(const Eigen::Matrix< var, R, C > &m)
stan::agrad::fvar< T > trace_quad_form(const Eigen::Matrix< stan::agrad::fvar< T >, RA, CA > &A, const Eigen::Matrix< stan::agrad::fvar< T >, RB, CB > &B)
Eigen::Matrix< fvar< T >, R1, C1 > multiply(const Eigen::Matrix< fvar< T >, R1, C1 > &m, const fvar< T > &c)
Definition: multiply.hpp:21
matrix_d tcrossprod(const matrix_d &M)
Returns the result of post-multiplying a matrix by its own transpose.
Definition: tcrossprod.hpp:17
boost::math::tools::promote_args< typename scalar_type< T1 >::type, typename scalar_type< T2 >::type, typename scalar_type< T3 >::type, typename scalar_type< T4 >::type, typename scalar_type< T5 >::type, typename scalar_type< T6 >::type >::type type
Definition: traits.hpp:406
Eigen::Matrix< T, C, R > transpose(const Eigen::Matrix< T, R, C > &m)
Definition: transpose.hpp:12
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
double trace_quad_form(const Eigen::Matrix< double, RA, CA > &A, const Eigen::Matrix< double, RB, CB > &B)
Compute trace(B^T A B).
fvar< T > dot_product(const Eigen::Matrix< fvar< T >, R1, C1 > &v1, const Eigen::Matrix< fvar< T >, R2, C2 > &v2)
Definition: dot_product.hpp:20
bool check_nonnegative(const char *function, const T_y &y, const char *name, T_result *result)
bool check_not_nan(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is nan.
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
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
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:15
boost::enable_if_c< boost::is_same< TA, var >::value||boost::is_same< TB, var >::value, Eigen::Matrix< var, CB, CB > >::type quad_form_sym(const Eigen::Matrix< TA, RA, CA > &A, const Eigen::Matrix< TB, RB, CB > &B)
Definition: quad_form.hpp:162
T log_determinant_spd(const Eigen::Matrix< T, R, C > &m)
Returns the log absolute determinant of the specified square matrix.
double dot_product(const Eigen::Matrix< double, R1, C1 > &v1, const Eigen::Matrix< double, R2, C2 > &v2)
Returns the dot product of the specified vectors.
Definition: dot_product.hpp:22

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