state-observation 1.7.0
Loading...
Searching...
No Matches
kalman-filter-base.hpp
Go to the documentation of this file.
1
22#ifndef KALMANFILTERBASEHPP
23#define KALMANFILTERBASEHPP
24
25#include <state-observation/api.h>
28
29namespace stateObservation
30{
31
51{
52public:
54 typedef Matrix Amatrix;
55
57 typedef Matrix Cmatrix;
58
60 typedef Matrix Qmatrix;
61
63 typedef Matrix Rmatrix;
64
66 typedef Matrix Pmatrix;
67
68 typedef Eigen::LLT<Pmatrix> LLTPMatrix;
69
72
75
78
84
100
102 virtual void setA(const Amatrix & A);
103
104 virtual Matrix getA() const;
105
107 virtual void clearA();
108
110 virtual void setC(const Cmatrix & C);
111
112 virtual Matrix getC() const;
113
115 virtual void clearC();
116
118 virtual void setR(const Rmatrix & R);
119 inline void setMeasurementCovariance(const Rmatrix & R)
120 {
121 setR(R);
122 }
123
125 virtual Matrix getR() const;
127 {
128 return getR();
129 }
130
132 virtual void clearR();
133
135 virtual void setQ(const Qmatrix & Q);
136 inline void setProcessCovariance(const Qmatrix & Q)
137 {
138 setQ(Q);
139 }
140
142 virtual Matrix getQ() const;
144 {
145 return getQ();
146 }
147
149 virtual void clearQ();
150
152 virtual void setStateCovariance(const Pmatrix & P);
153
156 virtual void clearStates();
157
160 virtual void clearStateCovariance();
161
164
166 virtual void reset();
167
170
173
176
179
181 bool checkAmatrix(const Amatrix &) const;
182
185
188
191
193 bool checkCmatrix(const Cmatrix &) const;
194
197
200
203
206
208 bool checkQmatrix(const Qmatrix &) const;
209
212
215
218
221
223 bool checkRmatrix(const Rmatrix &) const;
224
227
230
233
236
238 bool checkPmatrix(const Pmatrix &) const;
239
242
245
248
250 virtual bool checkStateTangentVector(const StateVectorTan & v) const;
251
254
257
260
262 virtual bool checkMeasureTangentVector(const MeasureVectorTan &) const;
263
267 virtual void setStateSize(Index n);
268
274 virtual void setStateSize(Index n, Index nt);
275
279 virtual void setMeasureSize(Index m);
280
286 virtual void setMeasureSize(Index m, Index mt);
287
290
292 virtual const StateVector & getInnovation();
293
298 inline const StateVector & updateStatePrediction();
299
304 inline MeasureVector updateStateAndMeasurementPrediction();
305
308
311
314
316 const Matrix & getLastGain() const;
317
321
322protected:
325
328
331
334
337
340
343
346
349
352
355
358
361
364
367
381
383
384public:
386};
387
389{
390 prediction_(this->x_.getTime() + 1);
391 return xbar_();
392}
393
400
401} // namespace stateObservation
402
403#endif // KALMANFILTERBASEHPP
TimeIndex getTime() const
Get the time index.
It mostly implements the equations of Kalman filtering It is suitablle by derivation to be used incas...
Definition kalman-filter-base.hpp:51
Pmatrix getPmatrixZero() const
Get a matrix having the size of the P matrix having zero values.
const StateVector & getLastPrediction() const
get the last predicted state
const MeasureVector & getLastMeasurement() const
get the last measurement
virtual StateVector prediction_(TimeIndex k)=0
The abstract method to overload to implement f(x,u)
Cmatrix getCmatrixZero() const
Get a matrix having the size of the C matrix having zero values.
Index mt_
the size of tangent space of the measurement space
Definition kalman-filter-base.hpp:327
KalmanFilterBase(Index n, Index m, Index p=0)
virtual MeasureVector predictSensor_(TimeIndex k)
Predicts the sensor measurement,.
const Matrix & getLastGain() const
get the last Kalman gain matrix
Matrix c_
Containers for the jacobian matrix of the measurement.
Definition kalman-filter-base.hpp:348
Rmatrix getRmatrixIdentity() const
Get an identity matrix having the size of the R matrix.
virtual void clearQ()
Clear the process noise covariance matrix.
virtual void clearC()
Clear the jacobian dh/dx.
virtual void setC(const Cmatrix &C)
Set the value of the Jacobian dh/dx.
Matrix Amatrix
The type of the jacobian df/dx.
Definition kalman-filter-base.hpp:54
Matrix getMeasurementCovariance() const
Definition kalman-filter-base.hpp:126
bool checkQmatrix(const Qmatrix &) const
checks whether or not a matrix has the dimensions of the Q matrix
Rmatrix getRmatrixConstant(double c) const
Get a matrix having the size of the R matrix having "c" values.
virtual void clearA()
Clear the jacobian df/dx.
Matrix Pmatrix
The type of the covariance matrix of the state estimation error.
Definition kalman-filter-base.hpp:66
Cmatrix getCmatrixRandom() const
Get a matrix having the size of the C matrix having random values.
bool checkRmatrix(const Rmatrix &) const
checks whether or not a matrix has the dimensions of the R matrix
void setMeasurementCovariance(const Rmatrix &R)
Definition kalman-filter-base.hpp:119
virtual MeasureVector simulateSensor_(const StateVector &x, TimeIndex k)=0
The abstract method to overload to implement h(x,u)
KalmanFilterBase()
Default constructor.
virtual void reset()
Resets all the observer.
KalmanFilterBase(Index n, Index nt, Index m, Index mt, Index p)
virtual StateVector oneStepEstimation_()
The Kalman filter loop.
virtual void setQ(const Qmatrix &Q)
Set the process noise covariance matrix.
Qmatrix getQmatrixRandom() const
Get a matrix having the size of the Q matrix having random values.
Matrix pr_
Container for the covariance matrix of the estimation error.
Definition kalman-filter-base.hpp:357
virtual MeasureVectorTan measureTangentVectorZero() const
Gives a vector of measurement tangent vector size having zero values.
virtual MeasureVector getSimulatedMeasurement(TimeIndex k)
Get simulation of the measurement y_k using the state estimation.
const StateVector & updateStatePrediction()
Definition kalman-filter-base.hpp:388
Qmatrix getQmatrixConstant(double c) const
Get a matrix having the size of the Q matrix having "c" values.
virtual void setStateSize(Index n)
Matrix r_
Container for the measurement noise covariance matrice.
Definition kalman-filter-base.hpp:354
Vector innovation_
Vector containing the inovation of the Kalman filter.
Definition kalman-filter-base.hpp:366
virtual void setMeasureSize(Index m, Index mt)
Amatrix getAmatrixRandom() const
Get a matrix having the size of the A matrix having random values.
virtual Pmatrix getStateCovariance() const
Get the covariance matrix of the current time state estimation.
StateVectorArithmetics * arithm_
Definition kalman-filter-base.hpp:382
virtual Matrix getR() const
Set the measurement noise covariance matrix.
Matrix getProcessCovariance() const
Definition kalman-filter-base.hpp:143
IndexedVector xbar_
container for the prediction
Definition kalman-filter-base.hpp:360
void setStateArithmetics(StateVectorArithmetics *arith)
virtual MeasureVectorTan measureTangentVectorConstant(double c) const
Gives a vector of measurement tangent vector size having duplicated "c" value.
Qmatrix getQmatrixIdentity() const
Get an identity matrix having the size of the Q matrix.
Pmatrix getPmatrixConstant(double c) const
Get a matrix having the size of the P matrix having "c" values.
Vector MeasureVectorTan
MeasureVector is the type of measurement tanegnt vector.
Definition kalman-filter-base.hpp:74
Pmatrix getPmatrixRandom() const
Get a matrix having the size of the P matrix having random values.
Pmatrix getPmatrixIdentity() const
Get an identity matrix having the size of the P matrix.
virtual void setStateSize(Index n, Index nt)
virtual Matrix getC() const
Rmatrix getRmatrixRandom() const
Get a matrix having the size of the R matrix having random values.
Qmatrix getQmatrixZero() const
Get a matrix having the size of the Q matrix having zero values.
Matrix Kmatrix
The type of Kalman gain matrix.
Definition kalman-filter-base.hpp:330
IndexedVector ybar_
container for the prediction of the sensor
Definition kalman-filter-base.hpp:363
Matrix Qmatrix
The type of the covariance matrix of the process noise v.
Definition kalman-filter-base.hpp:60
Amatrix getAmatrixIdentity() const
Get an identity matrix having the size of the A matrix.
virtual bool checkMeasureTangentVector(const MeasureVectorTan &) const
Tells whether or not the vector has the dimensions of a measurement tangent vector.
Index nt_
the size of tangent space of the state space
Definition kalman-filter-base.hpp:324
void setProcessCovariance(const Qmatrix &Q)
Definition kalman-filter-base.hpp:136
Matrix a_
Containers for the jacobian matrix of the process.
Definition kalman-filter-base.hpp:345
Rmatrix getRmatrixZero() const
Get a matrix having the size of the R matrix having zero values.
Matrix q_
Container for the process noise covariance matrice.
Definition kalman-filter-base.hpp:351
virtual MeasureVectorTan measureTangentVectorRandom() const
Gives a vector of measurement tangent vector size having random values.
const MeasureVector & getLastPredictedMeasurement() const
get the last predicted measurement
bool checkPmatrix(const Pmatrix &) const
Checks whether or not a matrix has the dimensions of the P matrix.
virtual Matrix getQ() const
Set the process noise covariance matrix.
virtual void setMeasureSize(Index m)
bool checkCmatrix(const Cmatrix &) const
checks whether or not a matrix has the dimensions of the C matrix
virtual bool checkStateTangentVector(const StateVectorTan &v) const
Tells whether or not the vector has the dimensions of a state tangent vector.
Amatrix getAmatrixConstant(double c) const
Get a matrix having the size of the A matrix having "c" values.
Vector StateVectorTan
StateVector is the type of state tangent vector.
Definition kalman-filter-base.hpp:71
Matrix Rmatrix
The type of the covariance matrix of the measurement noise w.
Definition kalman-filter-base.hpp:63
virtual StateVectorTan stateTangentVectorRandom() const
Gives a vector of state tangent vector size having random values.
virtual void clearStates()
clears the state and the state covariance
virtual const StateVector & getInnovation()
Get the last vector of innovation of the Kalman filter.
virtual StateVectorTan stateTangentVectorConstant(double c) const
Gives a vector of state tangent vector size having duplicated "c" value.
bool checkAmatrix(const Amatrix &) const
checks whether or not a matrix has the dimensions of the A matrix
MeasureVector updateStateAndMeasurementPrediction()
Definition kalman-filter-base.hpp:394
virtual void clearR()
Clear the measurement noise covariance matrix.
Cmatrix getCmatrixConstant(double c) const
Get a matrix having the size of the C matrix having "c" values.
virtual void setR(const Rmatrix &R)
Set the measurement noise covariance matrix.
virtual StateVectorTan stateTangentVectorZero() const
Gives a vector of state tangent vector size having zero values.
Matrix Cmatrix
The type of the jacobian dh/dx.
Definition kalman-filter-base.hpp:57
virtual void setStateCovariance(const Pmatrix &P)
Set the covariance matrix of the current time state estimation error.
Eigen::LLT< Pmatrix > LLTPMatrix
Definition kalman-filter-base.hpp:68
virtual Matrix getA() const
Amatrix getAmatrixZero() const
Get a matrix having the size of the A matrix having zero values.
virtual void setA(const Amatrix &A)
Set the value of the jacobian df/dx.
Vector StateVector
StateVector is the type of state vector.
Definition observer-base.hpp:42
Vector MeasureVector
MeasureVector is the type of measurements vector.
Definition observer-base.hpp:45
This class is used to customize the way the difference between measurements, the state update functio...
Definition state-vector-arithmetics.hpp:28
Defines the base class of online zero delay observers. Zero delay observers are the classical state o...
Definition zero-delay-observer.hpp:44
IndexedVector x_
while the measurements and iputs are put in lists
Definition zero-delay-observer.hpp:194
Definition bidim-elastic-inv-pendulum-dyn-sys.hpp:21
long int TimeIndex
Definition definitions.hpp:139
Eigen::MatrixXd Matrix
Dynamic sized Matrix.
Definition definitions.hpp:100
Eigen::Index Index
Definition definitions.hpp:138
Eigen::VectorXd Vector
Dynamic sized scalar vector.
Definition definitions.hpp:76
LLTPMatrix inoMeasCovLLT
Definition kalman-filter-base.hpp:376
Matrix inoMeasCovInverse
Definition kalman-filter-base.hpp:375
Vector xhat
Definition kalman-filter-base.hpp:372
Matrix kGain
Definition kalman-filter-base.hpp:377
Matrix t
Definition kalman-filter-base.hpp:378
Matrix inoMeasCov
Definition kalman-filter-base.hpp:374
EIGEN_MAKE_ALIGNED_OPERATOR_NEW Matrix pbar
Definition kalman-filter-base.hpp:371
Matrix mKc
Definition kalman-filter-base.hpp:379
Vector inoMeas
Definition kalman-filter-base.hpp:373
Defines the base class of online zero delay observers. Zero delay observers are the classical state o...