pf
|
Another class template for Gamma filtering, but this time. More...
#include <cf_filters.h>
Public Types | |
using | psv = Eigen::Matrix< float_t, dim_pred, 1 > |
"predictor size vector" | |
using | bsm = Eigen::Matrix< float_t, dim_obs, dim_pred > |
"beta size matrix" | |
using | tsv = Eigen::Matrix< float_t, 2, 1 > |
"two by 1 vector" to store size and shapes of gamma distributions | |
using | osv = Eigen::Matrix< float_t, dim_obs, 1 > |
"observation size vector" | |
using | osm = Eigen::Matrix< float_t, dim_obs, dim_obs > |
"observation size matrix" | |
![]() | |
using | obs_sized_vec = Eigen::Matrix< float_t, dimobs, 1 > |
using | state_sized_vec = Eigen::Matrix< float_t, dimstate, 1 > |
Public Member Functions | |
multivGamFilter (const float_t &nOneTilde, const float_t &dOneTilde) | |
Constructor. More... | |
virtual | ~multivGamFilter () |
The (virtual) desuctor. | |
float_t | getLogCondLike () const |
Get the latest conditional likelihood. More... | |
tsv | getFilterVec () const |
Get the current filter vector. More... | |
void | update (const osv &yt, const psv &xt, const bsm &B, const osm &Sigma, const float_t &delta) |
Perform a filtering update. More... | |
osv | getFcastMean (const psv &xtp1, const bsm &B, const osm &Sigma, const float_t &delta) |
Get the forecast mean (assuming filtering has been performed already) More... | |
osm | getFcastCov (const psv &xtp1, const bsm &B, const osm &Sigma, const float_t &delta) |
Get the forecast covariance matrix (assuming filtering has been performed already) More... | |
![]() | |
virtual | ~cf_filter () |
The (virtual) destructor. | |
virtual float_t | getLogCondLike () const=0 |
returns the log of the most recent conditional likelihood More... | |
Private Attributes | |
tsv | m_filtVec |
filter vector (shape and rate) | |
float_t | m_lastLogCondLike |
last log of the conditional likelihood | |
bool | m_fresh |
has data been observed? | |
Another class template for Gamma filtering, but this time.
pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::multivGamFilter | ( | const float_t & | nOneTilde, |
const float_t & | dOneTilde | ||
) |
Constructor.
nOneTilde | degrees of freedom for time 1 prior. |
dOneTilde | rate parameter for time 1 prior. |
auto pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::getFcastCov | ( | const psv & | xtp1, |
const bsm & | B, | ||
const osm & | Sigma, | ||
const float_t & | delta | ||
) |
Get the forecast covariance matrix (assuming filtering has been performed already)
gets the forecast covariance matrix!
xtp1 | the next time period's predictor vector |
B | the loadings matrix |
Sigma | the observation "shape" matrix |
delta | between 0 and 1 the discount parameter |
auto pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::getFcastMean | ( | const psv & | xtp1, |
const bsm & | B, | ||
const osm & | Sigma, | ||
const float_t & | delta | ||
) |
Get the forecast mean (assuming filtering has been performed already)
gets the forecast mean!
xtp1 | the next time period's predictor vector |
B | the loadings matrix |
Sigma | the observation "shape" matrix |
delta | between 0 and 1 the discount parameter |
auto pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::getFilterVec |
Get the current filter vector.
get the current filtering distribution. First element is the shape, second is the rate.
auto pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::getLogCondLike |
Get the latest conditional likelihood.
void pf::filters::multivGamFilter< dim_obs, dim_pred, float_t >::update | ( | const osv & | yt, |
const psv & | xt, | ||
const bsm & | B, | ||
const osm & | Sigma, | ||
const float_t & | delta | ||
) |
Perform a filtering update.
Perform a Gamma filter update.
yt | the most recent dependent random variable |
xt | the most recent predictor vector |
B | the loadings matrix |
Sigma | the observation "shape" matrix. |
delta | between 0 and 1 the discount parameter |