1 #ifndef STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__GAUSSIAN_DLM_OBS_HPP
2 #define STAN__PROB__DISTRIBUTIONS__MULTIVARIATE__CONTINUOUS__GAUSSIAN_DLM_OBS_HPP
4 #include <boost/random/normal_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
71 template <
bool propto,
73 typename T_F,
typename T_G,
74 typename T_V,
typename T_W,
75 typename T_m0,
typename T_C0
79 typename return_type<T_F,T_G,T_V,T_W,T_m0,T_C0>::type >::type
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%)";
116 F.cols(),
"columns of F",
117 y.rows(),
"rows of y",
120 F.rows(),
"rows of F",
121 G.rows(),
"rows of G",
127 G.rows(),
"rows of G",
128 G.cols(),
"columns of G",
134 V.rows(),
"rows of V",
135 y.rows(),
"rows of y",
143 W.rows(),
"rows of W",
144 G.rows(),
"rows of G",
152 m0.size(),
"size of m0",
153 G.rows(),
"rows of G",
159 C0.rows(),
"rows of C0",
160 G.rows(),
"rows of G",
166 if (y.cols() == 0 || y.rows() == 0)
170 lp -= 0.5 * LOG_TWO_PI * r * T;
174 Eigen::Matrix<T_lp,Eigen::Dynamic, 1> m(n);
175 Eigen::Matrix<T_lp,Eigen::Dynamic, Eigen::Dynamic> C(n, n);
178 for (
int i = 0; i < m0.size(); i ++ ) {
181 for (
int i = 0; i < C0.rows(); i ++ ) {
182 for (
int j = 0; j < C0.cols(); j ++ ) {
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);
197 for (
int i = 0; i < y.cols(); i++) {
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
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);
279 template <
bool propto,
281 typename T_F,
typename T_G,
282 typename T_V,
typename T_W,
283 typename T_m0,
typename T_C0
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%)";
325 F.cols(),
"columns of F",
326 y.rows(),
"rows of y",
329 F.rows(),
"rows of F",
330 G.rows(),
"rows of G",
336 G.rows(),
"rows of G",
337 G.cols(),
"columns of G",
344 V.size(),
"size of V",
345 y.rows(),
"rows of y",
353 W.rows(),
"rows of W",
354 G.rows(),
"rows of G",
361 m0.size(),
"size of m0",
362 G.rows(),
"rows of G",
369 C0.rows(),
"rows of C0",
370 G.rows(),
"rows of G",
375 if (y.cols() == 0 || y.rows() == 0)
379 lp += 0.5 * NEG_LOG_TWO_PI * r * T;
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);
393 for (
int i = 0; i < m0.size(); i ++ ) {
396 for (
int i = 0; i < C0.rows(); i ++ ) {
397 for (
int j = 0; j < C0.cols(); j ++ ) {
402 for (
int i = 0; i < y.cols(); i++) {
407 for (
int j = 0; j < y.rows(); ++j) {
411 for (
int k = 0; k < F.rows(); ++k) {
430 lp -= 0.5 * (
log(Q) +
pow(e, 2) * Q_inv);
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
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);
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.
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)
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.
Eigen::Matrix< T, CB, CB > quad_form_sym(const Eigen::Matrix< T, RA, CA > &A, const Eigen::Matrix< T, RB, CB > &B)
Metaprogram to calculate the base scalar return type resulting from promoting all the scalar types of...
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)
matrix_d tcrossprod(const matrix_d &M)
Returns the result of post-multiplying a matrix by its own transpose.
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
Eigen::Matrix< T, C, R > transpose(const Eigen::Matrix< T, R, C > &m)
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
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...
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)
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.
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.
fvar< T > log(const fvar< T > &x)
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)
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.