ssmkit  master-68aed98
Particle< Process, Resampler > Class Template Reference

Particle Filter. More...

Public Types

using CompeleteState = std::tuple< arma::mat, arma::vec >
 Type of the state posterior. More...
 

Public Member Functions

 Particle (Process process, Resampler resampler, unsigned long particles_num)
 Constructor. More...
 
template<class... Args>
void predict (const Args &...args)
 Prediction. More...
 
template<class Measurement , class... TArgs>
CompeleteState correct (const Measurement &measurement, const TArgs &...args)
 Correction. More...
 
CompeleteState initialize ()
 Initialization. More...
 
const arma::vec & getWeights (void) const
 
const arma::mat & getStateParticles (void) const
 
- Public Member Functions inherited from RecursiveBayesianBase< Particle< Process, Resampler > >
decltype(auto) step (const TMeasurement &measurement, const TArgs &...args)
 
decltype(auto) filter (const std::vector< TMeasurement > &measurements, const std::vector< TControl > &controls=std::vector< std::tuple<>>())
 

Detailed Description

template<class Process, class Resampler>
class ssmkit::filter::Particle< Process, Resampler >

Particle Filter.

Definition at line 30 of file particle.hpp.

Member Typedef Documentation

using CompeleteState = std::tuple<arma::mat, arma::vec>

Type of the state posterior.

\( \{\mathbf{x}^{(i)}_t,\omega^{(i)}\}_{i=1}^{M}\)

Definition at line 36 of file particle.hpp.

Constructor & Destructor Documentation

Particle ( Process  process,
Resampler  resampler,
unsigned long  particles_num 
)
inline

Constructor.

returns a Particle filter object.

Parameters
processThe process model object that the PF is defined for
resamplerThe resampling algorithm
particles_numNumber of particles \(M\)

Definition at line 63 of file particle.hpp.

Member Function Documentation

void predict ( const Args &...  args)
inline

Prediction.

Performs prediction step.

\begin{equation}\mathbf{x}^{(i)}_t \sim p(\mathbf{x}_t| \tilde{\mathbf{x}}^{(i)}_{t-1}, y^d_1, \cdots, y^d_{N_d}) \quad \text{for} \quad i=1,\cdots,M\end{equation}

Parameters
args...Control variables \(y^d_1, \cdots, y^d_{N_d}\) of the dynamic process, if any.

Definition at line 80 of file particle.hpp.

CompeleteState correct ( const Measurement &  measurement,
const TArgs &...  args 
)
inline

Correction.

Performs correction step.

\begin{equation} \omega^{(i)} = \tilde{\omega}^{(i)} p(\mathbf{z}_t| \mathbf{x}^{(i)}_t, y^m_1, \cdots, y^m_{N_m}) \end{equation}

\begin{equation}\{\mathbf{x}^{(i)}_t,\omega^{(i)}\}_{i=1}^{M} \overset{\mbox{resample}}{\longrightarrow} \{\tilde{\mathbf{x}}^{(i)}_t,\tilde{\omega}^{(i)}\}_{i=1}^{M}\end{equation}

Parameters
measurementMeasurement vector \(\mathbf{z}_t\).
args...Control variables \(y^m_1, \cdots, y^m_{N_m}\) of the measurement process, if any.
Returns
Estimated state \(\{\tilde{\mathbf{x}}^{(i)}_t,\tilde{\omega}^{(i)}\}_{i=1}^{M}\)

Definition at line 99 of file particle.hpp.

CompeleteState initialize ( )
inline

Initialization.

Returns
Estimated state \(\{\tilde{\mathbf{x}}^{(i)}_0,\tilde{\omega}^{(i)}\}_{i=1}^{M}\)

Definition at line 117 of file particle.hpp.

const arma::vec& getWeights ( void  ) const
inline
Returns
Estimated state \(\{\tilde{\omega}^{(i)}\}_{i=1}^{M}\)

Definition at line 133 of file particle.hpp.

const arma::mat& getStateParticles ( void  ) const
inline
Returns
Estimated state \(\{\tilde{\mathbf{x}}^{(i)}_t\}_{i=1}^{M}\)

Definition at line 135 of file particle.hpp.


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