-
Notifications
You must be signed in to change notification settings - Fork 7
/
Copy pathmainpage.dox
383 lines (316 loc) · 12.7 KB
/
mainpage.dox
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
/*!
\mainpage kalman-cpp
\section linear-kalman Kalman Filter for Linear Systems
\subsection definition_kf Definition
A linear system is described as follows.
\f[x_k = Ax_{k-1} + Bu_{k-1} + v_{k-1}\f]
\f[z_k = Hx_k + w_k\f]
where:\n
\f$v\f$ is the process noise (Gaussian with covariance Q)\n
\f$w\f$ is the measurement noise (Gaussian with covariance R)\n
\f$A\f$ is the system matrix\n
\f$B\f$ is the input matrix\n
\f$H\f$ is the output matrix\n
\f$x\f$ is the state vector\n
\f$z\f$ is the output vector\n
\f$u\f$ is the input vector\n\n
The noise covariance matrices must follow these conditions:\n
\f[Q = Q^T\f]
\f[R = R^T\f]
\subsection example_kf Examples
An example is provided in main1.cpp solves <a href="http://www.mathworks.com/matlabcentral/fileexchange/18465-learning-the-kalman-filter-in-simulink-v2-1/content/html/runkalmanfilter.html">the car voltmeter problem</a> using this library.
The system model for this problem is governed by:
\n
\f[x_k = 12 + v_{k-1}\f]
\f[z_k = x_k + w_ḳ\f]
\n
Let us extract the model parameters from the equation above: \f$A = 0\f$, \f$B = 1\f$, \f$H = 1\f$, \f$u = 12\f$, \f$x_0 = 12\f$, \f$Q=4\f$, and \f$R=4\f$.
Now that we know all the parameter values, we can then set up the following piece of code:
\n
\code{.cpp}
mat A(1,1), B(1,1), H(1,1), Q(1,1), R(1,1);
colvec u(1);
A << 0;
B << 1;
Q << 4;
H << 1;
R << 4;
u << 12.0;
\endcode
The first and second line of the code above are used to create all the necessary matrices with correct dimensions. In the following lines, we apply the correct value for each matrix.
Our next step is to create an instance of class class and initialize it with the previously created
\f$A\f$, \f$B\f$, \f$Q\f$, \f$H\f$, and \f$R\f$ matrices.
\n
\code{.cpp}
KF kalman;
kalman.InitSystem(A, B, H, Q, R);
\endcode
And or our last step, we run the Kalman procedure repetitively as shown in the code below.
\n
\code{.cpp}
for (int k = 0; k < 100; k++) {
kalman.Kalmanf(u);
// The rest of the code
}
\endcode
\n
@image HTML ../images/ex1.png
\section non-linear-kalman Extended Kalman Filter (EKF) for nonlinear systems
\subsection definition_ekf Definition
A non-linear system is described by:
\f[x_k = f(x_{k-1}, u_{k-1}) + v_{k-1}\f]
\f[z_k = h(x_k) + w_k\f]
where:\n
\f$f\f$ is the dynamic model of the system\n
\f$h\f$ is the measurement model of the system\n
\f$v\f$ is the process noise (Gaussian with covariance Q)\n
\f$w\f$ is the measurement noise (Gaussian with covariance R)\n
\f$x\f$ is the state vector\n
\f$z\f$ is the output vector\n
\f$u\f$ is the input vector\n\n
The noise covariance matrices must follow these conditions:\n
\f[Q = Q^T\f]
\f[R = R^T\f]
\subsection example_ekf Example
Here, let us take main4.cpp as an example. This content of this file is adapted from this
<a href="http://ch.mathworks.com/matlabcentral/fileexchange/38302-kalman-filter-package/content//Kalman%20Filter%20Package/Examples/ExtendedKalmanFilterDemo.m">MATLAB Central page.</a>
\n\n
Let us define a nonlinear system as follows.
\n\n
\f[ f = \begin{bmatrix} \sin(x_2(k-1))(k-1) \\ x_2(k-1) \end{bmatrix} \f]
\f[ h = \begin{bmatrix} x_1(k) \\ x_2(k) \end{bmatrix} \f]
\f[ x_0 = \begin{bmatrix} 0 \\ \frac{1 \pi}{500} \end{bmatrix}\f]
\n
A nonlinear function does not have a strict tempate as in a linear function. Therefore, we need to let the user to derive the EKF class for flexible implementation of a nonlinear model for both the process and the the measurement.
\n
\code{.cpp}
class MyEKF: public EKF
{
public:
virtual colvec f(const colvec& x, const colvec& u) {
colvec xk(nOutputs_);
xk(0) = sin(x(1) * u(0));
xk(1) = x(1);
return xk;
}
virtual colvec h(const colvec& x) {
colvec zk(nOutputs_);
zk(0) = x(0);
zk(1) = x(1);
return zk;
}
};
\endcode
The next steps below are very similiar to the previous example in the linear Kalman filter section.
First, we start by creating all necessary matrices and set their values accordingly. Afterward, we continue with initializing. the covariance matrices as well as the state vector. In the final step, we run the EKF procedurerepetitively while logging necessary data.
\n\n
\code{.cpp}
int main(int argc, char** argv)
{
//
// Log the result into a tab delimitted file, later we can open
// it with Matlab. Use: plot_data4.m to plot the results.
//
ofstream log_file;
log_file.open("log_file4.txt");
int n_states = 2;
int n_outputs = 2;
mat Q(2, 2);
mat R(2, 2);
Q << 0.001 << 0 << endr
<< 0 << 0 << endr;
R << 0.1 << 0 << endr
<< 0 << 0.01 << endr;
colvec x0(2);
x0 << 0 << 1 * M_PI / 500;
colvec u(1);
MyEKF myekf;
myekf.InitSystem(n_states, n_outputs, Q, R);
myekf.InitSystemState(x0);
for (int k = 0; k < 1000; k ++) {
u(0) = k;
myekf.EKalmanf(u);
colvec *x = myekf.GetCurrentState();
colvec *x_m = myekf.GetCurrentEstimatedState();
colvec *z = myekf.GetCurrentOutput();
log_file << k
<< '\t' << z->at(0,0) << '\t' << x->at(0,0) << '\t' << x_m->at(0,0)
<< '\t' << z->at(1,0) << '\t' << x->at(1,0) << '\t' << x_m->at(1,0)
<< '\n';
}
log_file.close();
return 0;
}
\endcode
@image HTML ../images/ex4.png
\section unscented-kalman Unscented Kalman Filter (UKF) for a nonlinear system
Another type of Kalman Filter for a nonlinear system is the Unscented Kalman Filter. It addresses the accuracy problem which arises during linearization process of an Extended Kalman filter when Jacobian is used. The Unscented Kalman Filter is implemented in class UKF and the steps for using the provided UKF class is exactly similiar as in using the EKF class.
\n\n
The implementation of the UKF in this C++ library is based on the following paper:
\n\n
Wan, E. A., & Van Der Merwe, R. (2006). The unscented Kalman filter for nonlinear estimation. Proceedings of the IEEE 2000 Adaptive Systems for Signal Processing, Communications, and Control Symposium (Cat. No.00EX373), 31(2), 153–158.
\n\n
The paper above can be downloaded from <a href="https://doi.org/10.1109/ASSPCC.2000.882463">here</a>.
\subsection example_ukf Example
Here, we take main10.cpp as an example. This example is adapted from this
\n
This example is taken from<a href="https://www.mathworks.com/matlabcentral/fileexchange/18217-learning-the-unscented-kalman-filter">here</a>. The nonlinear system is described as follows.
\n
\f[ f = \begin{bmatrix} x_2(k) \\ x_3(k) \\ 0.005 \, x_1(k) \bigg(x_2(k) + x_3(k) \bigg) \end{bmatrix} \f]
\f[ h = \begin{bmatrix} x_1(k) \end{bmatrix} \f]
\f[ x_0 = \begin{bmatrix} 0 \\ 0 \\ 1 \end{bmatrix}\f]
\n
The following lines of codes show how to use the UKF class. As mentioned before, we need to first derive the UKF class and define the virtual functions f (process function) and h (output function).
\n\n
\code{.cpp}
class MyUKF : public UKF
{
public:
virtual colvec f(const colvec& x, const colvec& u) {
colvec xk(nStates_);
xk.at(0) = x(1);
xk.at(1) = x(2);
xk.at(2) = 0.05*x(0)*(x(1)+x(2));
return xk;
}
virtual colvec h(const colvec& x) {
colvec zk(nOutputs_);
zk(0) = x(0);
return zk;
}
};
/////////////////////////////////////////////////////////////////////////////
int main(int argc, char** argv)
{
mat Q(3, 3);
mat R(1, 1);
MyUKF myukf;
double r = 0.1;
double q = 0.1;
Q = eye(3,3)*q*q;
R << r*r << endr;
colvec x0(3);
x0 << 0 << 0 << 1;
mat P0 = eye<mat>(3, 3);
P0 = P0;
colvec u;
// No inputs
u = u.zeros();
myukf.InitSystem(3, 1, Q, R);
myukf.InitSystemState(x0);
myukf.InitSystemStateCovariance(P0);
for (int k = 0; k < 20; k++) {
myukf.UKalmanf(u);
colvec *x = myukf.GetCurrentState();
colvec *x_m = myukf.GetCurrentEstimatedState();
colvec *z = myukf.GetCurrentOutput();
colvec *z_m = myukf.GetCurrentEstimatedOutput();
}
return 0;
}
\endcode
\n\n
These are the plots that show comparisons between measurements and estimations of the three states: \f$x_1\f$, \f$x_2\f$ and \f$x_3\f$.
\n
@image HTML ../images/ex10.png
\section practical-application Practical application: Kalman filter for noisy measurements
The examples we have so far are theoretical. Very often, what we would like to do is to reduce noise from pre-acquired measurement data. There are several reasons why we want to use Kalman filter. For example, noise has a vast spectrum. Thus, using a frequency-based filter hurts the data.
\n\n
Principally, there are two scenarios of using the Kalman filtering techniquie here. The first scenario is by first simulating the system, as shown in the figure below.
In this scenario, we only need to supply \f$u_k\f$ to the Kalman filter function.
The Kalman procedure will give us four outputs as a result: \f$x_k\f$, \f$z_k\f$, \f$\hat{x}_k\f$, and \f$\hat{z}_k\f$.
\f$x_k\f$ and \f$z_k\f$ are called the true states and the true outputs, respectvely. They are noisy.
\f$\hat{x}_k\f$, and \f$\hat{z}_k\f$ are called the estimated states and the estimated outputs, respectively.
They are filtered, look smoother and closer to the noiseless true states.
\n\n
The function prototype for the scenario above can be written as:
\n
[\f$x_k\f$, \f$z_k\f$, \f$\hat{x}_k\f$, \f$\hat{z}_k\f$] = function kalmanf( \f$u_k\f$)
\n
@image HTML ../images/Kalman_concept1.jpg
\n\n
The second scenario is used when the measurements are available. Thus, simulating the system becomes unnecessary. In such a scenario, we need to supply both \f$z_k\f$ and \f$u_k\f$ to the kalman filter function.
The Kalman filter will give us 2 outputs: \f$\hat{x}_k\f$ (the estimated system sates) and \f$\hat{z}_k\f$ (the estimated system outputs).
\n\n
The function prototype for this scenario can be written as:
\n
[\f$\hat{x}_k\f$, \f$\hat{z}_k\f$] = function kalmanf( \f$z_k\f$, \f$u_k\f$)
\n
@image HTML ../images/Kalman_concept2.jpg
\n
The second scenario is useful for smoothing noisy measurment data. However, both scenarios are availabe in this library.
\n
\subsection practical_example Example
As one example, let us see the problem provided in main1.cpp.
This time, however, we assume we have already had several noisy data points from measurements.
Therfore, we will not need the Kalman procedure to simulate the system.
\n
The system model can be formulated into:
\n
\f[x_k = x_{k-1} + v_{k-1}\f]
\f[z_k = x_k + w_ḳ\f]
\n
We then generate random number as the prerecorded noisy data (**z_measurement**), as follows.
\n
\code{.cpp}
// Assume we have a noisy signal that is acquired from a measurement tool
double w = 2; // Stdev of the noise, in reality we don't know this
colvec z_measurement(1);
z_measurement.randn(1);
z_measurement = z_measurement * w + 12.0;
\endcode
After than, **z_measurement** is sent to the Kalman filter function, as the following:
\n\n
\code{.cpp}
// Put the z_measurment to the Kalman filter
kalman.Kalmanf(z_measurement, u);
colvec *z_m = kalman.GetCurrentEstimatedOutput();
colvec *x_m = kalman.GetCurrentEstimatedState();
\endcode
\n\n
The **kalman.Kalmanf(z_measurement, u)** function above has different parameters with what we can find in main1.cpp.
Since **z_measurement** is available, we then send it as a parameter for the **kalmanf** function.
The processes above are done in several iterations according to the number of the availabe data points.
The final codes will look like as the following:
\n\n
\code{.cpp}
int main(int argc, char** argv)
{
// Log the result into a tab delimitted file, later we can open
// it with Matlab. Use: plot_data7.m to plot the results.
ofstream log_file;
log_file.open("log_file7.txt");
// Define the system and initialize the Kalman filter
mat A(1,1), B(1,1), H(1,1), Q(1,1), R(1,1);
colvec u(1);
colvec x0(1);
A << 1;
B << 1;
Q << 0.01; // Heuristic tuning parameter
H << 1;
R << 1; // Heuristic tuning parameter
u << 0;
KF kalman;
kalman.InitSystem(A, B, H, Q, R);
int N = 500;
for (int k = 0; k < N; k++) {
// Assume we have a noisy signal that is acquired from a measurement tool
double w = 2; // Stdev of the noise, in reality we don't know this
colvec z_measurement(1);
z_measurement.randn(1);
z_measurement = z_measurement * w + 12.0;
// Put the z_measurment to the Kalman filter
kalman.Kalmanf(z_measurement, u);
colvec *z_m = kalman.GetCurrentEstimatedOutput();
colvec *x_m = kalman.GetCurrentEstimatedState();
log_file << k << '\t' << z_measurement.at(0,0) << '\t' << z_m->at(0,0)
<< '\n';
}
return 0;
}
\endcode
\n\n
The values that are given to **Q** and **R** are defined **heuristically** since we do not know the actual variances of the noise of the data.
\n
@image HTML ../images/ex7.png
*/