state-observation 1.7.0
Loading...
Searching...
No Matches
lipm-dcm-estimator.hpp
Go to the documentation of this file.
1
12#ifndef LIPMDCMBIASESTIMATOR_HPP
13#define LIPMDCMBIASESTIMATOR_HPP
14
15#include <state-observation/api.h>
19
20namespace stateObservation
21{
22
70
71class STATE_OBSERVATION_DLLAPI LipmDcmEstimator
72{
73private:
74 constexpr static double defaultDt_ = 0.005;
75 constexpr static double defaultOmega_ = tools::sqrt(cst::gravityConstant);
76
77public:
79 constexpr static double defaultBiasDriftSecond = 0.002;
80
82 constexpr static double defaultZmpErrorStd = 0.005;
83 constexpr static double defaultDcmErrorStd = 0.01;
84
86 constexpr static double defaultDCMUncertainty = 0.01;
87 constexpr static double defaultBiasUncertainty = 0.01;
88
90 constexpr static double defaultBiasLimit = -1;
91
108 LipmDcmEstimator(double dt = defaultDt_,
109 double omega_0 = defaultOmega_,
110 double biasDriftPerSecondStd = defaultBiasDriftSecond,
111 double dcmMeasureErrorStd = defaultDcmErrorStd,
112 double zmpMeasureErrorStd = defaultZmpErrorStd,
113 const Vector2 & biasLimit = Vector2::Constant(defaultBiasLimit),
114 const Vector2 & initZMP = Vector2::Zero(),
115 const Vector2 & initDcm = Vector2::Zero(),
116 const Vector2 & initBias = Vector2::Zero(),
117 const Vector2 & initDcmUncertainty = Vector2::Constant(defaultDCMUncertainty),
118 const Vector2 & initBiasUncertainty = Vector2::Constant(defaultBiasUncertainty));
119
131 void resetWithMeasurements(const Vector2 & measuredDcm,
132 const Vector2 & measuredZMP,
133 const Matrix2 & yaw = Matrix2::Identity(),
134 bool measurementIsWithBias = true,
135 const Vector2 & initBias = Vector2::Constant(0),
136 const Vector2 & initBiasuncertainty = Vector2::Constant(defaultBiasUncertainty));
137
149 inline void resetWithMeasurements(const Vector2 & measuredDcm,
150 const Vector2 & measuredZMP,
151 double yaw,
152 bool measurementIsWithBias = true,
153 const Vector2 & initBias = Vector2::Constant(0),
154 const Vector2 & initBiasuncertainty = Vector2::Constant(defaultBiasUncertainty))
155 {
156 resetWithMeasurements(measuredDcm, measuredZMP, Rotation2D(yaw).toRotationMatrix(), measurementIsWithBias, initBias,
157 initBiasuncertainty);
158 }
159
172 inline void resetWithMeasurements(const Vector2 & measuredDcm,
173 const Vector2 & measuredZMP,
174 const Matrix3 & rotation,
175 bool measurementIsWithBias = true,
176 const Vector2 & initBias = Vector2::Constant(0),
177 const Vector2 & initBiasuncertainty = Vector2::Constant(defaultBiasUncertainty))
178 {
179 resetWithMeasurements(measuredDcm, measuredZMP, kine::rotationMatrixToYawAxisAgnostic(rotation),
180 measurementIsWithBias, initBias, initBiasuncertainty);
181 }
182
185
189 void setLipmNaturalFrequency(double omega_0);
190
194 inline double getLipmNaturalFrequency() const
195 {
196 return omega0_;
197 }
198
205 void setUnbiasedCoMOffset(const Vector2 & gamma);
206
212 {
213 return gamma_;
214 }
215
218 void setZMPCoef(double kappa);
219
222 double getZMPCoef() const
223 {
224 return kappa_;
225 }
226
230 void setSamplingTime(double dt);
231
235 void setBias(const Vector2 & bias);
236
240 void setBias(const Vector2 & bias, const Vector2 & uncertainty);
241
245 void setBiasDriftPerSecond(double driftPerSecond);
246
251 void setBiasLimit(const Vector2 & biasLimit);
252
256 void setUnbiasedDCM(const Vector2 & dcm);
257
261 void setUnbiasedDCM(const Vector2 & dcm, const Vector2 & uncertainty);
262
266
270
281 inline void setInputs(const Vector2 & dcm,
282 const Vector2 & zmp,
283 const Matrix3 & orientation,
284 const Vector2 & CoMOffset_gamma = Vector2::Zero(),
285 const double ZMPCoef_kappa = 1)
286 {
287 setInputs(dcm, zmp, kine::rotationMatrixToYawAxisAgnostic(orientation), CoMOffset_gamma, ZMPCoef_kappa);
288 }
289
298 inline void setInputs(const Vector2 & dcm,
299 const Vector2 & zmp,
300 double yaw,
301 const Vector2 & CoMOffset_gamma = Vector2::Zero(),
302 const double ZMPCoef_kappa = 1)
303 {
304 setInputs(dcm, zmp, Rotation2D(yaw).toRotationMatrix(), CoMOffset_gamma, ZMPCoef_kappa);
305 }
306
314 void setInputs(const Vector2 & dcm,
315 const Vector2 & zmp,
316 const Matrix2 & R = Matrix2::Identity(),
317 const Vector2 & CoMOffset_gamma = Vector2::Zero(),
318 const double ZMPCoef_kappa = 1);
319
323 void update();
324
330
335
339 inline Vector2 getLocalBias() const
340 {
341 return previousOrientation_.transpose() * getBias();
342 }
343
348 {
349 return filter_;
350 }
351
354 inline const LinearKalmanFilter & getFilter() const
355 {
356 return filter_;
357 }
358
359protected:
360 typedef Eigen::Matrix<double, 4, 2> Matrix42;
361 typedef Eigen::Matrix<double, 2, 4> Matrix24;
362
365
366 double omega0_;
367 Vector2 gamma_ = Vector2::Zero();
368 double kappa_ = 1;
369 double dt_;
372
373 bool needUpdateMatrices_ = true;
374
376
378
388
390
392 inline static Matrix2 Vec2ToSqDiag_(const Vector2 & v)
393 {
394 return Vector2(v.array().square()).asDiagonal();
395 }
396
398 inline static Matrix2 dblToDiag_(const double & d)
399 {
400 return Vector2::Constant(d).asDiagonal();
401 }
402
404 inline static Matrix2 dblToSqDiag_(const double & d)
405 {
406 return dblToDiag_(d * d);
407 }
408
409public:
411};
412
413} // namespace stateObservation
414
415#endif
The class of a Linear Kalman filter.
Definition linear-kalman-filter.hpp:48
Filtering of divergent component of motion (DCM) and estimation of a bias betweeen the DCM and the co...
Definition lipm-dcm-estimator.hpp:72
void setZMPCoef(double kappa)
Set the ZMP oefficient.
void setBias(const Vector2 &bias)
Set the Bias object from a guess.
void setZmpMeasureErrorStd(double)
Set the Zmp Measurement Error Stamdard devbiation.
void setLipmNaturalFrequency(double omega_0)
Set the Lipm Natural Frequency.
double getZMPCoef() const
get the ZMP Coefficient
Definition lipm-dcm-estimator.hpp:222
~LipmDcmEstimator()
Destroy the Lipm Dcm Bias Estimator object.
Eigen::Matrix< double, 2, 4 > Matrix24
Definition lipm-dcm-estimator.hpp:361
Matrix2 previousOrientation_
Definition lipm-dcm-estimator.hpp:389
static Matrix2 dblToSqDiag_(const double &d)
builds a constant 2x2 diagonal from a square of a double
Definition lipm-dcm-estimator.hpp:404
double zmpErrorStd_
Definition lipm-dcm-estimator.hpp:371
Matrix24 C_
this needs to be transposed
Definition lipm-dcm-estimator.hpp:383
void setBiasLimit(const Vector2 &biasLimit)
Set the Bias Limit.
Matrix4 A_
Definition lipm-dcm-estimator.hpp:380
double getLipmNaturalFrequency() const
Get the Lipm Natural Frequency.
Definition lipm-dcm-estimator.hpp:194
Matrix4 Q_
process noise
Definition lipm-dcm-estimator.hpp:387
Vector2 biasLimit_
Definition lipm-dcm-estimator.hpp:377
Vector2 getLocalBias() const
Get the estimated Bias expressed in the local frame of the robot.
Definition lipm-dcm-estimator.hpp:339
LipmDcmEstimator(double dt=defaultDt_, double omega_0=defaultOmega_, double biasDriftPerSecondStd=defaultBiasDriftSecond, double dcmMeasureErrorStd=defaultDcmErrorStd, double zmpMeasureErrorStd=defaultZmpErrorStd, const Vector2 &biasLimit=Vector2::Constant(defaultBiasLimit), const Vector2 &initZMP=Vector2::Zero(), const Vector2 &initDcm=Vector2::Zero(), const Vector2 &initBias=Vector2::Zero(), const Vector2 &initDcmUncertainty=Vector2::Constant(defaultDCMUncertainty), const Vector2 &initBiasUncertainty=Vector2::Constant(defaultBiasUncertainty))
Construct a new Lipm Dcm Estimator object.
static Matrix2 dblToDiag_(const double &d)
builds a constant 2x2 diagonal from a double
Definition lipm-dcm-estimator.hpp:398
double dt_
Definition lipm-dcm-estimator.hpp:369
Vector2 getUnbiasedDCM() const
Get the Unbiased DCM filtered by the estimator.
void setUnbiasedDCM(const Vector2 &dcm)
set the unbiased DCM position from a guess
LinearKalmanFilter filter_
Definition lipm-dcm-estimator.hpp:379
void setUnbiasedCoMOffset(const Vector2 &gamma)
Set an offset on the CoM dynamics.
void update()
Runs the estimation. Needs to be called every timestep.
void setInputs(const Vector2 &dcm, const Vector2 &zmp, const Matrix3 &orientation, const Vector2 &CoMOffset_gamma=Vector2::Zero(), const double ZMPCoef_kappa=1)
Set the Inputs of the estimator.
Definition lipm-dcm-estimator.hpp:281
Eigen::Matrix< double, 4, 2 > Matrix42
Definition lipm-dcm-estimator.hpp:360
static Matrix2 Vec2ToSqDiag_(const Vector2 &v)
builds a diagonal out of the square valued of the Vec2
Definition lipm-dcm-estimator.hpp:392
void setBias(const Vector2 &bias, const Vector2 &uncertainty)
Set the Bias object from a guess.
void setInputs(const Vector2 &dcm, const Vector2 &zmp, const Matrix2 &R=Matrix2::Identity(), const Vector2 &CoMOffset_gamma=Vector2::Zero(), const double ZMPCoef_kappa=1)
Set the Inputs of the estimator.
double omega0_
Definition lipm-dcm-estimator.hpp:366
void updateMatricesABQ_()
set Matrices: A, B, Q
Vector2 getBias() const
Get the estimated Bias.
Vector2 previousZmp_
Definition lipm-dcm-estimator.hpp:375
Matrix2 R_
measurement noise
Definition lipm-dcm-estimator.hpp:385
void resetWithMeasurements(const Vector2 &measuredDcm, const Vector2 &measuredZMP, double yaw, bool measurementIsWithBias=true, const Vector2 &initBias=Vector2::Constant(0), const Vector2 &initBiasuncertainty=Vector2::Constant(defaultBiasUncertainty))
Resets the estimator with first measurements.
Definition lipm-dcm-estimator.hpp:149
Vector2 getUnbiasedCoMOffset() const
Get the ubiased offset on the CoM dynamics.
Definition lipm-dcm-estimator.hpp:211
void setUnbiasedDCM(const Vector2 &dcm, const Vector2 &uncertainty)
void resetWithMeasurements(const Vector2 &measuredDcm, const Vector2 &measuredZMP, const Matrix2 &yaw=Matrix2::Identity(), bool measurementIsWithBias=true, const Vector2 &initBias=Vector2::Constant(0), const Vector2 &initBiasuncertainty=Vector2::Constant(defaultBiasUncertainty))
Resets the estimator with first measurements.
void setInputs(const Vector2 &dcm, const Vector2 &zmp, double yaw, const Vector2 &CoMOffset_gamma=Vector2::Zero(), const double ZMPCoef_kappa=1)
Set the Inputs of the estimator.
Definition lipm-dcm-estimator.hpp:298
const LinearKalmanFilter & getFilter() const
Get the Kalman Filter object This can be used to run specific Advanced Kalman filter related funcions...
Definition lipm-dcm-estimator.hpp:354
void setDcmMeasureErrorStd(double)
Set the Dcm Measurement Error Standard.
LinearKalmanFilter & getFilter()
Get the Kalman Filter object This can be used to run specific Advanced Kalman filter related funcions...
Definition lipm-dcm-estimator.hpp:347
void resetWithMeasurements(const Vector2 &measuredDcm, const Vector2 &measuredZMP, const Matrix3 &rotation, bool measurementIsWithBias=true, const Vector2 &initBias=Vector2::Constant(0), const Vector2 &initBiasuncertainty=Vector2::Constant(defaultBiasUncertainty))
Resets the estimator with first measurements.
Definition lipm-dcm-estimator.hpp:172
double biasDriftStd_
Definition lipm-dcm-estimator.hpp:370
void setSamplingTime(double dt)
Set the Sampling Time.
void setBiasDriftPerSecond(double driftPerSecond)
Set the Bias Drift Per Second.
Matrix4 B_
Definition lipm-dcm-estimator.hpp:381
Defines the class of a Linear Kalman filter.
Gathers many kinds of algorithms.
Definition bidim-elastic-inv-pendulum-dyn-sys.hpp:21
Eigen::Matrix4d Matrix4
4x4 Scalar Matrix
Definition definitions.hpp:115
Eigen::Rotation2D< double > Rotation2D
2D rotations
Definition definitions.hpp:136
Eigen::Matrix2d Matrix2
2D scalar Matrix
Definition definitions.hpp:106
Eigen::Matrix3d Matrix3
3x3 Scalar Matrix
Definition definitions.hpp:109
Eigen::MatrixXd Matrix
Dynamic sized Matrix.
Definition definitions.hpp:100
Eigen::Matrix< double, 2, 1 > Vector2
2d Vector
Definition definitions.hpp:82
Implements integrators for the kinematics, in terms or rotations and translations.