state-observation 1.7.0
Loading...
Searching...
No Matches
unidim-lipm-dcm-estimator.hpp
Go to the documentation of this file.
1
12#ifndef UNIDIMLIPMDCMBIASESTIMATOR_HPP
13#define UNIDIMLIPMDCMBIASESTIMATOR_HPP
14
15#include <state-observation/api.h>
18
19namespace stateObservation
20{
21
33
35{
36private:
37 constexpr static double defaultDt_ = 0.005;
38 constexpr static double defaultOmega_ = tools::sqrt(cst::gravityConstant);
39
40public:
42 constexpr static double defaultBiasDriftSecond = 0.002;
43
45 constexpr static double defaultZmpErrorStd = 0.005;
46 constexpr static double defaultDcmErrorStd = 0.01;
47
49 constexpr static double defaultDCMUncertainty = 0.01;
50 constexpr static double defaultBiasUncertainty = 0.01;
51
64 UnidimLipmDcmEstimator(double dt = defaultDt_,
65 double omega_0 = defaultOmega_,
66 double biasDriftPerSecondStd = defaultBiasDriftSecond,
67 double initDcm = 0,
68 double initZMP = 0,
69 double initBias = 0,
70 double dcmMeasureErrorStd = defaultDcmErrorStd,
71 double zmpMeasureErrorStd = defaultZmpErrorStd,
72 double initDcmUncertainty = defaultDCMUncertainty,
73 double initBiasUncertainty = defaultBiasUncertainty);
74
89 double measuredZMP,
90 bool measurementIsWithBias = true,
91 double biasDriftPerSecondStd = defaultBiasDriftSecond,
92 double dcmMeasureErrorStd = defaultDcmErrorStd,
93 double zmpMeasureErrorStd = defaultZmpErrorStd,
94 double initBias = 0,
95 double initBiasuncertainty = defaultBiasUncertainty);
96
99
104
108 void setSamplingTime(double dt);
109
113 void setBias(double bias);
114
118 void setBias(double bias, double uncertainty);
119
124
128 void setUnbiasedDCM(double dcm);
129
134 void setUnbiasedDCM(double dcm, double uncertainty);
135
139
143
148 void setInputs(double dcm, double zmp);
149
154 {
155 return filter_.getEstimatedState(filter_.getMeasurementTime());
156 }
157
162 double getUnbiasedDCM() const;
163
167 double getBias() const;
168
173 {
174 return filter_;
175 }
176
179 inline const LinearKalmanFilter & getFilter() const
180 {
181 return filter_;
182 }
183
184private:
186 void updateMatricesABQ_();
187
188 double omega0_;
189 double dt_;
190 double biasDriftStd_;
191 double zmpErrorStd_;
192
193 double previousZmp_;
194
195 LinearKalmanFilter filter_;
196 Matrix2 A_;
197 Vector2 B_;
199 Vector2 C_;
201 Matrix1 R_;
202
204 Matrix2 Q_;
205
206public:
208};
209
210} // namespace stateObservation
211
212#endif
The class of a Linear Kalman filter.
Definition linear-kalman-filter.hpp:48
1D version of the estimation of a bias betweeen the divergent component of motion and the correspondi...
Definition unidim-lipm-dcm-estimator.hpp:35
~UnidimLipmDcmEstimator()
Destroy the Lipm Dcm Bias Estimator object.
Definition unidim-lipm-dcm-estimator.hpp:98
void setBias(double bias)
Set the Bias object from a guess.
double getUnbiasedDCM() const
Get the Unbiased DCM filtered by the estimator.
LinearKalmanFilter & getFilter()
Get the Kalman Filter object This can be used to run specific Advanced Kalman filter related funcions...
Definition unidim-lipm-dcm-estimator.hpp:172
double getBias() const
Get the estimated Bias.
void setUnbiasedDCM(double dcm)
set the real DCM position from a guess
UnidimLipmDcmEstimator(double dt=defaultDt_, double omega_0=defaultOmega_, double biasDriftPerSecondStd=defaultBiasDriftSecond, double initDcm=0, double initZMP=0, double initBias=0, double dcmMeasureErrorStd=defaultDcmErrorStd, double zmpMeasureErrorStd=defaultZmpErrorStd, double initDcmUncertainty=defaultDCMUncertainty, double initBiasUncertainty=defaultBiasUncertainty)
Construct a new Unidimensional Lipm Dcm Estimator.
void setSamplingTime(double dt)
Set the Sampling Time.
void setDcmMeasureErrorStd(double)
Set the Dcm Measurement Error Standard.
void setUnbiasedDCM(double dcm, double uncertainty)
set the real DCM position from a guess
void setInputs(double dcm, double zmp)
Set the Inputs of the estimator.
void resetWithInputs(double measuredDcm, double measuredZMP, bool measurementIsWithBias=true, double biasDriftPerSecondStd=defaultBiasDriftSecond, double dcmMeasureErrorStd=defaultDcmErrorStd, double zmpMeasureErrorStd=defaultZmpErrorStd, double initBias=0, double initBiasuncertainty=defaultBiasUncertainty)
Construct a new Lipm Dcm Bias Estimator object.
void setBias(double bias, double uncertainty)
Set the Bias object from a guess.
Vector2 update()
Runs the estimation. Needs to be called every timestep.
Definition unidim-lipm-dcm-estimator.hpp:153
const LinearKalmanFilter & getFilter() const
Get the Kalman Filter object This can be used to run specific Advanced Kalman filter related funcions...
Definition unidim-lipm-dcm-estimator.hpp:179
void setBiasDriftPerSecond(double driftPerSecond)
Set the Bias Drift Per Second.
void setLipmNaturalFrequency(double omega_0)
Set the Lipm Natural Frequency.
void setZmpMeasureErrorStd(double)
Set the Zmp Measurement Error Stamdard devbiation.
Defines the class of a Linear Kalman filter.
Gathers many kinds of algorithms.
constexpr double gravityConstant
Definition definitions.hpp:663
double constexpr sqrt(double x)
Constexpr version of the square root.
Definition miscellaneous-algorithms.hpp:186
Definition bidim-elastic-inv-pendulum-dyn-sys.hpp:21
Eigen::Matrix2d Matrix2
2D scalar Matrix
Definition definitions.hpp:106
Eigen::MatrixXd Matrix
Dynamic sized Matrix.
Definition definitions.hpp:100
Eigen::Matrix< double, 2, 1 > Vector2
2d Vector
Definition definitions.hpp:82
Eigen::Matrix< double, 1, 1 > Matrix1
1D scalar Matrix
Definition definitions.hpp:103