Bayesian Filtering Library  Generated from SVN r
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
KalmanFilter Class Referenceabstract

Class representing the family of all Kalman Filters (EKF, IEKF, ...) More...

#include <kalmanfilter.h>

Inheritance diagram for KalmanFilter:
Filter< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > ExtendedKalmanFilter IteratedExtendedKalmanFilter NonminimalKalmanFilter SRIteratedExtendedKalmanFilter

Public Member Functions

 KalmanFilter (Gaussian *prior)
 Constructor. More...
 
virtual ~KalmanFilter ()
 Destructor.
 
virtual GaussianPostGet ()
 Get Posterior density. More...
 
void AllocateMeasModel (const vector< unsigned int > &meas_dimensions)
 Function to allocate memory needed during the measurement update,.
 
void AllocateMeasModel (const unsigned int &meas_dimensions)
 Function to allocate memory needed during the measurement update.
 
virtual void Reset (Pdf< MatrixWrapper::ColumnVector > *prior)
 Reset Filter.
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u, MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)
 Full Update (system with inputs/sensing params) More...
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)
 Full Update (system without inputs, with sensing params) More...
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z)
 Full Update (system without inputs/sensing params) More...
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u, MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z)
 Full Update (system with inputs, without sensing params) More...
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u)
 System Update (system with inputs) More...
 
virtual bool Update (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel)
 System Update (system without inputs) More...
 
virtual bool Update (MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)
 Measurement Update (system with "sensing params") More...
 
virtual bool Update (MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z)
 Measurement Update (system without "sensing params") More...
 
int TimeStepGet () const
 Get current time. More...
 

Protected Member Functions

void PostSigmaSet (const MatrixWrapper::SymmetricMatrix &s)
 Set covariance of posterior estimate.
 
void PostMuSet (const MatrixWrapper::ColumnVector &c)
 Set expected value of posterior estimate.
 
void CalculateSysUpdate (const MatrixWrapper::ColumnVector &J, const MatrixWrapper::Matrix &F, const MatrixWrapper::SymmetricMatrix &Q)
 
void CalculateMeasUpdate (const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &Z, const MatrixWrapper::Matrix &H, const MatrixWrapper::SymmetricMatrix &R)
 
virtual void SysUpdate (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u)=0
 System Update. More...
 
virtual void MeasUpdate (MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)=0
 Measurement Update (overloaded) More...
 
virtual bool UpdateInternal (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u, MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)
 Actual implementation of Update, varies along filters. More...
 
void CalculateNis (const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &Z, const MatrixWrapper::Matrix &H, const MatrixWrapper::SymmetricMatrix &R)
 

Protected Attributes

ColumnVector _Mu_new
 
SymmetricMatrix _Sigma_new
 
Matrix _Sigma_temp
 
Matrix _Sigma_temp_par
 
Matrix _SMatrix
 
Matrix _K
 
double _Nis
 
std::map< unsigned int, MeasUpdateVariables > _mapMeasUpdateVariables
 
std::map< unsigned int, MeasUpdateVariables >::iterator _mapMeasUpdateVariables_it
 
Pdf< MatrixWrapper::ColumnVector > * _prior
 prior Pdf
 
Pdf< MatrixWrapper::ColumnVector > * _post
 Pointer to the Posterior Pdf. More...
 
int _timestep
 Represents the current timestep of the filter. More...
 

Friends

class NonminimalKalmanFilter
 

Detailed Description

Class representing the family of all Kalman Filters (EKF, IEKF, ...)

This is a class representing the family of all Kalman Filter (KF). Kalman filters are filters in which the Posterior density is represented by a Gaussian density. Kalman filters are only applicable to continuous systems.

The system of updating the Posterior density is implemented in this base class. However, the parameters used for this update differ for different KFs (Simple KF,EKF,IEKF): that's why the xUpdate members are still pure virtual functions.

This class is the base class for all sorts of KFs.

See also
Gaussian
LinearAnalyticSystemModelGaussianUncertainty

Definition at line 49 of file kalmanfilter.h.

Constructor & Destructor Documentation

§ KalmanFilter()

KalmanFilter ( Gaussian prior)

Constructor.

Precondition
you created the prior
Parameters
priorpointer to the Gaussian Pdf prior density

Member Function Documentation

§ CalculateMeasUpdate()

void CalculateMeasUpdate ( const MatrixWrapper::ColumnVector z,
const MatrixWrapper::ColumnVector Z,
const MatrixWrapper::Matrix H,
const MatrixWrapper::SymmetricMatrix R 
)
protected

Calculate Kalman filter Measurement Update

\[ x_k = x_{k-} + K.(z - Z) \]

\[ P_k = (I-K.H).P_{k-} \]

with

\[ K = P_{k-}.H'.(H.P_{k-}.H'+R)^{-1} \]

§ CalculateSysUpdate()

void CalculateSysUpdate ( const MatrixWrapper::ColumnVector J,
const MatrixWrapper::Matrix F,
const MatrixWrapper::SymmetricMatrix Q 
)
protected

Calculate Kalman filter System Update

\[ x_k = J \]

\[ P_k = F.P_{k-}.F' + Q \]

§ MeasUpdate()

virtual void MeasUpdate ( MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const  measmodel,
const MatrixWrapper::ColumnVector z,
const MatrixWrapper::ColumnVector s 
)
protectedpure virtual

Measurement Update (overloaded)

Update the filter's Posterior density using the sensor measurements, an input and the measurement model. This method is used when the measurements depend on the inputs too (doesn't happen very often, does it?) BEWARE: the first time the measurment update is called with a new size of measurement, new allocations are done

