pf
pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug > Class Template Referenceabstract

Rao-Blackwellized/Marginal Particle Filter with inner Kalman Filter objectss. More...

#include <rbpf.h>

Inheritance diagram for pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >:
Collaboration diagram for pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >:

Public Types

using sssv = Eigen::Matrix< float_t, dimss, 1 >
 
using nsssv = Eigen::Matrix< float_t, dimnss, 1 >
 
using osv = Eigen::Matrix< float_t, dimy, 1 >
 
using Mat = Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic >
 
using nsssMat = Eigen::Matrix< float_t, dimnss, dimnss >
 
using cfModType = kalman< dimnss, dimy, 0, float_t, debug >
 
using arrayMod = std::array< cfModType, nparts >
 
using arrayVec = std::array< sssv, nparts >
 
using arrayfloat_t = std::array< float_t, nparts >
 
- Public Types inherited from pf::bases::rbpf_base< float_t, dimss, dimnss, dimy >
using obs_sized_vec = Eigen::Matrix< float_t, dimobs, 1 >
 
using sampled_state_sized_vec = Eigen::Matrix< float_t, dim_s_state, 1 >
 
using not_sampled_state_sized_vec = Eigen::Matrix< float_t, dim_ns_state, 1 >
 
using dynamic_matrix = Eigen::Matrix< float_t, Eigen::Dynamic, Eigen::Dynamic >
 
using func = std::function< const dynamic_matrix(const not_sampled_state_sized_vec &, const sampled_state_sized_vec &)>
 
using func_vec = std::vector< func >
 

Public Member Functions

 rbpf_kalman (const unsigned int &resamp_sched=1)
 The constructor. More...
 
void filter (const osv &data, const std::vector< std::function< const Mat(const nsssv &x1t, const sssv &x2t)> > &fs=std::vector< std::function< const Mat(const nsssv &x1t, const sssv &x2t)> >())
 Filter! More...
 
float_t getLogCondLike () const
 Get the latest log conditional likelihood. More...
 
std::vector< MatgetExpectations () const
 Get the latest filtered expectation E[h(x_1t, x_2t) | y_{1:t}]. More...
 
virtual float_t logMuEv (const sssv &x21)=0
 Evaluates the first time state density. More...
 
virtual sssv q1Samp (const osv &y1)=0
 Sample from the first time's proposal distribution. More...
 
virtual nsssv initKalmanMean (const sssv &x21)=0
 Provides the initial mean vector for each Kalman filter object. More...
 
virtual nsssMat initKalmanVar (const sssv &x21)=0
 Provides the initial covariance matrix for each Kalman filter object. More...
 
virtual sssv qSamp (const sssv &x2tm1, const osv &yt)=0
 Samples the time t second component. More...
 
virtual float_t logQ1Ev (const sssv &x21, const osv &y1)=0
 Evaluates the proposal density of the second state component at time 1. More...
 
virtual float_t logFEv (const sssv &x2t, const sssv &x2tm1)=0
 Evaluates the state transition density for the second state component. More...
 
virtual float_t logQEv (const sssv &x2t, const sssv &x2tm1, const osv &yt)=0
 Evaluates the proposal density at time t > 1. More...
 
virtual void updateKalman (cfModType &kMod, const osv &yt, const sssv &x2t)=0
 How to update your inner Kalman filter object at each time. More...
 
- Public Member Functions inherited from pf::bases::rbpf_base< float_t, dimss, dimnss, dimy >
virtual void filter (const obs_sized_vec &data, const func_vec &fs=func_vec())=0
 the filtering function that must be defined More...
 
virtual ~rbpf_base ()
 virtual destructor
 

Private Attributes

unsigned int m_rs
 
arrayMod m_p_innerMods
 
arrayVec m_p_samps
 
arrayfloat_t m_logUnNormWeights
 
unsigned int m_now
 
float_t m_lastLogCondLike
 
resamp_t m_resampler
 
std::vector< Matm_expectations
 

Additional Inherited Members

- Static Public Attributes inherited from pf::bases::rbpf_base< float_t, dimss, dimnss, dimy >
static constexpr unsigned int dim_sampled_state
 
static constexpr unsigned int dim_not_sampled_state
 
static constexpr unsigned int dim_obs
 

Detailed Description

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t, typename float_t, bool debug = false>
class pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >

Rao-Blackwellized/Marginal Particle Filter with inner Kalman Filter objectss.

Author
t

Member Typedef Documentation

◆ arrayfloat_t

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::arrayfloat_t = std::array<float_t,nparts>

array of weights

◆ arrayMod

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::arrayMod = std::array<cfModType,nparts>

array of model objects

◆ arrayVec

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::arrayVec = std::array<sssv,nparts>

array of samples

◆ cfModType

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::cfModType = kalman<dimnss,dimy,0,float_t,debug>

closed-form model type

◆ Mat

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::Mat = Eigen::Matrix<float_t,Eigen::Dynamic,Eigen::Dynamic>

dynamic size matrices

◆ nsssMat

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::nsssMat = Eigen::Matrix<float_t,dimnss,dimnss>

"not sampled state size matrix"

◆ nsssv

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::nsssv = Eigen::Matrix<float_t,dimnss,1>

