pf
auxiliary_pf.h
Go to the documentation of this file.
1 #ifndef AUXILIARY_PF_H
2 #define AUXILIARY_PF_H
3 
4 #include <array> //array
5 #include <vector> // vector
6 #include <functional> // function
7 
8 #ifdef DROPPINGTHISINRPACKAGE
9  #include <RcppEigen.h>
10  // [[Rcpp::depends(RcppEigen)]]
11 #else
12  #include <Eigen/Dense>
13 #endif
14 
15 #include <cmath>
16 
17 #include "pf_base.h"
18 #include "rv_samp.h" // for k_generator
19 
20 
21 namespace pf {
22 
23 namespace filters {
24 
26 
38 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug=false>
39 class APF : public bases::pf_base<float_t, dimy, dimx>
40 {
41 private:
42 
44  using ssv = Eigen::Matrix<float_t,dimx,1>;
46  using osv = Eigen::Matrix<float_t,dimy,1>;
48  using Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic>;
50  using arrayfloat_t = std::array<float_t, nparts>;
52  using arrayVec = std::array<ssv, nparts>;
54  using arrayUInt = std::array<unsigned int, nparts>;
55 
56 public:
57 
62  APF(const unsigned int &rs=1);
63 
64 
68  virtual ~APF();
69 
74  float_t getLogCondLike () const;
75 
76 
81  std::vector<Mat> getExpectations () const;
82 
83 
89  void filter(const osv &data, const std::vector<std::function<const Mat(const ssv&)> >& fs = std::vector<std::function<const Mat(const ssv&)> >());
90 
91 
97  virtual float_t logMuEv (const ssv &x1 ) = 0;
98 
99 
105  virtual ssv propMu (const ssv &xtm1 ) = 0;
106 
107 
113  virtual ssv q1Samp (const osv &y1) = 0;
114 
115 
121  virtual ssv fSamp (const ssv &xtm1) = 0;
122 
123 
130  virtual float_t logQ1Ev (const ssv &x1, const osv &y1) = 0;
131 
132 
139  virtual float_t logGEv (const osv &yt, const ssv &xt) = 0;
140 
141 
142 protected:
144  std::array<ssv,nparts> m_particles;
145 
147  std::array<float_t,nparts> m_logUnNormWeights;
148 
150  unsigned int m_now;
151 
154 
156  unsigned int m_rs;
157 
159  resamp_t m_resampler;
160 
163 
165  std::vector<Mat> m_expectations;
166 
167 };
168 
169 
170 
171 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
173  : m_now(0)
174  , m_logLastCondLike(0.0)
175  , m_rs(rs)
176 {
177  std::fill(m_logUnNormWeights.begin(), m_logUnNormWeights.end(), 0.0);
178 }
179 
180 
181 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
183 
184 
185 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
186 void APF<nparts, dimx, dimy, resamp_t, float_t, debug>::filter(const osv &data, const std::vector<std::function<const Mat(const ssv&)> >& fs)
187 {
188 
189  if(m_now > 0)
190  {
191 
192  // set up "first stage weights" to make k index sampler
193  arrayfloat_t logFirstStageUnNormWeights = m_logUnNormWeights;
194  arrayVec oldPartics = m_particles;
195  float_t m3(-std::numeric_limits<float_t>::infinity());
196  float_t m2(-std::numeric_limits<float_t>::infinity());
197  for(size_t ii = 0; ii < nparts; ++ii)
198  {
199  // update m3
200  if(m_logUnNormWeights[ii] > m3)
201  m3 = m_logUnNormWeights[ii];
202 
203  // sample
204  ssv xtm1 = oldPartics[ii];
205  logFirstStageUnNormWeights[ii] += logGEv(data, propMu(xtm1));
206 
207  // accumulate things
208  if(logFirstStageUnNormWeights[ii] > m2)
209  m2 = logFirstStageUnNormWeights[ii];
210 
211  // print stuff if debug mode is on
212  #ifndef DROPPINGTHISINRPACKAGE
213  if constexpr(debug) {
214  std::cout << "time: " << m_now
215  << ", first stage log unnorm weight: " << logFirstStageUnNormWeights[ii]
216  << "\n";
217  }
218  #endif
219 
220 
221  }
222 
223  // draw ks (indexes) (handles underflow issues)
224  arrayUInt myKs = m_kGen.sample(logFirstStageUnNormWeights);
225 
226  // now draw xts
227  float_t m1(-std::numeric_limits<float_t>::infinity());
228  float_t first_cll_sum(0.0);
229  float_t second_cll_sum(0.0);
230  float_t third_cll_sum(0.0);
231  ssv xtm1k;
232  ssv muT;
233  for(size_t ii = 0; ii < nparts; ++ii)
234  {
235  // calclations for log p(y_t|y_{1:t-1}) (using log-sum-exp trick)
236  second_cll_sum += std::exp( logFirstStageUnNormWeights[ii] - m2 );
237  third_cll_sum += std::exp( m_logUnNormWeights[ii] - m3 );
238 
239  // sampling and unnormalized weight update
240  xtm1k = oldPartics[myKs[ii]];
241  m_particles[ii] = fSamp(xtm1k);
242  muT = propMu(xtm1k);
243  m_logUnNormWeights[ii] += logGEv(data, m_particles[ii]) - logGEv(data, muT);
244 
245 
246  #ifndef DROPPINGTHISINRPACKAGE
247  if constexpr(debug){
248  std::cout << "time: " << m_now
249  << ", transposed sample: " << m_particles[ii].transpose()
250  << ", log unnorm weight: " << m_logUnNormWeights[ii] << "\n";
251  }
252  #endif
253 
254  // update m1
255  if(m_logUnNormWeights[ii] > m1)
256  m1 = m_logUnNormWeights[ii];
257  }
258 
259  // calculate estimate for log of last conditonal likelihood
260  for(size_t p = 0; p < nparts; ++p)
261  first_cll_sum += std::exp( m_logUnNormWeights[p] - m1 );
262  m_logLastCondLike = m1 + std::log(first_cll_sum) + m2 + std::log(second_cll_sum) - 2*m3 - 2*std::log(third_cll_sum);
263 
264  #ifndef DROPPINGTHISINRPACKAGE
265  if constexpr(debug)
266  std::cout << "time: " << m_now << ", log cond like: " << m_logLastCondLike << "\n";
267  #endif
268 
269  // calculate expectations before you resample
270  unsigned int fId(0);
271  for(auto & h : fs){
272 
273  Mat testOutput = h(m_particles[0]);
274  unsigned int rows = testOutput.rows();
275  unsigned int cols = testOutput.cols();
276  Mat numer = Mat::Zero(rows,cols);
277  float_t denom(0.0);
278 
279  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){
280  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - m1);
281  denom += std::exp(m_logUnNormWeights[prtcl] - m1);
282  }
283  m_expectations[fId] = numer/denom;
284 
285  #ifndef DROPPINGTHISINRPACKAGE
286  if constexpr(debug)
287  std::cout << "transposed expectation " << fId << "; " << m_expectations[fId] << "\n";
288  #endif
289 
290  fId++;
291  }
292 
293  // if you have to resample
294  if( (m_now+1)%m_rs == 0)
295  m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
296 
297  // advance time
298  m_now += 1;
299 
300  } else { // (m_now == 0)
301 
302  float_t max(-std::numeric_limits<float_t>::infinity());
303  for(size_t ii = 0; ii < nparts; ++ii)
304  {
305  // sample particles
306  m_particles[ii] = q1Samp(data);
307  m_logUnNormWeights[ii] = logMuEv(m_particles[ii]);
308  m_logUnNormWeights[ii] += logGEv(data, m_particles[ii]);
309  m_logUnNormWeights[ii] -= logQ1Ev(m_particles[ii], data);
310 
311  // print stuff if debug mode is on
312  #ifndef DROPPINGTHISINRPACKAGE
313  if constexpr(debug) {
314  std::cout << "time: " << m_now
315  << ", log unnorm weight: " << m_logUnNormWeights[ii]
316  << ", transposed sample: " << m_particles[ii].transpose()
317  << "\n";
318  }
319  #endif
320 
321  // update maximum
322  if( m_logUnNormWeights[ii] > max)
323  max = m_logUnNormWeights[ii];
324  }
325 
326  // calculate log-likelihood with log-exp-sum trick
327  float_t sumExp(0.0);
328  for( size_t i = 0; i < nparts; ++i){
329  sumExp += std::exp( m_logUnNormWeights[i] - max );
330  }
331  m_logLastCondLike = - std::log( static_cast<float_t>(nparts) ) + max + std::log(sumExp);
332 
333  // calculate expectations before you resample
334  m_expectations.resize(fs.size());
335  unsigned int fId(0);
336  for(auto & h : fs){
337 
338  Mat testOutput = h(m_particles[0]);
339  unsigned int rows = testOutput.rows();
340  unsigned int cols = testOutput.cols();
341  Mat numer = Mat::Zero(rows,cols);
342  float_t denom(0.0);
343  for(size_t prtcl = 0; prtcl < nparts; ++prtcl){
344  numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - max);
345  denom += std::exp(m_logUnNormWeights[prtcl] - max);
346  }
347  m_expectations[fId] = numer/denom;
348 
349  #ifndef DROPPINGTHISINRPACKAGE
350  if constexpr(debug)
351  std::cout << "transposed expectation " << fId << "; " << m_expectations[fId] << "\n";
352  #endif
353 
354  fId++;
355  }
356 
357  // resample if you should (automatically normalizes)
358  if( (m_now+1) % m_rs == 0)
359  m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
360 
361  // advance time step
362  m_now += 1;
363  }
364 
365 }
366 
367 
368 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
370 {
371  return m_logLastCondLike;
372 }
373 
374 
375 template<size_t nparts, size_t dimx, size_t dimy, typename resamp_t, typename float_t, bool debug>
377 {
378  return m_expectations;
379 }
380 
381 } // namespace filters
382 } // namespace pf
383 #endif //APF_H
384 
pf::bases::pf_base
Definition: pf_base.h:35
pf::filters::APF::m_logUnNormWeights
std::array< float_t, nparts > m_logUnNormWeights
particle unnormalized weights
Definition: auxiliary_pf.h:147
pf::filters::APF::m_rs
unsigned int m_rs
the resampling schedule
Definition: auxiliary_pf.h:156
pf::filters::APF::filter
void filter(const osv &data, const std::vector< std::function< const Mat(const ssv &)> > &fs=std::vector< std::function< const Mat(const ssv &)> >())
Use a new datapoint to update the filtering distribution (or smoothing if pathLength > 0).
Definition: auxiliary_pf.h:186
pf::filters::APF::propMu
virtual ssv propMu(const ssv &xtm1)=0
Evaluates the proposal distribution taking a Eigen::Matrix<float_t,dimx,1> from the previous time's s...
pf::filters::APF::getExpectations
std::vector< Mat > getExpectations() const
return all stored expectations (taken with respect to $p(x_t|y_{1:t})$
Definition: auxiliary_pf.h:376
pf::filters::APF::ssv
Eigen::Matrix< float_t, dimx, 1 > ssv
Definition: auxiliary_pf.h:44
pf_base.h
All non Rao-Blackwellized particle filters without covariates inherit from this.
pf::filters::APF::APF
APF(const unsigned int &rs=1)
The constructor.
Definition: auxiliary_pf.h:172
pf::filters::APF::arrayUInt
std::array< unsigned int, nparts > arrayUInt
Definition: auxiliary_pf.h:54
pf::filters::APF::logQ1Ev
virtual float_t logQ1Ev(const ssv &x1, const osv &y1)=0
Evaluates the log of q1.
pf::filters::APF::q1Samp
virtual ssv q1Samp(const osv &y1)=0
Samples from q1.
rv_samp.h
all rv samplers must inherit from this.
pf::filters::APF::m_expectations
std::vector< Mat > m_expectations
expectations E[h(x_t) | y_{1:t}] for user defined "h"s
Definition: auxiliary_pf.h:165
pf::filters::APF::~APF
virtual ~APF()
The (virtual) destructor.
Definition: auxiliary_pf.h:182
pf::filters::APF::osv
Eigen::Matrix< float_t, dimy, 1 > osv
Definition: auxiliary_pf.h:46
pf::filters::APF::m_now
unsigned int m_now
curren time
Definition: auxiliary_pf.h:150
pf::rvsamp::k_gen< nparts, float_t >
pf::filters::APF
A base-class for Auxiliary Particle Filtering. Filtering only, no smoothing.
Definition: auxiliary_pf.h:39
pf::filters::APF::m_kGen
rvsamp::k_gen< nparts, float_t > m_kGen
k generator object (default ctor'd)
Definition: auxiliary_pf.h:162
pf::filters::APF::fSamp
virtual ssv fSamp(const ssv &xtm1)=0
Samples from f.
pf::filters::APF::m_particles
std::array< ssv, nparts > m_particles
particle samples
Definition: auxiliary_pf.h:144
pf::filters::APF::m_logLastCondLike
float_t m_logLastCondLike
log p(y_t|y_{1:t-1}) or log p(y1)
Definition: auxiliary_pf.h:153
pf::filters::APF::getLogCondLike
float_t getLogCondLike() const
Get the latest log conditional likelihood.
Definition: auxiliary_pf.h:369
pf::filters::APF::logGEv
virtual float_t logGEv(const osv &yt, const ssv &xt)=0
Evaluates the log of g.
pf::filters::APF::m_resampler
resamp_t m_resampler
resampler object (default ctor'd)
Definition: auxiliary_pf.h:159
pf::filters::APF::arrayfloat_t
std::array< float_t, nparts > arrayfloat_t
Definition: auxiliary_pf.h:50
pf::filters::APF::Mat
Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic > Mat
Definition: auxiliary_pf.h:48
pf::filters::APF::logMuEv
virtual float_t logMuEv(const ssv &x1)=0
Evaluates the log of mu.
pf::filters::APF::arrayVec
std::array< ssv, nparts > arrayVec
Definition: auxiliary_pf.h:52