Parameters
measmodelpointer to the measurement model the filter should use
zsensor measurement
sinput to the system (must be of the same type as u for now, since this was not yet implemented in ConditionalPdf

Implemented in ExtendedKalmanFilter, SRIteratedExtendedKalmanFilter, and IteratedExtendedKalmanFilter.

§ PostGet()

virtual Gaussian* PostGet ( )
virtual

Get Posterior density.

Get the current Posterior density

Returns
a pointer to the current posterior

Reimplemented from Filter< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector >.

§ SysUpdate()

virtual void SysUpdate ( SystemModel< MatrixWrapper::ColumnVector > *const  sysmodel,
const MatrixWrapper::ColumnVector u 
)
protectedpure virtual

System Update.

Update the filter's Posterior density using the deterministic inputs to the system and the system model

Parameters
sysmodelpointer to the system model the filter should use
uinput to the system

Implemented in ExtendedKalmanFilter, SRIteratedExtendedKalmanFilter, and IteratedExtendedKalmanFilter.

§ TimeStepGet()

int TimeStepGet ( ) const
inherited

Get current time.

Get the current time of the filter

Returns
the current timestep

Definition at line 50 of file filter.h.

§ Update() [1/8]

bool Update ( SystemModel< SVar > *const  sysmodel,
const SVar &  u,
MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z,
const SVar &  s 
)
virtualinherited

Full Update (system with inputs/sensing params)

Parameters
sysmodelpointer to the system model to use for update
uinput to the system
measmodelpointer to the measurement model to use for update
zmeasurement
s"sensing parameter"

Definition at line 56 of file filter.h.

§ Update() [2/8]

bool Update ( SystemModel< SVar > *const  sysmodel,
MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z,
const SVar &  s 
)
virtualinherited

Full Update (system without inputs, with sensing params)

Parameters
sysmodelpointer to the system model to use for update
measmodelpointer to the measurement model to use for update
zmeasurement
s"sensing parameter"

Definition at line 66 of file filter.h.

§ Update() [3/8]

bool Update ( SystemModel< SVar > *const  sysmodel,
MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z 
)
virtualinherited

Full Update (system without inputs/sensing params)

Parameters
sysmodelpointer to the system model to use for update
measmodelpointer to the measurement model to use for update
zmeasurement

Definition at line 76 of file filter.h.

§ Update() [4/8]

bool Update ( SystemModel< SVar > *const  sysmodel,
const SVar &  u,
MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z 
)
virtualinherited

Full Update (system with inputs, without sensing params)

Parameters
sysmodelpointer to the system model to use for update
uinput to the system
measmodelpointer to the measurement model to use for update
zmeasurement

Definition at line 85 of file filter.h.

§ Update() [5/8]

bool Update ( SystemModel< SVar > *const  sysmodel,
const SVar &  u 
)
virtualinherited

System Update (system with inputs)

Parameters
sysmodelpointer to the system model to use for update
uinput to the system

Definition at line 95 of file filter.h.

§ Update() [6/8]

bool Update ( SystemModel< SVar > *const  sysmodel)
virtualinherited

System Update (system without inputs)

Parameters
sysmodelpointer to the system model to use for update

Definition at line 103 of file filter.h.

§ Update() [7/8]

bool Update ( MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z,
const SVar &  s 
)
virtualinherited

Measurement Update (system with "sensing params")

Parameters
measmodelpointer to the measurement model to use for update
zmeasurement
s"sensing parameter"

Definition at line 110 of file filter.h.

§ Update() [8/8]

bool Update ( MeasurementModel< MVar, SVar > *const  measmodel,
const MVar &  z 
)
virtualinherited

Measurement Update (system without "sensing params")

Parameters
measmodelpointer to the measurement model to use for update
zmeasurement

Definition at line 119 of file filter.h.

§ UpdateInternal()

virtual bool UpdateInternal ( SystemModel< MatrixWrapper::ColumnVector > *const  sysmodel,
const MatrixWrapper::ColumnVector u,
MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const  measmodel,
const MatrixWrapper::ColumnVector z,
const MatrixWrapper::ColumnVector s 
)
protectedvirtual

Actual implementation of Update, varies along filters.

Parameters
sysmodelpointer to the used system model
uinput param for proposal density
measmodelpointer to the used measurementmodel
zmeasurement param for proposal density
ssensor param for proposal density

Implements Filter< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector >.

Friends And Related Function Documentation

§ NonminimalKalmanFilter

friend class NonminimalKalmanFilter
friend

Very dirty hack to avoid ugly methods PostSigmaSet and PostMuSet to be public! NonMinimalKalmanFilter should be redesigned though!

Definition at line 111 of file kalmanfilter.h.

Member Data Documentation

§ _post

Pdf<MatrixWrapper::ColumnVector >* _post
protectedinherited

Pointer to the Posterior Pdf.

The Posterior Pdf represents the subjective belief of the person applying the filter AFTER processing inputs and measurements. A filter does not maintain the beliefs at all timesteps t, since this leads to non-constant (or ever growing if you prefer) memory requirements. However, it is possible, to copy the Posterior density at all timesteps in your application by means of the PostGet() member function

See also
PostGet()

Definition at line 95 of file filter.h.

§ _timestep

int _timestep
protectedinherited

Represents the current timestep of the filter.

Todo:
Check wether this really belongs here

Definition at line 100 of file filter.h.


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