#include <KalmanFilter.h>
Public Types | |
| enum | Factorization { QR, CHOLESKY } |
| Specifies the factorization variant to use. More... | |
| typedef GaussianDensity::shared_ptr | State |
| The Kalman filter state, represented as a shared pointer to a GaussianDensity. More... | |
Public Member Functions | |
| size_t | dim () const |
| State | init (const Vector &x0, const Matrix &P0) const |
| State | init (const Vector &x0, const SharedDiagonal &P0) const |
| KalmanFilter (size_t n, Factorization method=KALMANFILTER_DEFAULT_FACTORIZATION) | |
| State | predict (const State &p, const Matrix &F, const Matrix &B, const Vector &u, const SharedDiagonal &modelQ) const |
| State | predict2 (const State &p, const Matrix &A0, const Matrix &A1, const Vector &b, const SharedDiagonal &model=nullptr) const |
| State | predictQ (const State &p, const Matrix &F, const Matrix &B, const Vector &u, const Matrix &Q) const |
| void | print (const std::string &s="") const |
| State | update (const State &p, const Matrix &H, const Vector &z, const SharedDiagonal &model) const |
| State | updateQ (const State &p, const Matrix &H, const Vector &z, const Matrix &R) const |
Static Public Member Functions | |
| static Key | step (const State &p) |
Private Member Functions | |
| State | fuse (const State &p, GaussianFactor::shared_ptr newFactor) const |
| State | solve (const GaussianFactorGraph &factorGraph) const |
Private Attributes | |
| const GaussianFactorGraph::Eliminate | function_ |
| Elimination algorithm used. More... | |
| const Matrix | I_ |
Identity matrix of size . More... | |
| const size_t | n_ |
| Dimensionality of the state. More... | |
Kalman Filter class
Maintains a Gaussian density under linear-Gaussian motion and measurement models using the square-root information form.
The filter is functional; it does not maintain internal state. Instead:
Definition at line 42 of file KalmanFilter.h.
The Kalman filter state, represented as a shared pointer to a GaussianDensity.
Definition at line 55 of file KalmanFilter.h.
Specifies the factorization variant to use.
| Enumerator | |
|---|---|
| QR | |
| CHOLESKY | |
Definition at line 48 of file KalmanFilter.h.
|
inline |
Constructor.
| n | Dimensionality of the state. |
| method | Factorization method (default: QR unless compile-flag set). |
Definition at line 84 of file KalmanFilter.h.
|
inline |
Return the dimensionality of the state.
Definition at line 214 of file KalmanFilter.h.
|
private |
Fuse two states.
| p | The prior state. |
| newFactor | The new factor to incorporate. |
Definition at line 70 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::init | ( | const Vector & | x0, |
| const Matrix & | P0 | ||
| ) | const |
Create the initial state with a full covariance matrix.
| x0 | Initial state estimate. |
| P0 | Full covariance matrix. |
Definition at line 92 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::init | ( | const Vector & | x0, |
| const SharedDiagonal & | P0 | ||
| ) | const |
Create the initial state (prior density at time
).
In Kalman Filter notation:
: Initial state estimate.
: Initial covariance matrix.| x0 | Estimate of the state at time 0 ( ). |
| P0 | Covariance matrix ( ), given as a diagonal Gaussian model. |
Definition at line 82 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::predict | ( | const State & | p, |
| const Matrix & | F, | ||
| const Matrix & | B, | ||
| const Vector & | u, | ||
| const SharedDiagonal & | modelQ | ||
| ) | const |
Predict the next state
.
In Kalman Filter notation:
: Predicted state.
: Predicted covariance.Motion model:
where
is zero-mean Gaussian noise with covariance
.
| p | Previous state ( ). |
| F | State transition matrix ( ). |
| B | Control input matrix ( ). |
| u | Control vector ( ). |
| modelQ | Noise model ( , diagonal Gaussian). |
). Definition at line 113 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::predict2 | ( | const State & | p, |
| const Matrix & | A0, | ||
| const Matrix & | A1, | ||
| const Vector & | b, | ||
| const SharedDiagonal & | model = nullptr |
||
| ) | const |
Predict the next state using a GaussianFactor motion model.
| p | Previous state. |
| A0 | Factor matrix. |
| A1 | Factor matrix. |
| b | Constant term vector. |
| model | Noise model (optional). |
Definition at line 149 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::predictQ | ( | const State & | p, |
| const Matrix & | F, | ||
| const Matrix & | B, | ||
| const Vector & | u, | ||
| const Matrix & | Q | ||
| ) | const |
Predict the next state with a full covariance matrix.
| p | Previous state. |
| F | State transition matrix. |
| B | Control input matrix. |
| u | Control vector. |
| Q | Full covariance matrix ( ). |
Definition at line 125 of file KalmanFilter.cpp.
| void gtsam::KalmanFilter::print | ( | const std::string & | s = "" | ) | const |
Print the Kalman filter details.
| s | Optional string prefix. |
Definition at line 108 of file KalmanFilter.cpp.
|
private |
Solve the factor graph.
| factorGraph | The Gaussian factor graph to solve. |
Definition at line 56 of file KalmanFilter.cpp.
Return the step index
(starts at 0, incremented at each predict step).
| p | The current state. |
Definition at line 126 of file KalmanFilter.h.
| KalmanFilter::State gtsam::KalmanFilter::update | ( | const State & | p, |
| const Matrix & | H, | ||
| const Vector & | z, | ||
| const SharedDiagonal & | model | ||
| ) | const |
Update the Kalman filter with a measurement.
Observation model:
where
is zero-mean Gaussian noise with covariance R. In this version, R is restricted to diagonal Gaussians (model parameter)
| p | Previous state. |
| H | Observation matrix. |
| z | Measurement vector. |
| model | Noise model (diagonal Gaussian). |
Definition at line 159 of file KalmanFilter.cpp.
| KalmanFilter::State gtsam::KalmanFilter::updateQ | ( | const State & | p, |
| const Matrix & | H, | ||
| const Vector & | z, | ||
| const Matrix & | R | ||
| ) | const |
Update the Kalman filter with a measurement using a full covariance matrix.
| p | Previous state. |
| H | Observation matrix. |
| z | Measurement vector. |
| R | Full covariance matrix. |
Definition at line 170 of file KalmanFilter.cpp.
|
private |
Elimination algorithm used.
Definition at line 61 of file KalmanFilter.h.
|
private |
Identity matrix of size
.
Definition at line 59 of file KalmanFilter.h.
|
private |
Dimensionality of the state.
Definition at line 58 of file KalmanFilter.h.