SigPack - the C++ signal processing library
 
Loading...
Searching...
No Matches
kalman_linear.cpp
Go to the documentation of this file.
1
8#include "sigpack.h"
9
10using namespace std;
11using namespace sp;
12
13int main()
14{
15 // Number of samples
16 arma::uword Nsamp = 120;
17 arma::uword N = 6; // Nr of states
18 arma::uword M = 2; // Nr of measurements
19 arma::uword L = 1; // Nr of inputs
20
21 // Instatiate a Kalman Filter
22 KF kalman( N, M, L );
23
24 // Initialisation and setup of system
25 double P0 = 100;
26 double Q0 = 0.00003;
27 double R0 = 25;
28
29 // Meas interval
30 double dT = 0.1;
31
32 arma::mat x = { 0, 10, 1, 50, -0.08, -9 };
33 kalman.set_state_vec( 0.9 * x.t() );
34 // [X,Y,Vx,Vy,Ax,Ay] use position, velocity and acceleration as states
35 arma::mat A = { { 1, 0, dT, 0, dT * dT / 2, 0 }, { 0, 1, 0, dT, 0, dT * dT / 2 }, { 0, 0, 1, 0, dT, 0 }, { 0, 0, 0, 1, 0, dT }, { 0, 0, 0, 0, 1, 0 }, { 0, 0, 0, 0, 0, 1 } };
36 kalman.set_trans_mat( A );
37
38 arma::mat H = { { 1, 0, 0, 0, 0, 0 }, { 0, 1, 0, 0, 0, 0 } };
39 kalman.set_meas_mat( H );
40
41 arma::mat P = P0 * arma::eye( N, N );
42 kalman.set_err_cov( P );
43 arma::mat Q = arma::zeros( N, N );
44 Q( N - 2, N - 2 ) = 1;
45 Q( N - 1, N - 1 ) = 1;
46 Q = Q0 * Q;
47 kalman.set_proc_noise( Q );
48
49 arma::mat R = R0 * arma::eye( M, M );
50 kalman.set_meas_noise( R );
51
52 // Create simulation data
53 arma::mat z( M, Nsamp, arma::fill::zeros );
54 arma::mat z0( M, Nsamp, arma::fill::zeros );
55
56 arma::mat xx( N, 1, arma::fill::zeros );
57 xx = x.t();
58 for( arma::uword n = 1; n < Nsamp; n++ )
59 {
60 xx = A * xx + 0.1 * Q * arma::randn( N, 1 );
61 z0.col( n ) = H * xx; //
62 }
63
64 z.row( 0 ) = z0.row( 0 ) + 0.001 * R0 * arma::randn( 1, Nsamp );
65 z.row( 1 ) = z0.row( 1 ) + 0.8 * R0 * arma::randn( 1, Nsamp );
66
67 arma::mat x_log( N, Nsamp );
68 arma::mat e_log( M, Nsamp );
69 arma::cube P_log( N, N, Nsamp );
70 arma::mat xs_log( M, Nsamp );
71 arma::cube Ps_log( N, N, Nsamp );
72
73 // Kalman filter loop
74 for( arma::uword n = 0; n < Nsamp; n++ )
75 {
76 kalman.predict();
77 kalman.update( z.col( n ) );
78
79 x_log.col( n ) = kalman.get_state_vec();
80 e_log.col( n ) = kalman.get_err();
81 P_log.slice( n ) = kalman.get_err_cov();
82 }
83
84 // RTS smoother
85 kalman.rts_smooth( x_log, P_log, xs_log, Ps_log );
86
87 // Display result
88 gplot gp0;
89 gp0.window( "Plot", 10, 10, 500, 500 );
90 gp0.set_term( "qt" );
91 gp0.plot_add( z0.row( 0 ), z0.row( 1 ), "True Y", "lines dashtype 2 linecolor \"black\"" );
92 gp0.plot_add( z.row( 0 ), z.row( 1 ), "Meas Y", "points" );
93 gp0.plot_add( x_log.row( 0 ), x_log.row( 1 ), "Kalman" );
94 gp0.plot_add( xs_log.row( 0 ), xs_log.row( 1 ), "RTS smooth" );
95 gp0.plot_show();
96
97 return 0;
98}
Kalman filter class.
Definition kalman.h:125
Gnuplot class.
Definition gplot.h:27
void plot_show(void)
Show plots.
Definition gplot.h:401
void set_term(const char *ttype)
Set output terminal.
Definition gplot.h:757
void plot_add(const T1 &x, const T2 &y, const std::string lb, const std::string ls="lines")
Push plot y vs. x with label and linespec.
Definition gplot.h:325
void window(const int fig, const char *name, const int x, const int y, const int width, const int height)
Configure the figure used Windows environment.
Definition gplot.h:190
arma::mat get_err_cov(void)
Definition kalman.h:261
void rts_smooth(const arma::mat &Xf, const arma::cube &Pf, arma::mat &Xs, arma::cube &Ps)
Rauch-Tung-Striebel smoother. See https://users.aalto.fi/~ssarkka/course_k2012/full_course_booklet_20...
Definition kalman.h:307
arma::mat get_state_vec(void)
Definition kalman.h:246
void set_meas_noise(const arma::mat &_R)
Definition kalman.h:224
void set_state_vec(const arma::mat &_x)
Definition kalman.h:194
arma::mat get_err(void)
Definition kalman.h:251
void set_meas_mat(const arma::mat &_H)
Definition kalman.h:209
void set_err_cov(const arma::mat &_P)
Definition kalman.h:214
void set_trans_mat(const arma::mat &_A)
Definition kalman.h:199
void set_proc_noise(const arma::mat &_Q)
Definition kalman.h:219
void predict(const arma::mat u)
Predict the internal states using a control input.
Definition kalman.h:270
void update(const arma::mat z)
Correct and update the internal states.
Definition kalman.h:288
int main()
Definition base.h:8