pf
bootstrap_filter.h
Go to the documentation of this file.
1 #ifndef BOOTSTRAP_FILTER_H
2 #define BOOTSTRAP_FILTER_H
3 
4 #include <array>
5 #include <vector>
6 
7 
8 #ifdef DROPPINGTHISINRPACKAGE
9  #include <RcppEigen.h>
10  // [[Rcpp::depends(RcppEigen)]]
11 #else
12  #include <Eigen/Dense>
13 #endif
14 
15 
16 #include "pf_base.h"
17 
18 
19 namespace pf {
20 
21 namespace filters {
22 
24 
34 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug=false>
35 class BSFilter : public bases::pf_base<float_t, dimy, dimx>
36 {
37 private:
38 
40  using ssv = Eigen::Matrix<float_t, dimx, 1>;
42  using osv = Eigen::Matrix<float_t, dimy, 1>; // obs size vec
44  using Mat = Eigen::Matrix<float_t, Eigen::Dynamic, Eigen::Dynamic>;
46  using arrayStates = std::array<ssv, nparts>;
48  using arrayFloat = std::array<float_t, nparts>;
49 
50 public:
55  BSFilter(const unsigned int &rs = 1);
56 
57 
61  virtual ~BSFilter();
62 
63 
68  float_t getLogCondLike() const;
69 
70 
77  void filter(const osv &data, const std::vector<std::function<const Mat(const ssv&)> >& fs = std::vector<std::function<const Mat(const ssv&)> >());
78 
79 
84  auto getExpectations () const -> std::vector<Mat>;
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 
126  virtual ssv fSamp (const ssv &xtm1) = 0;
127 
128 protected:
131 
134 
136  unsigned int m_now;
137 
140 
142  resamp_t m_resampler;
143 
145  std::vector<Mat> m_expectations;
146 
148  unsigned int m_resampSched;
149 };
150 
151 
152 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
153 BSFilter<nparts, dimx, dimy, resamp_t, float_t, debug>::BSFilter(const unsigned int &rs)
154  : m_now(0)
155  , m_logLastCondLike(0.0)
156  , m_resampSched(rs)
157 
158 {
159  std::fill(m_logUnNormWeights.begin(), m_logUnNormWeights.end(), 0.0);
160 }
161 
162 
163 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
165 
166 
167 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
168 void BSFilter<nparts, dimx, dimy, resamp_t, float_t, debug>::filter(const osv &dat, const std::vector<std::function<const Mat(const ssv&)> >& fs)
169 {
170 
171  if( m_now > 0)
172  {
173 
174  // try to iterate over particles all at once
175  ssv newSamp;
176  float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
177  arrayFloat oldLogUnNormWts = m_logUnNormWeights;
178  for(size_t ii = 0; ii < nparts; ++ii)
179  {
180  // update max of old logUnNormWts
181  if (m_logUnNormWeights[ii] > maxOldLogUnNormWts)
182  maxOldLogUnNormWts = m_logUnNormWeights[ii];
183 
184  // sample and get weight adjustments
185  newSamp = fSamp(m_particles[ii]);
186  m_logUnNormWeights[ii] += logGEv(dat, newSamp);
187 
188  // overwrite stuff
189  m_particles[ii] = newSamp;
190 
191  // print stuff if debug mode is on
192  #ifndef DROPPINGTHISINRPACKAGE
193  if constexpr(debug)
194  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
195  #endif
196  }
197 
198  // compute estimate of log p(y_t|y_{1:t-1}) with log-exp-sum trick
199  float_t maxNumer = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end()); //because you added log adjustments
200  float_t sumExp1(0.0);
201  float_t sumExp2(0.0);
202  for(size_t i = 0; i < nparts; ++i){
203  sumExp1 += std::exp(m_logUnNormWeights[i] - maxNumer);
204  sumExp2 += std::exp(oldLogUnNormWts[i] - maxOldLogUnNormWts); //1
205  }
206  m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
207 
208  // calculate expectations before you resample
209  int fId(0);
210  for(auto & h : fs){
211 
212  Mat testOutput = h(m_particles[0]);
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){
218  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - maxNumer);
219  weightNormConst += std::exp(m_logUnNormWeights[prtcl] - maxNumer);
220  }
221  m_expectations[fId] = numer/weightNormConst;
222 
223  // print stuff if debug mode is on
224  #ifndef DROPPINGTHISINRPACKAGE
225  if constexpr(debug)
226  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
227  #endif
228 
229  fId++;
230  }
231 
232  // resample if you should
233  if ( (m_now+1) % m_resampSched == 0)
235 
236  // advance time
237  m_now += 1;
238  }
239  else // (m_now == 0) //time 1
240  {
241  // only need to iterate over particles once
242  for(size_t ii = 0; ii < nparts; ++ii)
243  {
244  // sample particles
245  m_particles[ii] = q1Samp(dat);
247  m_logUnNormWeights[ii] += logGEv(dat, m_particles[ii]);
248  m_logUnNormWeights[ii] -= logQ1Ev(m_particles[ii], dat);
249 
250  // print stuff if debug mode is on
251  #ifndef DROPPINGTHISINRPACKAGE
252  if constexpr(debug)
253  std::cout << "time: " << m_now << ", transposed sample: " << m_particles[ii].transpose() << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
254  #endif
255  }
256 
257  // calculate log cond likelihood with log-exp-sum trick
258  float_t max = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
259  float_t sumExp(0.0);
260  for(size_t i = 0; i < nparts; ++i){
261  sumExp += std::exp(m_logUnNormWeights[i] - max);
262  }
263  m_logLastCondLike = -std::log(nparts) + max + std::log(sumExp);
264 
265  // calculate expectations before you resample
266  // paying mind to underflow
267  m_expectations.resize(fs.size());
268  unsigned int fId(0);
269  for(auto & h : fs){
270 
271  Mat testOutput = h(m_particles[0]);
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){
277  numer += h(m_particles[prtcl]) * std::exp( m_logUnNormWeights[prtcl] - (max) );
278  weightNormConst += std::exp( m_logUnNormWeights[prtcl] - (max) );
279  }
280  m_expectations[fId] = numer/weightNormConst;
281 
282  // print stuff if debug mode is on
283  #ifndef DROPPINGTHISINRPACKAGE
284  if constexpr(debug)
285  std::cout << "transposed expectation " << fId << ": " << m_expectations[fId].transpose() << "\n";
286  #endif
287 
288  fId++;
289  }
290 
291  // resample if you should
292  if ( (m_now+1) % m_resampSched == 0){
294  }
295 
296  // advance time step
297  m_now += 1;
298  }
299 
300 }
301 
302 
303 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
305 {
306  return m_logLastCondLike;
307 }
308 
309 
310 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
312 {
313  return m_expectations;
314 }
315 
316 
317 } // namespace filters
318 } // namespace pf
319 
320 #endif // BOOTSTRAP_FILTER_H
pf::bases::pf_base
Definition: pf_base.h:35
pf::filters::BSFilter::arrayStates
std::array< ssv, nparts > arrayStates
Definition: bootstrap_filter.h:46
pf::filters::BSFilter::fSamp
virtual ssv fSamp(const ssv &xtm1)=0
Sample from the state transition distribution.
pf::filters::BSFilter::m_now
unsigned int m_now
time point
Definition: bootstrap_filter.h:136
pf::filters::BSFilter::m_logUnNormWeights
arrayFloat m_logUnNormWeights
particle unnormalized weights
Definition: bootstrap_filter.h:133
pf_base.h
All non Rao-Blackwellized particle filters without covariates inherit from this.
pf::filters::BSFilter::logQ1Ev
virtual float_t logQ1Ev(const ssv &x1, const osv &y1)=0
Calculate q1Ev or log q1Ev.
pf::filters::BSFilter::Mat
Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > Mat
Definition: bootstrap_filter.h:44
pf::filters::BSFilter::BSFilter
BSFilter(const unsigned int &rs=1)
The constructor.
Definition: bootstrap_filter.h:153
pf::filters::BSFilter::m_resampSched
unsigned int m_resampSched
resampling schedule (e.g. resample every __ time points)
Definition: bootstrap_filter.h:148
pf::filters::BSFilter::m_resampler
resamp_t m_resampler
resampler object
Definition: bootstrap_filter.h:142
pf::filters::BSFilter::arrayFloat
std::array< float_t, nparts > arrayFloat
Definition: bootstrap_filter.h:48
pf::filters::BSFilter::~BSFilter
virtual ~BSFilter()
The (virtual) destructor.
Definition: bootstrap_filter.h:164
pf::filters::BSFilter::logGEv
virtual float_t logGEv(const osv &yt, const ssv &xt)=0
Calculate gEv or logGEv.
pf::filters::BSFilter
A base class for the bootstrap particle filter.
Definition: bootstrap_filter.h:35
pf::filters::BSFilter::m_expectations
std::vector< Mat > m_expectations
expectations E[h(x_t) | y_{1:t}] for user defined "h"s
Definition: bootstrap_filter.h:145
pf::filters::BSFilter::ssv
Eigen::Matrix< float_t, dimx, 1 > ssv
Definition: bootstrap_filter.h:40
pf::filters::BSFilter::getExpectations
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
pf::filters::BSFilter::osv
Eigen::Matrix< float_t, dimy, 1 > osv
Definition: bootstrap_filter.h:42
pf::filters::BSFilter::q1Samp
virtual ssv q1Samp(const osv &y1)=0
Samples from time 1 proposal.
pf::filters::BSFilter::getLogCondLike
float_t getLogCondLike() const
Returns the most recent (log-) conditional likelihood.
Definition: bootstrap_filter.h:304
pf::filters::BSFilter::m_particles
arrayStates m_particles
particle samples
Definition: bootstrap_filter.h:130
pf::filters::BSFilter::m_logLastCondLike
float_t m_logLastCondLike
log p(y_t|y_{1:t-1}) or log p(y1)
Definition: bootstrap_filter.h:139
pf::filters::BSFilter::logMuEv
virtual float_t logMuEv(const ssv &x1)=0
Calculate muEv or logmuEv.
pf::filters::BSFilter::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: bootstrap_filter.h:168