Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Friends
BFL::KalmanFilter Class Reference

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

#include <kalmanfilter.h>

Inheritance diagram for BFL::KalmanFilter:
Inheritance graph
[legend]

List of all members.

Classes

struct  MeasUpdateVariables

Public Member Functions

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.
 KalmanFilter (Gaussian *prior)
 Constructor.
virtual GaussianPostGet ()
 Get Posterior density.
virtual ~KalmanFilter ()
 Destructor.

Protected Member Functions

void CalculateMeasUpdate (const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &Z, const MatrixWrapper::Matrix &H, const MatrixWrapper::SymmetricMatrix &R)
void CalculateSysUpdate (const MatrixWrapper::ColumnVector &J, const MatrixWrapper::Matrix &F, const MatrixWrapper::SymmetricMatrix &Q)
virtual void MeasUpdate (MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const measmodel, const MatrixWrapper::ColumnVector &z, const MatrixWrapper::ColumnVector &s)=0
 Measurement Update (overloaded)
void PostMuSet (const MatrixWrapper::ColumnVector &c)
 Set expected value of posterior estimate.
void PostSigmaSet (const MatrixWrapper::SymmetricMatrix &s)
 Set covariance of posterior estimate.
virtual void SysUpdate (SystemModel< MatrixWrapper::ColumnVector > *const sysmodel, const MatrixWrapper::ColumnVector &u)=0
 System Update.
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.

Protected Attributes

Matrix _K
std::map< unsigned int,
MeasUpdateVariables
_mapMeasUpdateVariables
std::map< unsigned int,
MeasUpdateVariables >
::iterator 
_mapMeasUpdateVariables_it
ColumnVector _Mu_new
SymmetricMatrix _Sigma_new
Matrix _Sigma_temp
Matrix _Sigma_temp_par
Matrix _SMatrix

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

Constructor.

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

Definition at line 28 of file kalmanfilter.cpp.

Destructor.

Definition at line 39 of file kalmanfilter.cpp.


Member Function Documentation

void BFL::KalmanFilter::AllocateMeasModel ( const vector< unsigned int > &  meas_dimensions)

Function to allocate memory needed during the measurement update,.

Definition at line 45 of file kalmanfilter.cpp.

void BFL::KalmanFilter::AllocateMeasModel ( const unsigned int &  meas_dimensions)

Function to allocate memory needed during the measurement update.

Definition at line 63 of file kalmanfilter.cpp.

void BFL::KalmanFilter::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} \]

Definition at line 88 of file kalmanfilter.cpp.

void BFL::KalmanFilter::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 \]

Definition at line 76 of file kalmanfilter.cpp.

virtual void BFL::KalmanFilter::MeasUpdate ( MeasurementModel< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector > *const  measmodel,
const MatrixWrapper::ColumnVector &  z,
const MatrixWrapper::ColumnVector &  s 
) [protected, pure 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 BFL::ExtendedKalmanFilter, BFL::SRIteratedExtendedKalmanFilter, and BFL::IteratedExtendedKalmanFilter.

Get Posterior density.

Get the current Posterior density

Returns:
a pointer to the current posterior

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

Definition at line 160 of file kalmanfilter.cpp.

void BFL::KalmanFilter::PostMuSet ( const MatrixWrapper::ColumnVector &  c) [protected]

Set expected value of posterior estimate.

Definition at line 153 of file kalmanfilter.cpp.

void BFL::KalmanFilter::PostSigmaSet ( const MatrixWrapper::SymmetricMatrix &  s) [protected]

Set covariance of posterior estimate.

Definition at line 147 of file kalmanfilter.cpp.

virtual void BFL::KalmanFilter::SysUpdate ( SystemModel< MatrixWrapper::ColumnVector > *const  sysmodel,
const MatrixWrapper::ColumnVector &  u 
) [protected, pure 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 BFL::ExtendedKalmanFilter, BFL::SRIteratedExtendedKalmanFilter, and BFL::IteratedExtendedKalmanFilter.

bool BFL::KalmanFilter::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 
) [protected, virtual]

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 BFL::Filter< MatrixWrapper::ColumnVector, MatrixWrapper::ColumnVector >.

Definition at line 130 of file kalmanfilter.cpp.


Friends And Related Function Documentation

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 110 of file kalmanfilter.h.


Member Data Documentation

Matrix BFL::KalmanFilter::_K [protected]

Definition at line 101 of file kalmanfilter.h.

Definition at line 102 of file kalmanfilter.h.

std::map<unsigned int, MeasUpdateVariables>::iterator BFL::KalmanFilter::_mapMeasUpdateVariables_it [protected]

Definition at line 103 of file kalmanfilter.h.

ColumnVector BFL::KalmanFilter::_Mu_new [protected]

Definition at line 96 of file kalmanfilter.h.

SymmetricMatrix BFL::KalmanFilter::_Sigma_new [protected]

Definition at line 97 of file kalmanfilter.h.

Matrix BFL::KalmanFilter::_Sigma_temp [protected]

Definition at line 98 of file kalmanfilter.h.

Definition at line 99 of file kalmanfilter.h.

Matrix BFL::KalmanFilter::_SMatrix [protected]

Definition at line 100 of file kalmanfilter.h.


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


bfl
Author(s): Klaas Gadeyne, Wim Meeussen, Tinne Delaet and many others. See web page for a full contributor list. ROS package maintained by Wim Meeussen.
autogenerated on Sun Oct 5 2014 22:29:53