WPILibC++ 2025.2.1
Loading...
Searching...
No Matches
SteadyStateKalmanFilter.h
Go to the documentation of this file.
1// Copyright (c) FIRST and other WPILib contributors.
2// Open Source Software; you can modify and/or share it under the terms of
3// the WPILib BSD license file in the root directory of this project.
4
5#pragma once
6
7#include <cmath>
8#include <stdexcept>
9#include <string>
10
11#include <Eigen/Cholesky>
12#include <wpi/SymbolExports.h>
13#include <wpi/array.h>
14
15#include "frc/DARE.h"
16#include "frc/EigenCore.h"
17#include "frc/StateSpaceUtil.h"
18#include "frc/fmt/Eigen.h"
21#include "units/time.h"
22#include "wpimath/MathShared.h"
23
24namespace frc {
25
26/**
27 * A Kalman filter combines predictions from a model and measurements to give an
28 * estimate of the true system state. This is useful because many states cannot
29 * be measured directly as a result of sensor noise, or because the state is
30 * "hidden".
31 *
32 * Kalman filters use a K gain matrix to determine whether to trust the model or
33 * measurements more. Kalman filter theory uses statistics to compute an optimal
34 * K gain which minimizes the sum of squares error in the state estimate. This K
35 * gain is used to correct the state estimate by some amount of the difference
36 * between the actual measurements and the measurements predicted by the model.
37 *
38 * This class assumes predict() and correct() are called in pairs, so the Kalman
39 * gain converges to a steady-state value. If they aren't, use KalmanFilter
40 * instead.
41 *
42 * For more on the underlying math, read
43 * https://file.tavsys.net/control/controls-engineering-in-frc.pdf chapter 9
44 * "Stochastic control theory".
45 *
46 * @tparam States Number of states.
47 * @tparam Inputs Number of inputs.
48 * @tparam Outputs Number of outputs.
49 */
50template <int States, int Inputs, int Outputs>
52 public:
56
59
60 /**
61 * Constructs a steady-state Kalman filter with the given plant.
62 *
63 * See
64 * https://docs.wpilib.org/en/stable/docs/software/advanced-controls/state-space/state-space-observers.html#process-and-measurement-noise-covariance-matrices
65 * for how to select the standard deviations.
66 *
67 * @param plant The plant used for the prediction step.
68 * @param stateStdDevs Standard deviations of model states.
69 * @param measurementStdDevs Standard deviations of measurements.
70 * @param dt Nominal discretization timestep.
71 * @throws std::invalid_argument If the system is undetectable.
72 */
74 const StateArray& stateStdDevs,
75 const OutputArray& measurementStdDevs,
76 units::second_t dt) {
77 m_plant = &plant;
78
79 auto contQ = MakeCovMatrix(stateStdDevs);
80 auto contR = MakeCovMatrix(measurementStdDevs);
81
84 DiscretizeAQ<States>(plant.A(), contQ, dt, &discA, &discQ);
85
86 auto discR = DiscretizeR<Outputs>(contR, dt);
87
88 const auto& C = plant.C();
89
90 if (!IsDetectable<States, Outputs>(discA, C)) {
91 std::string msg = fmt::format(
92 "The system passed to the Kalman filter is undetectable!\n\n"
93 "A =\n{}\nC =\n{}\n",
94 discA, C);
95
97 throw std::invalid_argument(msg);
98 }
99
100 if (auto P = DARE<States, Outputs>(discA.transpose(), C.transpose(), discQ,
101 discR)) {
102 // S = CPCᵀ + R
103 Matrixd<Outputs, Outputs> S = C * P.value() * C.transpose() + discR;
104
105 // We want to put K = PCᵀS⁻¹ into Ax = b form so we can solve it more
106 // efficiently.
107 //
108 // K = PCᵀS⁻¹
109 // KS = PCᵀ
110 // (KS)ᵀ = (PCᵀ)ᵀ
111 // SᵀKᵀ = CPᵀ
112 //
113 // The solution of Ax = b can be found via x = A.solve(b).
114 //
115 // Kᵀ = Sᵀ.solve(CPᵀ)
116 // K = (Sᵀ.solve(CPᵀ))ᵀ
117 m_K = S.transpose().ldlt().solve(C * P.value().transpose()).transpose();
118 } else if (P.error() == DAREError::QNotSymmetric ||
120 std::string msg =
121 fmt::format("{}\n\nQ =\n{}\n", to_string(P.error()), discQ);
122
124 throw std::invalid_argument(msg);
125 } else if (P.error() == DAREError::RNotSymmetric ||
126 P.error() == DAREError::RNotPositiveDefinite) {
127 std::string msg =
128 fmt::format("{}\n\nR =\n{}\n", to_string(P.error()), discR);
129
131 throw std::invalid_argument(msg);
132 } else if (P.error() == DAREError::ABNotStabilizable) {
133 std::string msg = fmt::format(
134 "The (A, C) pair is not detectable.\n\nA =\n{}\nC =\n{}\n",
135 to_string(P.error()), discA, C);
136
138 throw std::invalid_argument(msg);
139 } else if (P.error() == DAREError::ACNotDetectable) {
140 std::string msg = fmt::format("{}\n\nA =\n{}\nQ =\n{}\n",
141 to_string(P.error()), discA, discQ);
142
144 throw std::invalid_argument(msg);
145 }
146
147 Reset();
148 }
149
152
153 /**
154 * Returns the steady-state Kalman gain matrix K.
155 */
156 const Matrixd<States, Outputs>& K() const { return m_K; }
157
158 /**
159 * Returns an element of the steady-state Kalman gain matrix K.
160 *
161 * @param i Row of K.
162 * @param j Column of K.
163 */
164 double K(int i, int j) const { return m_K(i, j); }
165
166 /**
167 * Returns the state estimate x-hat.
168 */
169 const StateVector& Xhat() const { return m_xHat; }
170
171 /**
172 * Returns an element of the state estimate x-hat.
173 *
174 * @param i Row of x-hat.
175 */
176 double Xhat(int i) const { return m_xHat(i); }
177
178 /**
179 * Set initial state estimate x-hat.
180 *
181 * @param xHat The state estimate x-hat.
182 */
183 void SetXhat(const StateVector& xHat) { m_xHat = xHat; }
184
185 /**
186 * Set an element of the initial state estimate x-hat.
187 *
188 * @param i Row of x-hat.
189 * @param value Value for element of x-hat.
190 */
191 void SetXhat(int i, double value) { m_xHat(i) = value; }
192
193 /**
194 * Resets the observer.
195 */
196 void Reset() { m_xHat.setZero(); }
197
198 /**
199 * Project the model into the future with a new control input u.
200 *
201 * @param u New control input from controller.
202 * @param dt Timestep for prediction.
203 */
204 void Predict(const InputVector& u, units::second_t dt) {
205 m_xHat = m_plant->CalculateX(m_xHat, u, dt);
206 }
207
208 /**
209 * Correct the state estimate x-hat using the measurements in y.
210 *
211 * @param u Same control input used in the last predict step.
212 * @param y Measurement vector.
213 */
214 void Correct(const InputVector& u, const OutputVector& y) {
215 const auto& C = m_plant->C();
216 const auto& D = m_plant->D();
217
218 // x̂ₖ₊₁⁺ = x̂ₖ₊₁⁻ + K(y − (Cx̂ₖ₊₁⁻ + Duₖ₊₁))
219 m_xHat += m_K * (y - (C * m_xHat + D * u));
220 }
221
222 private:
224
225 /**
226 * The steady-state Kalman gain matrix.
227 */
229
230 /**
231 * The state estimate.
232 */
233 StateVector m_xHat;
234};
235
236extern template class EXPORT_TEMPLATE_DECLARE(WPILIB_DLLEXPORT)
237 SteadyStateKalmanFilter<1, 1, 1>;
238extern template class EXPORT_TEMPLATE_DECLARE(WPILIB_DLLEXPORT)
239 SteadyStateKalmanFilter<2, 1, 1>;
240
241} // namespace frc
#define WPILIB_DLLEXPORT
Definition SymbolExports.h:36
#define EXPORT_TEMPLATE_DECLARE(export)
Definition SymbolExports.hpp:92
A plant defined using state-space notation.
Definition LinearSystem.h:35
constexpr const Matrixd< Outputs, States > & C() const
Returns the output matrix C.
Definition LinearSystem.h:130
constexpr const Matrixd< States, States > & A() const
Returns the system matrix A.
Definition LinearSystem.h:104
A Kalman filter combines predictions from a model and measurements to give an estimate of the true sy...
Definition SteadyStateKalmanFilter.h:51
void SetXhat(const StateVector &xHat)
Set initial state estimate x-hat.
Definition SteadyStateKalmanFilter.h:183
SteadyStateKalmanFilter & operator=(SteadyStateKalmanFilter &&)=default
SteadyStateKalmanFilter(LinearSystem< States, Inputs, Outputs > &plant, const StateArray &stateStdDevs, const OutputArray &measurementStdDevs, units::second_t dt)
Constructs a steady-state Kalman filter with the given plant.
Definition SteadyStateKalmanFilter.h:73
void Correct(const InputVector &u, const OutputVector &y)
Correct the state estimate x-hat using the measurements in y.
Definition SteadyStateKalmanFilter.h:214
void Predict(const InputVector &u, units::second_t dt)
Project the model into the future with a new control input u.
Definition SteadyStateKalmanFilter.h:204
void Reset()
Resets the observer.
Definition SteadyStateKalmanFilter.h:196
void SetXhat(int i, double value)
Set an element of the initial state estimate x-hat.
Definition SteadyStateKalmanFilter.h:191
double K(int i, int j) const
Returns an element of the steady-state Kalman gain matrix K.
Definition SteadyStateKalmanFilter.h:164
double Xhat(int i) const
Returns an element of the state estimate x-hat.
Definition SteadyStateKalmanFilter.h:176
const Matrixd< States, Outputs > & K() const
Returns the steady-state Kalman gain matrix K.
Definition SteadyStateKalmanFilter.h:156
Vectord< Outputs > OutputVector
Definition SteadyStateKalmanFilter.h:55
SteadyStateKalmanFilter(SteadyStateKalmanFilter &&)=default
Vectord< Inputs > InputVector
Definition SteadyStateKalmanFilter.h:54
const StateVector & Xhat() const
Returns the state estimate x-hat.
Definition SteadyStateKalmanFilter.h:169
Vectord< States > StateVector
Definition SteadyStateKalmanFilter.h:53
This class is a wrapper around std::array that does compile time size checking.
Definition array.h:26
static void ReportError(const S &format, Args &&... args)
Definition MathShared.h:62
Definition CAN.h:11
void DiscretizeAQ(const Matrixd< States, States > &contA, const Matrixd< States, States > &contQ, units::second_t dt, Matrixd< States, States > *discA, Matrixd< States, States > *discQ)
Discretizes the given continuous A and Q matrices.
Definition Discretization.h:71
Eigen::Matrix< double, Rows, Cols, Options, MaxRows, MaxCols > Matrixd
Definition EigenCore.h:21
Eigen::Vector< double, Size > Vectord
Definition EigenCore.h:12
Matrixd< Outputs, Outputs > DiscretizeR(const Matrixd< Outputs, Outputs > &R, units::second_t dt)
Returns a discretized version of the provided continuous measurement noise covariance matrix.
Definition Discretization.h:113
constexpr Matrixd< sizeof...(Ts), sizeof...(Ts)> MakeCovMatrix(Ts... stdDevs)
Creates a covariance matrix from the given vector for use with Kalman filters.
Definition StateSpaceUtil.h:75
wpi::expected< Eigen::Matrix< double, States, States >, DAREError > DARE(const Eigen::Matrix< double, States, States > &A, const Eigen::Matrix< double, States, Inputs > &B, const Eigen::Matrix< double, States, States > &Q, const Eigen::Matrix< double, Inputs, Inputs > &R, bool checkPreconditions=true)
Computes the unique stabilizing solution X to the discrete-time algebraic Riccati equation:
Definition DARE.h:165
constexpr std::string_view to_string(const DAREError &error)
Converts the given DAREError enum to a string.
Definition DARE.h:39
@ QNotPositiveSemidefinite
Q was not positive semidefinite.
@ RNotSymmetric
R was not symmetric.
@ QNotSymmetric
Q was not symmetric.
@ ACNotDetectable
(A, C) pair where Q = CᵀC was not detectable.
@ RNotPositiveDefinite
R was not positive definite.
@ ABNotStabilizable
(A, B) pair was not stabilizable.
bool IsDetectable(const Matrixd< States, States > &A, const Matrixd< Outputs, States > &C)
Returns true if (A, C) is a detectable pair.
Definition StateSpaceUtil.h:309
#define S(label, offset, message)
Definition Errors.h:113