6 #ifdef DROPPINGTHISINRPACKAGE
10 #include <Eigen/Dense>
33 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug=false>
39 using ssv = Eigen::Matrix<float_t, dimx, 1>;
41 using osv = Eigen::Matrix<float_t, dimy, 1>;
43 using Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic>;
84 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;
127 virtual float_t
logFEv (
const ssv &xt,
const ssv &xtm1 ) = 0;
146 virtual float_t
logQEv (
const ssv &xt,
const ssv &xtm1,
const osv &yt ) = 0;
180 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug>
183 , m_logLastCondLike(0.0)
190 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug>
194 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug>
197 return m_logLastCondLike;
201 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug>
204 return m_expectations;
208 template<
size_t nparts,
size_t dimx,
size_t dimy,
typename resamp_t,
typename float_t,
bool debug>
218 float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
219 for(
size_t ii = 0; ii < nparts; ++ii)
223 if (m_logUnNormWeights[ii] > maxOldLogUnNormWts)
224 maxOldLogUnNormWts = m_logUnNormWeights[ii];
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);
233 m_particles[ii] = newSamp;
235 #ifndef DROPPINGTHISINRPACKAGE
237 std::cout <<
"time: " << m_now <<
", transposed sample: " << m_particles[ii].transpose() <<
", log unnorm weight: " << m_logUnNormWeights[ii] <<
"\n";
243 float_t maxNumer = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
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);
250 m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
254 float_t weightNormConst(0.0);
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);
263 for(
size_t prtcl = 0; prtcl < nparts; ++prtcl){
264 numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - maxNumer );
265 denom += std::exp(m_logUnNormWeights[prtcl] - maxNumer);
267 m_expectations[fId] = numer/denom;
270 #ifndef DROPPINGTHISINRPACKAGE
272 std::cout <<
"transposed expectation " << fId <<
": " << m_expectations[fId].transpose() <<
"\n";
279 if( (m_now + 1) % m_resampSched == 0)
280 m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
290 for(
size_t ii = 0; ii < nparts; ++ii)
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);
298 #ifndef DROPPINGTHISINRPACKAGE
300 std::cout <<
"time: " << m_now <<
", transposed sample: " << m_particles[ii].transpose() <<
", log unnorm weight: " << m_logUnNormWeights[ii] <<
"\n";
307 float_t max = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
309 for(
size_t i = 0; i < nparts; ++i){
310 sumExp += std::exp(m_logUnNormWeights[i] - max);
312 m_logLastCondLike = -std::log(nparts) + max + std::log(sumExp);
315 m_expectations.resize(fs.size());
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);
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);
329 m_expectations[fId] = numer/denom;
332 #ifndef DROPPINGTHISINRPACKAGE
334 std::cout <<
"transposed expectation " << fId <<
": " << m_expectations[fId].transpose() <<
"\n";
341 if( (m_now + 1) % m_resampSched == 0)
342 m_resampler.resampLogWts(m_particles, m_logUnNormWeights);
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>
372 using ssv = Eigen::Matrix<float_t, dimx, 1>;
374 using osv = Eigen::Matrix<float_t, dimy, 1>;
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>;
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&)> >());
460 virtual float_t
logGEv (
const osv &yt,
const ssv &xt ) = 0;
469 virtual float_t
logFEv (
const ssv &xt,
const ssv &xtm1 ) = 0;
489 virtual float_t
logQEv (
const ssv &xt,
const ssv &xtm1,
const osv &yt ) = 0;
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>
521 , m_logLastCondLike(0.0)
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>
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>
535 return m_logLastCondLike;
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>
542 return m_expectations;
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)
556 float_t maxOldLogUnNormWts(-std::numeric_limits<float_t>::infinity());
557 for(
size_t ii = 0; ii < nparts; ++ii)
561 if (m_logUnNormWeights[ii] > maxOldLogUnNormWts)
562 maxOldLogUnNormWts = m_logUnNormWeights[ii];
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);
571 m_particles[ii] = newSamp;
573 #ifndef DROPPINGTHISINRPACKAGE
575 std::cout <<
"time: " << m_now <<
", transposed sample: " << m_particles[ii].transpose() <<
", log unnorm weight: " << m_logUnNormWeights[ii] <<
"\n";
580 float_t maxNumer = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
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);
587 m_logLastCondLike = maxNumer + std::log(sumExp1) - maxOldLogUnNormWts - std::log(sumExp2);
591 float_t weightNormConst(0.0);
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);
600 for(
size_t prtcl = 0; prtcl < nparts; ++prtcl){
601 numer += h(m_particles[prtcl]) * std::exp(m_logUnNormWeights[prtcl] - maxNumer );
602 denom += std::exp(m_logUnNormWeights[prtcl] - maxNumer);
604 m_expectations[fId] = numer/denom;
607 #ifndef DROPPINGTHISINRPACKAGE
609 std::cout <<
"transposed expectation " << fId <<
": " << m_expectations[fId].transpose() <<
"\n";
616 if( (m_now + 1) % m_resampSched == 0)
617 m_resampler.resampLogWts(m_particles, m_logUnNormWeights, Uresamp);
627 for(
size_t ii = 0; ii < nparts; ++ii)
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);
635 #ifndef DROPPINGTHISINRPACKAGE
637 std::cout <<
"time: " << m_now <<
", transposed sample: " << m_particles[ii].transpose() <<
", log unnorm weight: " << m_logUnNormWeights[ii] <<
"\n";
643 float_t max = *std::max_element(m_logUnNormWeights.begin(), m_logUnNormWeights.end());
645 for(
size_t i = 0; i < nparts; ++i){
646 sumExp += std::exp(m_logUnNormWeights[i] - max);
648 m_logLastCondLike = -std::log(nparts) + max + std::log(sumExp);
651 m_expectations.resize(fs.size());
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);
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);
665 m_expectations[fId] = numer/denom;
668 #ifndef DROPPINGTHISINRPACKAGE
670 std::cout <<
"transposed expectation " << fId <<
": " << m_expectations[fId].transpose() <<
"\n";
677 if( (m_now + 1) % m_resampSched == 0)
678 m_resampler.resampLogWts(m_particles, m_logUnNormWeights, Uresamp);
689 #endif //SISR_FILTER_H