"not sampled state size vector"

◆ osv

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::osv = Eigen::Matrix<float_t,dimy,1>

"observation size vector"

◆ sssv

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
using pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::sssv = Eigen::Matrix<float_t,dimss,1>

"sampled state size vector"

Constructor & Destructor Documentation

◆ rbpf_kalman()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug>
pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::rbpf_kalman ( const unsigned int &  resamp_sched = 1)

The constructor.

Parameters
resamp_schedhow often you want to resample (e.g once every resamp_sched time points)

Member Function Documentation

◆ filter()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug>
void pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::filter ( const osv data,
const std::vector< std::function< const Mat(const nsssv &x1t, const sssv &x2t)> > &  fs = std::vector<std::function<const Mat(const nsssv &x1t, const sssv &x2t)> >() 
)

Filter!

The workhorse function

Parameters
datathe most recent observable portion of the time series.
fsa vector of functions computing E[h(x_1t, x_2t^i)| x_2t^i,y_1:t]. to be averaged to yield E[h(x_1t, x_2t)|,y_1:t]

◆ getExpectations()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug>
auto pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::getExpectations

Get the latest filtered expectation E[h(x_1t, x_2t) | y_{1:t}].

Get the expectations you're keeping track of.

Returns
a vector of Mats

◆ getLogCondLike()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug>
float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::getLogCondLike

Get the latest log conditional likelihood.

Returns
the latest log conditional likelihood.

◆ initKalmanMean()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual nsssv pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::initKalmanMean ( const sssv x21)
pure virtual

Provides the initial mean vector for each Kalman filter object.

provides the initial mean vector for each Kalman filter object.

Parameters
x21the second state componenent at time 1.
Returns
a nsssv representing the unconditional mean.

◆ initKalmanVar()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual nsssMat pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::initKalmanVar ( const sssv x21)
pure virtual

Provides the initial covariance matrix for each Kalman filter object.

provides the initial covariance matrix for each Kalman filter object.

Parameters
x21the second state component at time 1.
Returns
a covariance matrix.

◆ logFEv()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::logFEv ( const sssv x2t,
const sssv x2tm1 
)
pure virtual

Evaluates the state transition density for the second state component.

Evaluates the state transition density for the second state component.

Parameters
x2tthe current second state component.
x2tm1the previous second state component.
Returns
a float_t evaluation.

◆ logMuEv()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::logMuEv ( const sssv x21)
pure virtual

Evaluates the first time state density.

evaluates log mu(x21).

Parameters
x21component two at time 1
Returns
a float_t evaluation

◆ logQ1Ev()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::logQ1Ev ( const sssv x21,
const osv y1 
)
pure virtual

Evaluates the proposal density of the second state component at time 1.

Evaluates the proposal density of the second state component at time 1.

Parameters
x21the second state component at time 1 you sampled.
y1time 1 observation.
Returns
a float_t evaluation of the density.

◆ logQEv()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::logQEv ( const sssv x2t,
const sssv x2tm1,
const osv yt 
)
pure virtual

Evaluates the proposal density at time t > 1.

Evaluates the proposal density at time t > 1.

Parameters
x2tthe current second state component.
x2tm1the previous second state component.
ytthe current time series observation.
Returns
a float_t evaluation.

◆ q1Samp()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual sssv pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::q1Samp ( const osv y1)
pure virtual

Sample from the first time's proposal distribution.

samples the second component of the state at time 1.

Parameters
y1most recent datum.
Returns
a Vec sample for x21.

◆ qSamp()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual sssv pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::qSamp ( const sssv x2tm1,
const osv yt 
)
pure virtual

Samples the time t second component.

Samples the time t second component.

Parameters
x2tm1the previous time's second state component.
ytthe current observation.
Returns
a sssv sample of the second state component at the current time.

◆ updateKalman()

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
virtual void pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::updateKalman ( cfModType kMod,
const osv yt,
const sssv x2t 
)
pure virtual

How to update your inner Kalman filter object at each time.

How to update your inner Kalman filter object at each time.

Parameters
kModa Kalman filter object describing the conditional closed-form model.
ytthe current time series observation.
x2tthe current second state component.

Member Data Documentation

◆ m_expectations

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
std::vector<Mat> pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_expectations
private

expectations

◆ m_lastLogCondLike

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
float_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_lastLogCondLike
private

log p(y_t|y_{1:t-1}) or log p(y1)

◆ m_logUnNormWeights

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
arrayfloat_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_logUnNormWeights
private

the array of the (log of) unnormalized weights

◆ m_now

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
unsigned int pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_now
private

the current time period

◆ m_p_innerMods

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
arrayMod pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_p_innerMods
private

the array of inner Kalman filter objects

◆ m_p_samps

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
arrayVec pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_p_samps
private

the array of particle samples

◆ m_resampler

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
resamp_t pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_resampler
private

resampler object

◆ m_rs

template<size_t nparts, size_t dimnss, size_t dimss, size_t dimy, typename resamp_t , typename float_t , bool debug = false>
unsigned int pf::filters::rbpf_kalman< nparts, dimnss, dimss, dimy, resamp_t, float_t, debug >::m_rs
private

the resamplign schedule


The documentation for this class was generated from the following file: