A base class for the bootstrap particle filter with covariates.
More...
#include <bootstrap_filter_with_covariates.h>
|
using | ssv = Eigen::Matrix< float_t, dimx, 1 > |
|
using | osv = Eigen::Matrix< float_t, dimy, 1 > |
|
using | cvsv = Eigen::Matrix< float_t, dimcov, 1 > |
|
using | Mat = Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > |
|
using | arrayStates = std::array< ssv, nparts > |
|
using | arrayfloat_t = std::array< float_t, nparts > |
|
using | Funcs = std::vector< std::function< const Mat(const ssv &, const cvsv &)> > |
|
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t, typename float_t>
class BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >
A base class for the bootstrap particle filter with covariates.
- Author
- taylor
◆ arrayfloat_t
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
type alias for array of float_ts
◆ arrayStates
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
type alias for linear algebra stuff
◆ cvsv
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
using BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::cvsv = Eigen::Matrix<float_t,dimcov,1> |
covariate size vector" type alias for linear algebra stuff
◆ Funcs
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
using BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::Funcs = std::vector<std::function<const Mat(const ssv&, const cvsv&)> > |
◆ Mat
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
using BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic> |
type alias for dynamically sized matrix
◆ osv
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
using BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::osv = Eigen::Matrix<float_t, dimy, 1> |
"obs size vector" type alias for linear algebra stuff
◆ ssv
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
using BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::ssv = Eigen::Matrix<float_t, dimx, 1> |
"state size vector" type alias for linear algebra stuff
◆ BSFilterWC()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::BSFilterWC |
( |
const unsigned int & |
rs = 1 | ) |
|
The constructor.
- Parameters
-
rs | the resampling schedule (e.g. every rs time point) |
◆ filter()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
void BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::filter |
( |
const osv & |
ydata, |
|
|
const cvsv & |
covdata, |
|
|
const Funcs & |
fs = Funcs() |
|
) |
| |
updates filtering distribution on a new datapoint. Optionally stores expectations of functionals.
- Parameters
-
data | the most recent data point |
covData | covariate data |
fs | a vector of functions if you want to calculate expectations. |
◆ fSamp()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
virtual ssv BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::fSamp |
( |
const ssv & |
xtm1, |
|
|
const cvsv & |
zt |
|
) |
| |
|
pure virtual |
Sample from the state transition distribution.
- Parameters
-
xtm1 | is a const Vec& describing the time t-1 state |
zt | is a const Vec& describing the time t covariate |
- Returns
- the sample as a Vec
◆ getExpectations()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
auto BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::getExpectations |
( |
| ) |
const -> std::vector<Mat> |
return all stored expectations (taken with respect to $p(x_t|y_{1:t})$
- Returns
- return a std::vector<Mat> of expectations. How many depends on how many callbacks you gave to
◆ getLogCondLike()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
float_t BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::getLogCondLike |
( |
| ) |
const |
Returns the most recent (log-) conditiona likelihood.
- Returns
- log p(y_t | y_{1:t-1})
◆ logGEv()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
virtual float_t BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::logGEv |
( |
const osv & |
yt, |
|
|
const ssv & |
xt, |
|
|
const cvsv & |
zt |
|
) |
| |
|
pure virtual |
Calculate gEv or logGEv.
- Parameters
-
yt | is a const Vec& describing the time t datum |
xt | is a const Vec& describing the time t state |
zt | is a const Vec& describing the time t covariate |
- Returns
- the density or log-density evaluation as a float_t
◆ logMuEv()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
virtual float_t BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::logMuEv |
( |
const ssv & |
x1, |
|
|
const cvsv & |
z1 |
|
) |
| |
|
pure virtual |
Calculate muEv or logmuEv.
- Parameters
-
x1 | is a const Vec& describing the state sample |
z1 | is a const Vec& describing the covariate sample |
- Returns
- the density or log-density evaluation as a float_t
◆ logQ1Ev()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
virtual float_t BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::logQ1Ev |
( |
const ssv & |
x1, |
|
|
const osv & |
y1, |
|
|
const cvsv & |
z1 |
|
) |
| |
|
pure virtual |
Calculate q1Ev or log q1Ev.
- Parameters
-
x1 | is a const Vec& describing the time 1 state sample |
y1 | is a const Vec& describing the time 1 datum |
z1 | is a const Vec& describing the time 1 covariate |
- Returns
- the density or log-density evaluation as a float_t
◆ q1Samp()
template<size_t nparts, size_t dimx, size_t dimy, size_t dimcov, typename resamp_t , typename float_t >
virtual ssv BSFilterWC< nparts, dimx, dimy, dimcov, resamp_t, float_t >::q1Samp |
( |
const osv & |
y1, |
|
|
const cvsv & |
z1 |
|
) |
| |
|
pure virtual |
Samples from time 1 proposal.
- Parameters
-
y1 | is a const Vec& representing the first observed datum |
z1 | is the const Vec& representing the first covariate |
- Returns
- the sample as a Vec
The documentation for this class was generated from the following file: