ssmkit  master-68aed98
Kalman< STA_MAP, OBS_MAP > Class Template Reference

Kalman filter. More...

Public Types

using TProcess = Hierarchical< Markov< Gaussian, STA_MAP, Gaussian >, Memoryless< Gaussian, OBS_MAP >>
 Type of process object. More...
 
using TCompeleteState = std::tuple< arma::vec, arma::mat >
 Type of the posterior state \((\hat{\mathbf{x}}, \hat{\mathbf{P}})\). More...
 

Public Member Functions

 Kalman (const TProcess &process)
 Construct a Kalman filter. More...
 
template<class... TArgs>
void predict (const TArgs &...args)
 Prediction. More...
 
template<class... TArgs>
TCompeleteState correct (const arma::vec &measurement, const TArgs &...args)
 Correction. More...
 
TCompeleteState initialize ()
 Initialization. More...
 
- Public Member Functions inherited from RecursiveBayesianBase< Kalman< STA_MAP, OBS_MAP > >
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 STA_MAP, class OBS_MAP>
class ssmkit::filter::Kalman< STA_MAP, OBS_MAP >

Kalman filter.

Definition at line 30 of file kalman.hpp.

Member Typedef Documentation

using TProcess = Hierarchical<Markov<Gaussian, STA_MAP, Gaussian>, Memoryless<Gaussian, OBS_MAP>>

Type of process object.

Definition at line 37 of file kalman.hpp.

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

Type of the posterior state \((\hat{\mathbf{x}}, \hat{\mathbf{P}})\).

Definition at line 40 of file kalman.hpp.

Constructor & Destructor Documentation

Kalman ( const TProcess process)
inline

Construct a Kalman filter.

Construct a Kalman filter with parameters taken from process argument.

Definition at line 67 of file kalman.hpp.

Member Function Documentation

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

Prediction.

Performs the prediction step.

\begin{equation}\hat{\mathbf{x}}_{t|t-1} = \mathbf{F}\hat{\mathbf{x}}_{t-1|t-1}+u_d(y^d_1, \cdots, y^d_{N_d}) \end{equation}

\begin{equation}\mathbf{P}_{t|t-1} = \mathbf{F}\mathbf{P}_{t-1|t-1}\mathbf{F}^T+\mathbf{Q}\end{equation}

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

Definition at line 94 of file kalman.hpp.

TCompeleteState correct ( const arma::vec &  measurement,
const TArgs &...  args 
)
inline

Correction.

Performs correction step.

\begin{equation}\tilde{\mathbf{z}}_t=\mathbf{z}_t-\mathbf{H}\hat{\mathbf{x}}_{t|t-1}-u_m(y^m_1, \cdots, y^m_{N_m})\end{equation}

\begin{equation}\mathbf {S}_t=\mathbf{H}\mathbf{P}_{t|t-1}\mathbf{H}^T+\mathbf{R} \end{equation}

\begin{equation}\mathbf{K}_t=\mathbf{P}_{t|t-1}\mathbf{H}^T\mathbf{S}_t^{-1} \end{equation}

\begin{equation}\hat{\mathbf{x}}_{t|t}=\hat{\mathbf{x}}_{t|t-1}+\mathbf{K}_{t}\tilde{\mathbf{z}}_t \end{equation}

\begin{equation}\mathbf{P}_{t|t}=(I-\mathbf{K}_t\mathbf{H})\mathbf{P}_{t|t-1}\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 \((\hat{\mathbf{x}}_{t|t}, \mathbf{P}_{t|t})\)

Definition at line 118 of file kalman.hpp.

TCompeleteState initialize ( )
inline

Initialization.

Returns
Initial state \((\hat{\mathbf{x}}_{0|0}, \mathbf{P}_{0|0})\)

Definition at line 137 of file kalman.hpp.


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