pf
sisr_filter.h
Go to the documentation of this file.
1 #ifndef SISR_FILTER_H
2 #define SISR_FILTER_H
3 
4 #include <array>
5 
6 #ifdef DROPPINGTHISINRPACKAGE
7  #include <RcppEigen.h>
8  // [[Rcpp::depends(RcppEigen)]]
9 #else
10  #include <Eigen/Dense>
11 #endif
12 
13 #include "pf_base.h"
14 
15 
16 namespace pf {
17 
18 namespace filters {
19 
21 
33 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug=false>
34 class SISRFilter : public bases::pf_base<float_t, dimy, dimx>
35 {
36 private:
37 
39  using ssv = Eigen::Matrix<float_t, dimx, 1>;
41  using osv = Eigen::Matrix<float_t, dimy, 1>; // obs size vec
43  using Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic>;
45  using arrayStates = std::array<ssv, nparts>;
47  using arrayfloat_t = std::array<float_t, nparts>;
48 
49 public:
50 
55  SISRFilter(const unsigned int &rs=1);
56 
57 
61  virtual ~SISRFilter();
62 
63 
68  float_t getLogCondLike() const;
69 
70 
75  std::vector<Mat> getExpectations() const;
76 
77 
84  void filter(const osv &data, const std::vector<std::function<const Mat(const ssv&)> >& fs = std::vector<std::function<const Mat(const ssv&)> >());
85 
86 
92  virtual float_t logMuEv (const ssv &x1) = 0;
93 
94 
100  virtual ssv q1Samp (const osv &y1) = 0;
101 
102 
109  virtual float_t logQ1Ev (const ssv &x1, const osv &y1 ) = 0;
110 
111 
118  virtual float_t logGEv (const osv &yt, const ssv &xt ) = 0;
119 
120 
127  virtual float_t logFEv (const ssv &xt, const ssv &xtm1 ) = 0;
128 
129 
136  virtual ssv qSamp (const ssv &xtm1, const osv &yt ) = 0;
137 
138 
146  virtual float_t logQEv (const ssv &xt, const ssv &xtm1, const osv &yt ) = 0;
147 
148 protected:
149 
152 
155 
157  unsigned int m_now;
158 
161 
163  resamp_t m_resampler;
164 
166  std::vector<Mat> m_expectations; // stores any sample averages the user wants
167 
169  unsigned int m_resampSched;
170 
171 
176 };
177 
178 
179 
180 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
182  : m_now(0)
183  , m_logLastCondLike(0.0)
184  , m_resampSched(rs)
185 {
186  std::fill(m_logUnNormWeights.begin(), m_logUnNormWeights.end(), 0.0); // log(1) = 0
187 }
188 
189 
190 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
192 
193 
194 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
196 {
197  return m_logLastCondLike;
198 }
199 
200 
201 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
203 {
204  return m_expectations;
205 }
206 
207 
208 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
209 void SISRFilter<nparts,dimx,dimy,resamp_t,float_t, debug>::filter(const osv &data, const std::vector<std::function<const Mat(const ssv&)> >& fs)
210 {
211 
212  if(m_now > 0)
213  {
214 
215  // try to iterate over particles all at once
216  ssv newSamp;
217  arrayfloat_t oldLogUnNormWts = m_logUnNormWeights;
218  float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
219  for(size_t ii = 0; ii < nparts; ++ii)
220  {
221 
222  // update max of old logUnNormWts before you change the element
223  if (m_logUnNormWeights[ii] > maxOldLogUnNormWts)
224  maxOldLogUnNormWts = m_logUnNormWeights[ii];
225 
226  // sample and get weight adjustments
227  newSamp = qSamp(m_particles[ii], data);
228  m_logUnNormWeights[ii] += logFEv(newSamp, m_particles[ii]);
229  m_logUnNormWeights[ii] += logGEv(data, newSamp);
230  m_logUnNormWeights[ii] -= logQEv(newSamp, m_particles[ii], data);
231 
232  // overwrite stuff
233  m_particles[ii] = newSamp;
234 
235  #ifndef DROPPINGTHISINRPACKAGE
236  if constexpr(debug)
237  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
238  #endif
239 
240  }
241 
242  // compute estimate of log p(y_t|y_{1:t-1}) with log-exp-sum trick
243  float_t maxNumer = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end()); //because you added log adjustments
244  float_t sumExp1(0.0);
245  float_t sumExp2(0.0);
246  for(size_t i = 0; i < nparts; ++i){
247  sumExp1 += std::exp(m_logUnNormWeights[i] - maxNumer);
248  sumExp2 += std::exp(oldLogUnNormWts[i] - maxOldLogUnNormWts);
249  }
250  m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
251 
252  // calculate expectations before you resample
253  unsigned int fId(0);
254  float_t weightNormConst(0.0);
255  for(auto & h : fs){ // iterate over all functions
256 
257  Mat testOut = h(m_particles[0]);
258  unsigned int rows = testOut.rows();
259  unsigned int cols = testOut.cols();
260  Mat numer = Mat::Zero(rows,cols);
261  float_t denom(0.0);
262 
263  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){ // iterate over all particles
264  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - maxNumer );
265  denom += std::exp(m_logUnNormWeights[prtcl] - maxNumer);
266  }
267  m_expectations[fId] = numer/denom;
268 
269  // print stuff if debug mode is on
270  #ifndef DROPPINGTHISINRPACKAGE
271  if constexpr(debug)
272  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
273  #endif
274 
275  fId++;
276  }
277 
278  // resample if you should
279  if( (m_now + 1) % m_resampSched == 0)
280  m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
281 
282  // advance time
283  m_now += 1;
284  }
285  else // (m_now == 0) //time 1
286  {
287 
288  // only need to iterate over particles once
289  float_t sumWts(0.0);
290  for(size_t ii = 0; ii < nparts; ++ii)
291  {
292  // sample particles
293  m_particles[ii] = q1Samp(data);
294  m_logUnNormWeights[ii] += logMuEv(m_particles[ii]);
295  m_logUnNormWeights[ii] += logGEv(data, m_particles[ii]);
296  m_logUnNormWeights[ii] -= logQ1Ev(m_particles[ii], data);
297 
298  #ifndef DROPPINGTHISINRPACKAGE
299  if constexpr(debug)
300  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
301  #endif
302 
303 
304  }
305 
306  // calculate log cond likelihood with log-exp-sum trick
307  float_t max = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
308  float_t sumExp(0.0);
309  for(size_t i = 0; i < nparts; ++i){
310  sumExp += std::exp(m_logUnNormWeights[i] - max);
311  }
312  m_logLastCondLike = -std::log(nparts) + max + std::log(sumExp);
313 
314  // calculate expectations before you resample
315  m_expectations.resize(fs.size());
316  unsigned int fId(0);
317  for(auto & h : fs){
318 
319  Mat testOut = h(m_particles[0]);
320  unsigned int rows = testOut.rows();
321  unsigned int cols = testOut.cols();
322  Mat numer = Mat::Zero(rows,cols);
323  float_t denom(0.0);
324 
325  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){
326  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - max);
327  denom += std::exp(m_logUnNormWeights[prtcl] - max);
328  }
329  m_expectations[fId] = numer/denom;
330 
331  // print stuff if debug mode is on
332  #ifndef DROPPINGTHISINRPACKAGE
333  if constexpr(debug)
334  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
335  #endif
336 
337  fId++;
338  }
339 
340  // resample if you should
341  if( (m_now + 1) % m_resampSched == 0)
342  m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
343 
344  // advance time step
345  m_now += 1;
346  }
347 
348 }
349 
350 
352 
366 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug=false>
367 class SISRFilterCRN : public bases::pf_base_crn<float_t, dimy, dimx, dimu, dimur, nparts>
368 {
369 private:
370 
372  using ssv = Eigen::Matrix<float_t, dimx, 1>;
374  using osv = Eigen::Matrix<float_t, dimy, 1>; // obs size vec
376  using usv = Eigen::Matrix<float_t, dimu, 1>;
378  using usvr = Eigen::Matrix<float_t, dimur, 1>;
380  using Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic>;
382  using arrayStates = std::array<ssv, nparts>;
384  using arrayUs = std::array<usv, nparts>;
386  using arrayfloat_t = std::array<float_t, nparts>;
387 
388 public:
389 
394  SISRFilterCRN(const unsigned int &rs=1);
395 
396 
400  virtual ~SISRFilterCRN();
401 
402 
407  float_t getLogCondLike() const;
408 
409 
414  std::vector<Mat> getExpectations() const;
415 
416 
425  void filter(const osv &data, const arrayUs& Uarr, const usvr& Uresamp, const std::vector<std::function<const Mat(const ssv&)> >& fs = std::vector<std::function<const Mat(const ssv&)> >());
426 
427 
433  virtual float_t logMuEv (const ssv &x1) = 0;
434 
435 
442  virtual ssv Xi1 (const usv& U, const osv &y1) = 0;
443 
444 
451  virtual float_t logQ1Ev (const ssv &x1, const osv &y1 ) = 0;
452 
453 
460  virtual float_t logGEv (const osv &yt, const ssv &xt ) = 0;
461 
462 
469  virtual float_t logFEv (const ssv &xt, const ssv &xtm1 ) = 0;
470 
471 
479  virtual ssv Xit (const ssv &xtm1, const usv& U, const osv &yt) = 0;
480 
481 
489  virtual float_t logQEv (const ssv &xt, const ssv &xtm1, const osv &yt ) = 0;
490 
491 protected:
492 
495 
498 
500  unsigned int m_now;
501 
504 
506  resamp_t m_resampler;
507 
509  std::vector<Mat> m_expectations; // stores any sample averages the user wants
510 
512  unsigned int m_resampSched;
513 
514 
515 };
516 
517 
518 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug>
520  : m_now(0)
521  , m_logLastCondLike(0.0)
522  , m_resampSched(rs)
523 {
524  std::fill(m_logUnNormWeights.begin(), m_logUnNormWeights.end(), 0.0); // log(1) = 0
525 }
526 
527 
528 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug>
530 
531 
532 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug>
534 {
535  return m_logLastCondLike;
536 }
537 
538 
539 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug>
541 {
542  return m_expectations;
543 }
544 
545 
546 template<size_t nparts, size_t dimx, size_t dimy, size_t dimu, size_t dimur, typename resamp_t, typename float_t, bool debug>
547 void SISRFilterCRN<nparts,dimx,dimy,dimu,dimur,resamp_t,float_t, debug>::filter(const osv &data, const arrayUs& Uarr, const usvr& Uresamp, const std::vector<std::function<const Mat(const ssv&)> >& fs)
548 {
549 
550  if(m_now > 0)
551  {
552 
553  // try to iterate over particles all at once
554  ssv newSamp;
555  arrayfloat_t oldLogUnNormWts = m_logUnNormWeights;
556  float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
557  for(size_t ii = 0; ii < nparts; ++ii)
558  {
559 
560  // update max of old logUnNormWts before you change the element
561  if (m_logUnNormWeights[ii] > maxOldLogUnNormWts)
562  maxOldLogUnNormWts = m_logUnNormWeights[ii];
563 
564  // sample and get weight adjustments
565  newSamp = Xit(m_particles[ii], Uarr[ii], data);
566  m_logUnNormWeights[ii] += logFEv(newSamp, m_particles[ii]);
567  m_logUnNormWeights[ii] += logGEv(data, newSamp);
568  m_logUnNormWeights[ii] -= logQEv(newSamp, m_particles[ii], data);
569 
570  // overwrite stuff
571  m_particles[ii] = newSamp;
572 
573  #ifndef DROPPINGTHISINRPACKAGE
574  if constexpr(debug)
575  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
576  #endif
577  }
578 
579  // compute estimate of log p(y_t|y_{1:t-1}) with log-exp-sum trick
580  float_t maxNumer = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end()); //because you added log adjustments
581  float_t sumExp1(0.0);
582  float_t sumExp2(0.0);
583  for(size_t i = 0; i < nparts; ++i){
584  sumExp1 += std::exp(m_logUnNormWeights[i] - maxNumer);
585  sumExp2 += std::exp(oldLogUnNormWts[i] - maxOldLogUnNormWts);
586  }
587  m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
588 
589  // calculate expectations before you resample
590  unsigned int fId(0);
591  float_t weightNormConst(0.0);
592  for(auto & h : fs){ // iterate over all functions
593 
594  Mat testOut = h(m_particles[0]);
595  unsigned int rows = testOut.rows();
596  unsigned int cols = testOut.cols();
597  Mat numer = Mat::Zero(rows,cols);
598  float_t denom(0.0);
599 
600  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){ // iterate over all particles
601  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - maxNumer );
602  denom += std::exp(m_logUnNormWeights[prtcl] - maxNumer);
603  }
604  m_expectations[fId] = numer/denom;
605 
606  // print stuff if debug mode is on
607  #ifndef DROPPINGTHISINRPACKAGE
608  if constexpr(debug)
609  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
610  #endif
611 
612  fId++;
613  }
614 
615  // resample if you should
616  if( (m_now + 1) % m_resampSched == 0)
617  m_resampler.resampLogWts(m_particles, m_logUnNormWeights, Uresamp);
618 
619  // advance time
620  m_now += 1;
621  }
622  else // (m_now == 0) //time 1
623  {
624 
625  // only need to iterate over particles once
626  float_t sumWts(0.0);
627  for(size_t ii = 0; ii < nparts; ++ii)
628  {
629  // sample particles
630  m_particles[ii] = Xi1(Uarr[ii], data);
631  m_logUnNormWeights[ii] += logMuEv(m_particles[ii]);
632  m_logUnNormWeights[ii] += logGEv(data, m_particles[ii]);
633  m_logUnNormWeights[ii] -= logQ1Ev(m_particles[ii], data);
634 
635  #ifndef DROPPINGTHISINRPACKAGE
636  if constexpr(debug)
637  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
638  #endif
639 
640  }
641 
642  // calculate log cond likelihood with log-exp-sum trick
643  float_t max = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
644  float_t sumExp(0.0);
645  for(size_t i = 0; i < nparts; ++i){
646  sumExp += std::exp(m_logUnNormWeights[i] - max);
647  }
648  m_logLastCondLike = -std::log(nparts) + max + std::log(sumExp);
649 
650  // calculate expectations before you resample
651  m_expectations.resize(fs.size());
652  unsigned int fId(0);
653  for(auto & h : fs){
654 
655  Mat testOut = h(m_particles[0]);
656  unsigned int rows = testOut.rows();
657  unsigned int cols = testOut.cols();
658  Mat numer = Mat::Zero(rows,cols);
659  float_t denom(0.0);
660 
661  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){
662  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - max);
663  denom += std::exp(m_logUnNormWeights[prtcl] - max);
664  }
665  m_expectations[fId] = numer/denom;
666 
667  // print stuff if debug mode is on
668  #ifndef DROPPINGTHISINRPACKAGE
669  if constexpr(debug)
670  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
671  #endif
672 
673  fId++;
674  }
675 
676  // resample if you should
677  if( (m_now + 1) % m_resampSched == 0)
678  m_resampler.resampLogWts(m_particles, m_logUnNormWeights, Uresamp);
679 
680  // advance time step
681  m_now += 1;
682  }
683 
684 }
685 
686 } // namespace filters
687 } //namespace pf
688 
689 #endif //SISR_FILTER_H
pf::filters::SISRFilterCRN::m_particles
arrayStates m_particles
particle samples
Definition: sisr_filter.h:494
pf::bases::pf_base_crn
Definition: pf_base.h:724
pf::filters::SISRFilterCRN::Xit
virtual ssv Xit(const ssv &xtm1, const usv &U, const osv &yt)=0
"Samples" from the proposal/instrumental/importance density at time t. Really, it maps the normal ran...
pf::bases::pf_base
Definition: pf_base.h:35
pf::filters::SISRFilterCRN::filter
void filter(const osv &data, const arrayUs &Uarr, const usvr &Uresamp, 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: sisr_filter.h:547
pf::filters::SISRFilter::getExpectations
std::vector< Mat > getExpectations() const
return all stored expectations (taken with respect to $p(x_t|y_{1:t})$
Definition: sisr_filter.h:202
pf::filters::SISRFilterCRN::usvr
Eigen::Matrix< float_t, dimur, 1 > usvr
Definition: sisr_filter.h:378
pf::filters::SISRFilter
A base class for the Sequential Important Sampling with Resampling (SISR).
Definition: sisr_filter.h:34
pf::filters::SISRFilter::logMuEv
virtual float_t logMuEv(const ssv &x1)=0
Calculate muEv or logmuEv.
pf::filters::SISRFilter::qSamp
virtual ssv qSamp(const ssv &xtm1, const osv &yt)=0
Samples from the proposal/instrumental/importance density at time t.
pf::filters::SISRFilterCRN::Xi1
virtual ssv Xi1(const usv &U, const osv &y1)=0
"Samples" from time 1 proposal. Really, it maps the normal random vector into the sample.
pf::filters::SISRFilter::logFEv
virtual float_t logFEv(const ssv &xt, const ssv &xtm1)=0
Evaluates the state transition density.
pf::filters::SISRFilter::SISRFilter
SISRFilter(const unsigned int &rs=1)
The (one and only) constructor.
Definition: sisr_filter.h:181
pf::filters::SISRFilterCRN::ssv
Eigen::Matrix< float_t, dimx, 1 > ssv
Definition: sisr_filter.h:372
pf::filters::SISRFilterCRN::~SISRFilterCRN
virtual ~SISRFilterCRN()
The (virtual) destructor.
Definition: sisr_filter.h:529
pf::filters::SISRFilter::logQEv
virtual float_t logQEv(const ssv &xt, const ssv &xtm1, const osv &yt)=0
Evaluates the proposal/instrumental/importance density/pmf.
pf::filters::SISRFilter::logGEv
virtual float_t logGEv(const osv &yt, const ssv &xt)=0
Calculate gEv or logGEv.
pf::filters::SISRFilterCRN::logQ1Ev
virtual float_t logQ1Ev(const ssv &x1, const osv &y1)=0
Calculate q1Ev or log q1Ev.
pf::filters::SISRFilter::m_particles
arrayStates m_particles
particle samples
Definition: sisr_filter.h:151
pf_base.h
All non Rao-Blackwellized particle filters without covariates inherit from this.
pf::filters::SISRFilterCRN::m_logLastCondLike
float_t m_logLastCondLike
log p(y_t|y_{1:t-1}) or log p(y1)
Definition: sisr_filter.h:503
pf::filters::SISRFilterCRN::Mat
Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > Mat
Definition: sisr_filter.h:380
pf::filters::SISRFilterCRN::logQEv
virtual float_t logQEv(const ssv &xt, const ssv &xtm1, const osv &yt)=0
Evaluates the proposal/instrumental/importance density/pmf.
pf::filters::SISRFilter::arrayfloat_t
std::array< float_t, nparts > arrayfloat_t
Definition: sisr_filter.h:47
pf::filters::SISRFilter::m_logLastCondLike
float_t m_logLastCondLike
log p(y_t|y_{1:t-1}) or log p(y1)
Definition: sisr_filter.h:160
pf::filters::SISRFilter::m_resampSched
unsigned int m_resampSched
resampling schedule (e.g. resample every __ time points)
Definition: sisr_filter.h:169
pf::filters::SISRFilterCRN::arrayfloat_t
std::array< float_t, nparts > arrayfloat_t
Definition: sisr_filter.h:386
pf::filters::SISRFilter::arrayStates
std::array< ssv, nparts > arrayStates
Definition: sisr_filter.h:45
pf::filters::SISRFilter::filter
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: sisr_filter.h:209
pf::filters::SISRFilter::m_now
unsigned int m_now
current time point
Definition: sisr_filter.h:157
pf::filters::SISRFilterCRN::arrayStates
std::array< ssv, nparts > arrayStates
Definition: sisr_filter.h:382
pf::filters::SISRFilterCRN::m_resampSched
unsigned int m_resampSched
resampling schedule (e.g. resample every __ time points)
Definition: sisr_filter.h:512
pf::filters::SISRFilter::Mat
Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > Mat
Definition: sisr_filter.h:43
pf::filters::SISRFilter::getLogCondLike
float_t getLogCondLike() const
Returns the most recent (log-) conditional likelihood.
Definition: sisr_filter.h:195
pf::filters::SISRFilter::m_resampler
resamp_t m_resampler
resampling object
Definition: sisr_filter.h:163
pf::filters::SISRFilter::ssv
Eigen::Matrix< float_t, dimx, 1 > ssv
Definition: sisr_filter.h:39
pf::filters::SISRFilterCRN
A base class for the Sequential Important Sampling with Resampling (SISR). Uses normal common random ...
Definition: sisr_filter.h:367
pf::filters::SISRFilterCRN::getLogCondLike
float_t getLogCondLike() const
Returns the most recent (log-) conditional likelihood.
Definition: sisr_filter.h:533
pf::filters::SISRFilter::q1Samp
virtual ssv q1Samp(const osv &y1)=0
Samples from time 1 proposal.
pf::filters::SISRFilterCRN::osv
Eigen::Matrix< float_t, dimy, 1 > osv
Definition: sisr_filter.h:374
pf::filters::SISRFilterCRN::m_now
unsigned int m_now
current time point
Definition: sisr_filter.h:500
pf::filters::SISRFilterCRN::arrayUs
std::array< usv, nparts > arrayUs
Definition: sisr_filter.h:384
pf::filters::SISRFilter::m_expectations
std::vector< Mat > m_expectations
expectations E[h(x_t) | y_{1:t}] for user defined "h"s
Definition: sisr_filter.h:166
pf::filters::SISRFilterCRN::usv
Eigen::Matrix< float_t, dimu, 1 > usv
Definition: sisr_filter.h:376
pf::filters::SISRFilterCRN::SISRFilterCRN
SISRFilterCRN(const unsigned int &rs=1)
The (one and only) constructor.
Definition: sisr_filter.h:519
pf::filters::SISRFilter::~SISRFilter
virtual ~SISRFilter()
The (virtual) destructor.
Definition: sisr_filter.h:191
pf::filters::SISRFilterCRN::logMuEv
virtual float_t logMuEv(const ssv &x1)=0
Calculate muEv or logmuEv.
pf::filters::SISRFilter::m_logUnNormWeights
arrayfloat_t m_logUnNormWeights
particle weights
Definition: sisr_filter.h:154
pf::filters::SISRFilterCRN::m_resampler
resamp_t m_resampler
resampling object
Definition: sisr_filter.h:506
pf::filters::SISRFilterCRN::getExpectations
std::vector< Mat > getExpectations() const
return all stored expectations (taken with respect to $p(x_t|y_{1:t})$
Definition: sisr_filter.h:540
pf::filters::SISRFilterCRN::logGEv
virtual float_t logGEv(const osv &yt, const ssv &xt)=0
Calculate gEv or logGEv.
pf::filters::SISRFilter::osv
Eigen::Matrix< float_t, dimy, 1 > osv
Definition: sisr_filter.h:41
pf::filters::SISRFilterCRN::m_expectations
std::vector< Mat > m_expectations
expectations E[h(x_t) | y_{1:t}] for user defined "h"s
Definition: sisr_filter.h:509
pf::filters::SISRFilter::logQ1Ev
virtual float_t logQ1Ev(const ssv &x1, const osv &y1)=0
Calculate q1Ev or log q1Ev.
pf::filters::SISRFilterCRN::logFEv
virtual float_t logFEv(const ssv &xt, const ssv &xtm1)=0
Evaluates the state transition density.
pf::filters::SISRFilterCRN::m_logUnNormWeights
arrayfloat_t m_logUnNormWeights
particle weights
Definition: sisr_filter.h:497