Go to the documentation of this file.    1 #ifndef BOOTSTRAP_FILTER_H 
    2 #define BOOTSTRAP_FILTER_H 
    8 #ifdef DROPPINGTHISINRPACKAGE 
   12     #include <Eigen/Dense> 
   34 template<
size_t nparts, 
size_t dimx, 
size_t dimy, 
typename resamp_t, 
typename float_t, 
bool debug=false>
 
   40     using ssv         = Eigen::Matrix<float_t, dimx, 1>; 
 
   42     using osv         = Eigen::Matrix<float_t, dimy, 1>; 
 
   44     using Mat         = Eigen::Matrix<float_t, Eigen::Dynamic, Eigen::Dynamic>;
 
   55     BSFilter(
const unsigned int &rs = 1);
 
   77     void filter(
const osv &data, 
const std::vector<std::function<
const Mat(
const ssv&)> >& fs = std::vector<std::function<
const Mat(
const ssv&)> >()); 
 
  118     virtual float_t 
logGEv (const 
osv &yt, const 
ssv &xt ) = 0;
 
  152 template<
size_t nparts, 
size_t dimx, 
size_t dimy, typename resamp_t, typename float_t, 
bool debug>
 
  163 template<
size_t nparts, 
size_t dimx, 
size_t dimy, 
typename resamp_t, 
typename float_t, 
bool debug>
 
  167 template<
size_t nparts, 
size_t dimx, 
size_t dimy, 
typename resamp_t, 
typename float_t, 
bool debug>
 
  176         float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
 
  178         for(
size_t ii = 0; ii < nparts; ++ii)
 
  192             #ifndef DROPPINGTHISINRPACKAGE  
  200         float_t sumExp1(0.0);
 
  201         float_t sumExp2(0.0);
 
  202         for(
size_t i = 0; i < nparts; ++i){
 
  204             sumExp2 += std::exp(oldLogUnNormWts[i] - maxOldLogUnNormWts);  
 
  206         m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
 
  213             unsigned int rows = testOutput.rows();
 
  214             unsigned int cols = testOutput.cols();
 
  215             Mat numer = Mat::Zero(rows,cols);
 
  216             float_t weightNormConst (0.0);
 
  217             for(
size_t prtcl = 0; prtcl < nparts; ++prtcl){ 
 
  224             #ifndef DROPPINGTHISINRPACKAGE 
  226                 std::cout << 
"transposed expectation " << fId << 
": " << 
m_expectations[fId].transpose() << 
"\n";
 
  242         for(
size_t ii = 0; ii < nparts; ++ii)
 
  251             #ifndef DROPPINGTHISINRPACKAGE 
  260         for(
size_t i = 0; i < nparts; ++i){
 
  272             unsigned int rows = testOutput.rows();
 
  273             unsigned int cols = testOutput.cols();
 
  274             Mat numer = Mat::Zero(rows,cols);
 
  275             float_t weightNormConst (0.0);
 
  276             for(
size_t prtcl = 0; prtcl < nparts; ++prtcl){ 
 
  283             #ifndef DROPPINGTHISINRPACKAGE 
  285                 std::cout << 
"transposed expectation " << fId << 
": " << 
m_expectations[fId].transpose() << 
"\n";
 
  303 template<
size_t nparts, 
size_t dimx, 
size_t dimy, 
typename resamp_t, 
typename float_t, 
bool debug>
 
  310 template<
size_t nparts, 
size_t dimx, 
size_t dimy, 
typename resamp_t, 
typename float_t, 
bool debug>
 
  320 #endif // BOOTSTRAP_FILTER_H 
 
 
std::array< ssv, nparts > arrayStates
Definition: bootstrap_filter.h:46
 
virtual ssv fSamp(const ssv &xtm1)=0
Sample from the state transition distribution.
 
unsigned int m_now
time point
Definition: bootstrap_filter.h:136
 
arrayFloat m_logUnNormWeights
particle unnormalized weights
Definition: bootstrap_filter.h:133
 
All non Rao-Blackwellized particle filters without covariates inherit from this.
 
virtual float_t logQ1Ev(const ssv &x1, const osv &y1)=0
Calculate q1Ev or log q1Ev.
 
Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > Mat
Definition: bootstrap_filter.h:44
 
BSFilter(const unsigned int &rs=1)
The constructor.
Definition: bootstrap_filter.h:153
 
unsigned int m_resampSched
resampling schedule (e.g. resample every __ time points)
Definition: bootstrap_filter.h:148
 
resamp_t m_resampler
resampler object
Definition: bootstrap_filter.h:142
 
std::array< float_t, nparts > arrayFloat
Definition: bootstrap_filter.h:48
 
virtual ~BSFilter()
The (virtual) destructor.
Definition: bootstrap_filter.h:164
 
virtual float_t logGEv(const osv &yt, const ssv &xt)=0
Calculate gEv or logGEv.
 
A base class for the bootstrap particle filter.
Definition: bootstrap_filter.h:35
 
std::vector< Mat > m_expectations
expectations E[h(x_t) | y_{1:t}] for user defined "h"s
Definition: bootstrap_filter.h:145
 
Eigen::Matrix< float_t, dimx, 1 > ssv
Definition: bootstrap_filter.h:40
 
auto getExpectations() const -> std::vector< Mat >
return all stored expectations (taken with respect to $p(x_t|y_{1:t})$
Definition: bootstrap_filter.h:311
 
Eigen::Matrix< float_t, dimy, 1 > osv
Definition: bootstrap_filter.h:42
 
virtual ssv q1Samp(const osv &y1)=0
Samples from time 1 proposal.
 
float_t getLogCondLike() const
Returns the most recent (log-) conditional likelihood.
Definition: bootstrap_filter.h:304
 
arrayStates m_particles
particle samples
Definition: bootstrap_filter.h:130
 
float_t m_logLastCondLike
log p(y_t|y_{1:t-1}) or log p(y1)
Definition: bootstrap_filter.h:139
 
virtual float_t logMuEv(const ssv &x1)=0
Calculate muEv or logmuEv.
 
void filter(const osv &data, const std::vector< std::function< const Mat(const ssv &)> > &fs=std::vector< std::function< const Mat(const ssv &)> >())
updates filtering distribution on a new datapoint. Optionally stores expectations of functionals.
Definition: bootstrap_filter.h:168