1 | /*! |
---|
2 | \file |
---|
3 | \brief Bayesian Filtering for linear Gaussian models (Kalman Filter) and extensions |
---|
4 | \author Vaclav Smidl. |
---|
5 | |
---|
6 | ----------------------------------- |
---|
7 | BDM++ - C++ library for Bayesian Decision Making under Uncertainty |
---|
8 | |
---|
9 | Using IT++ for numerical operations |
---|
10 | ----------------------------------- |
---|
11 | */ |
---|
12 | |
---|
13 | #ifndef KF_H |
---|
14 | #define KF_H |
---|
15 | |
---|
16 | |
---|
17 | #include "../math/functions.h" |
---|
18 | #include "../stat/exp_family.h" |
---|
19 | #include "../math/chmat.h" |
---|
20 | #include "../base/user_info.h" |
---|
21 | //#include <../applications/pmsm/simulator_zdenek/ekf_example/pmsm_mod.h> |
---|
22 | |
---|
23 | namespace bdm { |
---|
24 | |
---|
25 | /*! |
---|
26 | * \brief Basic elements of linear state-space model |
---|
27 | |
---|
28 | Parameter evolution model:\f[ x_{t+1} = A x_{t} + B u_t + Q^{1/2} e_t \f] |
---|
29 | Observation model: \f[ y_t = C x_{t} + C u_t + R^{1/2} w_t. \f] |
---|
30 | Where $e_t$ and $w_t$ are mutually independent vectors of Normal(0,1)-distributed disturbances. |
---|
31 | */ |
---|
32 | template<class sq_T> |
---|
33 | class StateSpace { |
---|
34 | protected: |
---|
35 | //! Matrix A |
---|
36 | mat A; |
---|
37 | //! Matrix B |
---|
38 | mat B; |
---|
39 | //! Matrix C |
---|
40 | mat C; |
---|
41 | //! Matrix D |
---|
42 | mat D; |
---|
43 | //! Matrix Q in square-root form |
---|
44 | sq_T Q; |
---|
45 | //! Matrix R in square-root form |
---|
46 | sq_T R; |
---|
47 | public: |
---|
48 | StateSpace() : A(), B(), C(), D(), Q(), R() {} |
---|
49 | //!copy constructor |
---|
50 | StateSpace ( const StateSpace<sq_T> &S0 ) : A ( S0.A ), B ( S0.B ), C ( S0.C ), D ( S0.D ), Q ( S0.Q ), R ( S0.R ) {} |
---|
51 | //! set all matrix parameters |
---|
52 | void set_parameters ( const mat &A0, const mat &B0, const mat &C0, const mat &D0, const sq_T &Q0, const sq_T &R0 ); |
---|
53 | //! validation |
---|
54 | void validate(); |
---|
55 | //! not virtual in this case |
---|
56 | void from_setting ( const Setting &set ) { |
---|
57 | UI::get ( A, set, "A", UI::compulsory ); |
---|
58 | UI::get ( B, set, "B", UI::compulsory ); |
---|
59 | UI::get ( C, set, "C", UI::compulsory ); |
---|
60 | UI::get ( D, set, "D", UI::compulsory ); |
---|
61 | mat Qtm, Rtm; // full matrices |
---|
62 | if ( !UI::get ( Qtm, set, "Q", UI::optional ) ) { |
---|
63 | vec dq; |
---|
64 | UI::get ( dq, set, "dQ", UI::compulsory ); |
---|
65 | Qtm = diag ( dq ); |
---|
66 | } |
---|
67 | if ( !UI::get ( Rtm, set, "R", UI::optional ) ) { |
---|
68 | vec dr; |
---|
69 | UI::get ( dr, set, "dQ", UI::compulsory ); |
---|
70 | Rtm = diag ( dr ); |
---|
71 | } |
---|
72 | R = Rtm; // automatic conversion to square-root form |
---|
73 | Q = Qtm; |
---|
74 | |
---|
75 | validate(); |
---|
76 | } |
---|
77 | void to_setting ( Setting &set ) const { |
---|
78 | UI::save( A, set, "A" ); |
---|
79 | UI::save( B, set, "B" ); |
---|
80 | UI::save( C, set, "C" ); |
---|
81 | UI::save( D, set, "D" ); |
---|
82 | UI::save( Q.to_mat(), set, "Q" ); |
---|
83 | UI::save( R.to_mat(), set, "R" ); |
---|
84 | } |
---|
85 | //! access function |
---|
86 | const mat& _A() const { |
---|
87 | return A; |
---|
88 | } |
---|
89 | //! access function |
---|
90 | const mat& _B() const { |
---|
91 | return B; |
---|
92 | } |
---|
93 | //! access function |
---|
94 | const mat& _C() const { |
---|
95 | return C; |
---|
96 | } |
---|
97 | //! access function |
---|
98 | const mat& _D() const { |
---|
99 | return D; |
---|
100 | } |
---|
101 | //! access function |
---|
102 | const sq_T& _Q() const { |
---|
103 | return Q; |
---|
104 | } |
---|
105 | //! access function |
---|
106 | const sq_T& _R() const { |
---|
107 | return R; |
---|
108 | } |
---|
109 | }; |
---|
110 | |
---|
111 | //! Common abstract base for Kalman filters |
---|
112 | template<class sq_T> |
---|
113 | class Kalman: public BM, public StateSpace<sq_T> { |
---|
114 | protected: |
---|
115 | //! id of output |
---|
116 | RV yrv; |
---|
117 | //! Kalman gain |
---|
118 | mat _K; |
---|
119 | //!posterior |
---|
120 | enorm<sq_T> est; |
---|
121 | //!marginal on data f(y|y) |
---|
122 | enorm<sq_T> fy; |
---|
123 | public: |
---|
124 | Kalman<sq_T>() : BM(), StateSpace<sq_T>(), yrv(), _K(), est() {} |
---|
125 | //! Copy constructor |
---|
126 | Kalman<sq_T> ( const Kalman<sq_T> &K0 ) : BM ( K0 ), StateSpace<sq_T> ( K0 ), yrv ( K0.yrv ), _K ( K0._K ), est ( K0.est ), fy ( K0.fy ) {} |
---|
127 | //!set statistics of the posterior |
---|
128 | void set_statistics ( const vec &mu0, const mat &P0 ) { |
---|
129 | est.set_parameters ( mu0, P0 ); |
---|
130 | }; |
---|
131 | //!set statistics of the posterior |
---|
132 | void set_statistics ( const vec &mu0, const sq_T &P0 ) { |
---|
133 | est.set_parameters ( mu0, P0 ); |
---|
134 | }; |
---|
135 | //! return correctly typed posterior (covariant return) |
---|
136 | const enorm<sq_T>& posterior() const { |
---|
137 | return est; |
---|
138 | } |
---|
139 | |
---|
140 | /*! Create object from the following structure |
---|
141 | |
---|
142 | \code |
---|
143 | class = 'KalmanFull'; |
---|
144 | prior = configuration of bdm::epdf; % prior density represented by any offspring of epdf, bdm::epdf::from_setting - it will be converted to gaussian |
---|
145 | --- inherited fields --- |
---|
146 | bdm::StateSpace<sq_T>::from_setting |
---|
147 | bdm::BM::from_setting |
---|
148 | \endcode |
---|
149 | */ |
---|
150 | void from_setting ( const Setting &set ) { |
---|
151 | StateSpace<sq_T>::from_setting ( set ); |
---|
152 | BM::from_setting(set); |
---|
153 | |
---|
154 | shared_ptr<epdf> pri=UI::build<epdf>(set,"prior",UI::compulsory); |
---|
155 | //bdm_assert(pri->dimension()==); |
---|
156 | set_statistics ( pri->mean(), pri->covariance() ); |
---|
157 | } |
---|
158 | void to_setting ( Setting &set ) const { |
---|
159 | StateSpace<sq_T>::to_setting ( set ); |
---|
160 | BM::to_setting(set); |
---|
161 | |
---|
162 | UI::save(est, set, "prior"); |
---|
163 | } |
---|
164 | //! validate object |
---|
165 | void validate() { |
---|
166 | StateSpace<sq_T>::validate(); |
---|
167 | dimy = this->C.rows(); |
---|
168 | dimc = this->B.cols(); |
---|
169 | set_dim ( this->A.rows() ); |
---|
170 | |
---|
171 | bdm_assert ( est.dimension(), "Statistics and model parameters mismatch" ); |
---|
172 | } |
---|
173 | |
---|
174 | }; |
---|
175 | /*! |
---|
176 | * \brief Basic Kalman filter with full matrices |
---|
177 | */ |
---|
178 | |
---|
179 | class KalmanFull : public Kalman<fsqmat> { |
---|
180 | public: |
---|
181 | //! For EKFfull; |
---|
182 | KalmanFull() : Kalman<fsqmat>() {}; |
---|
183 | //! Here dt = [yt;ut] of appropriate dimensions |
---|
184 | void bayes ( const vec &yt, const vec &cond = empty_vec ); |
---|
185 | |
---|
186 | virtual KalmanFull* _copy() const { |
---|
187 | KalmanFull* K = new KalmanFull; |
---|
188 | K->set_parameters ( A, B, C, D, Q, R ); |
---|
189 | K->set_statistics ( est._mu(), est._R() ); |
---|
190 | return K; |
---|
191 | } |
---|
192 | }; |
---|
193 | UIREGISTER ( KalmanFull ); |
---|
194 | |
---|
195 | |
---|
196 | /*! \brief Kalman filter in square root form |
---|
197 | |
---|
198 | Trivial example: |
---|
199 | \include kalman_simple.cpp |
---|
200 | |
---|
201 | Complete constructor: |
---|
202 | */ |
---|
203 | class KalmanCh : public Kalman<chmat> { |
---|
204 | protected: |
---|
205 | //! @{ \name Internal storage - needs initialize() |
---|
206 | //! pre array (triangular matrix) |
---|
207 | mat preA; |
---|
208 | //! post array (triangular matrix) |
---|
209 | mat postA; |
---|
210 | //!@} |
---|
211 | public: |
---|
212 | //! copy constructor |
---|
213 | virtual KalmanCh* _copy() const { |
---|
214 | KalmanCh* K = new KalmanCh; |
---|
215 | K->set_parameters ( A, B, C, D, Q, R ); |
---|
216 | K->set_statistics ( est._mu(), est._R() ); |
---|
217 | K->validate(); |
---|
218 | return K; |
---|
219 | } |
---|
220 | //! set parameters for adapt from Kalman |
---|
221 | void set_parameters ( const mat &A0, const mat &B0, const mat &C0, const mat &D0, const chmat &Q0, const chmat &R0 ); |
---|
222 | //! initialize internal parametetrs |
---|
223 | void initialize(); |
---|
224 | |
---|
225 | /*!\brief Here dt = [yt;ut] of appropriate dimensions |
---|
226 | |
---|
227 | The following equality hold::\f[ |
---|
228 | \left[\begin{array}{cc} |
---|
229 | R^{0.5}\\ |
---|
230 | P_{t|t-1}^{0.5}C' & P_{t|t-1}^{0.5}CA'\\ |
---|
231 | & Q^{0.5}\end{array}\right]<\mathrm{orth.oper.}>=\left[\begin{array}{cc} |
---|
232 | R_{y}^{0.5} & KA'\\ |
---|
233 | & P_{t+1|t}^{0.5}\\ |
---|
234 | \\\end{array}\right]\f] |
---|
235 | |
---|
236 | Thus this object evaluates only predictors! Not filtering densities. |
---|
237 | */ |
---|
238 | void bayes ( const vec &yt, const vec &cond = empty_vec ); |
---|
239 | |
---|
240 | /*! Create object from the following structure |
---|
241 | |
---|
242 | \code |
---|
243 | class = 'KalmanCh'; |
---|
244 | --- inherited fields --- |
---|
245 | bdm::Kalman<chmat>::from_setting |
---|
246 | \endcode |
---|
247 | */ |
---|
248 | void from_setting ( const Setting &set ) { |
---|
249 | Kalman<chmat>::from_setting ( set ); |
---|
250 | } |
---|
251 | |
---|
252 | void validate() { |
---|
253 | Kalman<chmat>::validate(); |
---|
254 | initialize(); |
---|
255 | } |
---|
256 | }; |
---|
257 | UIREGISTER ( KalmanCh ); |
---|
258 | |
---|
259 | /*! |
---|
260 | \brief Extended Kalman Filter in full matrices |
---|
261 | |
---|
262 | An approximation of the exact Bayesian filter with Gaussian noices and non-linear evolutions of their mean. |
---|
263 | */ |
---|
264 | class EKFfull : public KalmanFull { |
---|
265 | protected: |
---|
266 | //! Internal Model f(x,u) |
---|
267 | shared_ptr<diffbifn> pfxu; |
---|
268 | |
---|
269 | //! Observation Model h(x,u) |
---|
270 | shared_ptr<diffbifn> phxu; |
---|
271 | |
---|
272 | public: |
---|
273 | //! Default constructor |
---|
274 | EKFfull (); |
---|
275 | |
---|
276 | //! Set nonlinear functions for mean values and covariance matrices. |
---|
277 | void set_parameters ( const shared_ptr<diffbifn> &pfxu, const shared_ptr<diffbifn> &phxu, const mat Q0, const mat R0 ); |
---|
278 | |
---|
279 | //! Here dt = [yt;ut] of appropriate dimensions |
---|
280 | void bayes ( const vec &yt, const vec &cond = empty_vec ); |
---|
281 | //! set estimates |
---|
282 | void set_statistics ( const vec &mu0, const mat &P0 ) { |
---|
283 | est.set_parameters ( mu0, P0 ); |
---|
284 | }; |
---|
285 | //! access function |
---|
286 | const mat _R() { |
---|
287 | return est._R().to_mat(); |
---|
288 | } |
---|
289 | |
---|
290 | /*! Create object from the following structure |
---|
291 | |
---|
292 | \code |
---|
293 | class = 'EKFfull'; |
---|
294 | |
---|
295 | OM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
296 | IM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
297 | dQ = [...]; % vector containing diagonal of Q |
---|
298 | dR = [...]; % vector containing diagonal of R |
---|
299 | --- optional fields --- |
---|
300 | mu0 = [...]; % vector of statistics mu0 |
---|
301 | dP0 = [...]; % vector containing diagonal of P0 |
---|
302 | -- or -- |
---|
303 | P0 = [...]; % full matrix P0 |
---|
304 | --- inherited fields --- |
---|
305 | bdm::BM::from_setting |
---|
306 | \endcode |
---|
307 | If the optional fields are not given, they will be filled as follows: |
---|
308 | \code |
---|
309 | mu0 = [0,0,0,....]; % empty statistics |
---|
310 | P0 = eye( dim ); |
---|
311 | \endcode |
---|
312 | */ |
---|
313 | void from_setting ( const Setting &set ) { |
---|
314 | BM::from_setting ( set ); |
---|
315 | shared_ptr<diffbifn> IM = UI::build<diffbifn> ( set, "IM", UI::compulsory ); |
---|
316 | shared_ptr<diffbifn> OM = UI::build<diffbifn> ( set, "OM", UI::compulsory ); |
---|
317 | |
---|
318 | //statistics |
---|
319 | int dim = IM->dimension(); |
---|
320 | vec mu0; |
---|
321 | if ( !UI::get ( mu0, set, "mu0" ) ) |
---|
322 | mu0 = zeros ( dim ); |
---|
323 | |
---|
324 | mat P0; |
---|
325 | vec dP0; |
---|
326 | if ( UI::get ( dP0, set, "dP0" ) ) |
---|
327 | P0 = diag ( dP0 ); |
---|
328 | else if ( !UI::get ( P0, set, "P0" ) ) |
---|
329 | P0 = eye ( dim ); |
---|
330 | |
---|
331 | set_statistics ( mu0, P0 ); |
---|
332 | |
---|
333 | //parameters |
---|
334 | vec dQ, dR; |
---|
335 | UI::get ( dQ, set, "dQ", UI::compulsory ); |
---|
336 | UI::get ( dR, set, "dR", UI::compulsory ); |
---|
337 | set_parameters ( IM, OM, diag ( dQ ), diag ( dR ) ); |
---|
338 | |
---|
339 | // pfxu = UI::build<diffbifn>(set, "IM", UI::compulsory); |
---|
340 | // phxu = UI::build<diffbifn>(set, "OM", UI::compulsory); |
---|
341 | // |
---|
342 | // mat R0; |
---|
343 | // UI::get(R0, set, "R",UI::compulsory); |
---|
344 | // mat Q0; |
---|
345 | // UI::get(Q0, set, "Q",UI::compulsory); |
---|
346 | // |
---|
347 | // |
---|
348 | // mat P0; vec mu0; |
---|
349 | // UI::get(mu0, set, "mu0", UI::optional); |
---|
350 | // UI::get(P0, set, "P0", UI::optional); |
---|
351 | // set_statistics(mu0,P0); |
---|
352 | // // Initial values |
---|
353 | // UI::get (yrv, set, "yrv", UI::optional); |
---|
354 | // UI::get (urv, set, "urv", UI::optional); |
---|
355 | // set_drv(concat(yrv,urv)); |
---|
356 | // |
---|
357 | // // setup StateSpace |
---|
358 | // pfxu->dfdu_cond(mu0, zeros(pfxu->_dimu()), A,true); |
---|
359 | // phxu->dfdu_cond(mu0, zeros(pfxu->_dimu()), C,true); |
---|
360 | // |
---|
361 | } |
---|
362 | |
---|
363 | void validate() { |
---|
364 | // KalmanFull::validate(); |
---|
365 | |
---|
366 | // check stats and IM and OM |
---|
367 | } |
---|
368 | }; |
---|
369 | UIREGISTER ( EKFfull ); |
---|
370 | |
---|
371 | |
---|
372 | /*! |
---|
373 | \brief Extended Kalman Filter in Square root |
---|
374 | |
---|
375 | An approximation of the exact Bayesian filter with Gaussian noices and non-linear evolutions of their mean. |
---|
376 | */ |
---|
377 | |
---|
378 | class EKFCh : public KalmanCh { |
---|
379 | protected: |
---|
380 | //! Internal Model f(x,u) |
---|
381 | shared_ptr<diffbifn> pfxu; |
---|
382 | |
---|
383 | //! Observation Model h(x,u) |
---|
384 | shared_ptr<diffbifn> phxu; |
---|
385 | public: |
---|
386 | //! copy constructor duplicated - calls different set_parameters |
---|
387 | EKFCh* _copy() const { |
---|
388 | return new EKFCh(*this); |
---|
389 | } |
---|
390 | //! Set nonlinear functions for mean values and covariance matrices. |
---|
391 | void set_parameters ( const shared_ptr<diffbifn> &pfxu, const shared_ptr<diffbifn> &phxu, const chmat Q0, const chmat R0 ); |
---|
392 | |
---|
393 | //! Here dt = [yt;ut] of appropriate dimensions |
---|
394 | void bayes ( const vec &yt, const vec &cond = empty_vec ); |
---|
395 | |
---|
396 | /*! Create object from the following structure |
---|
397 | |
---|
398 | \code |
---|
399 | class = 'EKFCh'; |
---|
400 | OM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
401 | IM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
402 | dQ = [...]; % vector containing diagonal of Q |
---|
403 | dR = [...]; % vector containing diagonal of R |
---|
404 | --- optional fields --- |
---|
405 | mu0 = [...]; % vector of statistics mu0 |
---|
406 | dP0 = [...]; % vector containing diagonal of P0 |
---|
407 | -- or -- |
---|
408 | P0 = [...]; % full matrix P0 |
---|
409 | --- inherited fields --- |
---|
410 | bdm::BM::from_setting |
---|
411 | \endcode |
---|
412 | If the optional fields are not given, they will be filled as follows: |
---|
413 | \code |
---|
414 | mu0 = [0,0,0,....]; % empty statistics |
---|
415 | P0 = eye( dim ); |
---|
416 | \endcode |
---|
417 | */ |
---|
418 | void from_setting ( const Setting &set ); |
---|
419 | |
---|
420 | void validate() {}; |
---|
421 | // TODO dodelat void to_setting( Setting &set ) const; |
---|
422 | |
---|
423 | }; |
---|
424 | |
---|
425 | UIREGISTER ( EKFCh ); |
---|
426 | SHAREDPTR ( EKFCh ); |
---|
427 | |
---|
428 | //! EKF using bierman and Thorton code |
---|
429 | class EKF_UD : public BM { |
---|
430 | protected: |
---|
431 | //! logger |
---|
432 | LOG_LEVEL(EKF_UD,logU, logG, logD); |
---|
433 | //! Internal Model f(x,u) |
---|
434 | shared_ptr<diffbifn> pfxu; |
---|
435 | |
---|
436 | //! Observation Model h(x,u) |
---|
437 | shared_ptr<diffbifn> phxu; |
---|
438 | |
---|
439 | //! U part |
---|
440 | mat U; |
---|
441 | //! D part |
---|
442 | vec D; |
---|
443 | |
---|
444 | mat A; |
---|
445 | mat C; |
---|
446 | mat Q; |
---|
447 | vec R; |
---|
448 | |
---|
449 | enorm<ldmat> est; |
---|
450 | public: |
---|
451 | |
---|
452 | //! copy constructor duplicated |
---|
453 | EKF_UD* _copy() const { |
---|
454 | return new EKF_UD(*this); |
---|
455 | } |
---|
456 | |
---|
457 | const enorm<ldmat>& posterior()const{return est;}; |
---|
458 | |
---|
459 | enorm<ldmat>& prior() { |
---|
460 | return const_cast<enorm<ldmat>&>(posterior()); |
---|
461 | } |
---|
462 | |
---|
463 | EKF_UD(){} |
---|
464 | |
---|
465 | |
---|
466 | EKF_UD(const EKF_UD &E0): pfxu(E0.pfxu),phxu(E0.phxu), U(E0.U), D(E0.D){} |
---|
467 | |
---|
468 | //! Set nonlinear functions for mean values and covariance matrices. |
---|
469 | void set_parameters ( const shared_ptr<diffbifn> &pfxu, const shared_ptr<diffbifn> &phxu, const mat Q0, const vec R0 ); |
---|
470 | |
---|
471 | //! Here dt = [yt;ut] of appropriate dimensions |
---|
472 | void bayes ( const vec &yt, const vec &cond = empty_vec ); |
---|
473 | |
---|
474 | void log_register ( bdm::logger& L, const string& prefix ){ |
---|
475 | BM::log_register ( L, prefix ); |
---|
476 | |
---|
477 | if ( log_level[logU] ) |
---|
478 | L.add_vector ( log_level, logU, RV ( dimension()*dimension() ), prefix ); |
---|
479 | if ( log_level[logG] ) |
---|
480 | L.add_vector ( log_level, logG, RV ( dimension()*dimension() ), prefix ); |
---|
481 | if ( log_level[logD] ) |
---|
482 | L.add_vector ( log_level, logD, RV ( dimension()), prefix ); |
---|
483 | |
---|
484 | } |
---|
485 | /*! Create object from the following structure |
---|
486 | |
---|
487 | \code |
---|
488 | class = 'EKF_UD'; |
---|
489 | OM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
490 | IM = configuration of bdm::diffbifn; % any offspring of diffbifn, bdm::diffbifn::from_setting |
---|
491 | dQ = [...]; % vector containing diagonal of Q |
---|
492 | dR = [...]; % vector containing diagonal of R |
---|
493 | --- optional fields --- |
---|
494 | mu0 = [...]; % vector of statistics mu0 |
---|
495 | dP0 = [...]; % vector containing diagonal of P0 |
---|
496 | -- or -- |
---|
497 | P0 = [...]; % full matrix P0 |
---|
498 | --- inherited fields --- |
---|
499 | bdm::BM::from_setting |
---|
500 | \endcode |
---|
501 | If the optional fields are not given, they will be filled as follows: |
---|
502 | \code |
---|
503 | mu0 = [0,0,0,....]; % empty statistics |
---|
504 | P0 = eye( dim ); |
---|
505 | \endcode |
---|
506 | */ |
---|
507 | void from_setting ( const Setting &set ); |
---|
508 | |
---|
509 | void validate() {}; |
---|
510 | // TODO dodelat void to_setting( Setting &set ) const; |
---|
511 | |
---|
512 | }; |
---|
513 | UIREGISTER(EKF_UD); |
---|
514 | |
---|
515 | class UKFCh : public EKFCh{ |
---|
516 | public: |
---|
517 | double kappa; |
---|
518 | |
---|
519 | void bayes ( const vec &yt, const vec &cond = empty_vec ){ |
---|
520 | |
---|
521 | vec &_mu = est._mu(); |
---|
522 | chmat &_P = est._R(); |
---|
523 | chmat &_Ry = fy._R(); |
---|
524 | vec &_yp = fy._mu(); |
---|
525 | |
---|
526 | int dim = dimension(); |
---|
527 | int dim2 = 1+dim+dim; |
---|
528 | |
---|
529 | double npk =dim+kappa;//n+kappa |
---|
530 | mat Xi(dim,dim2); |
---|
531 | vec w=ones(dim2)* 0.5/npk; |
---|
532 | w(0) = (npk-dim)/npk; // mean is special |
---|
533 | |
---|
534 | //step 1. |
---|
535 | int i; |
---|
536 | Xi.set_col(0,_mu); |
---|
537 | |
---|
538 | for ( i=0;i<dim; i++){ |
---|
539 | vec tmp=sqrt(npk)*_P._Ch().get_col(i); |
---|
540 | Xi.set_col(i+1, _mu+tmp); |
---|
541 | Xi.set_col(i+1+dim, _mu-tmp); |
---|
542 | } |
---|
543 | |
---|
544 | // step 2. |
---|
545 | mat Xik(dim,dim2); |
---|
546 | for (i=0; i<dim2; i++){ |
---|
547 | Xik.set_col(i, pfxu->eval(Xi.get_col(i), cond)); |
---|
548 | } |
---|
549 | |
---|
550 | //step 3 |
---|
551 | vec xp=zeros(dim); |
---|
552 | for (i=0;i<dim2;i++){ |
---|
553 | xp += w(i) * Xik.get_col(i); |
---|
554 | } |
---|
555 | |
---|
556 | //step 4 |
---|
557 | mat P4=Q.to_mat(); |
---|
558 | vec tmp; |
---|
559 | for (i=0;i<dim2;i++){ |
---|
560 | tmp = Xik.get_col(i)-xp; |
---|
561 | P4+=w(i)*outer_product(tmp,tmp); |
---|
562 | } |
---|
563 | |
---|
564 | //step 5 |
---|
565 | mat Yi(dimy,dim2); |
---|
566 | for (i=0; i<dim2; i++){ |
---|
567 | Yi.set_col(i, phxu->eval(Xik.get_col(i), cond)); |
---|
568 | } |
---|
569 | //step 6 |
---|
570 | _yp.clear(); |
---|
571 | for (i=0;i<dim2;i++){ |
---|
572 | _yp += w(i) * Yi.get_col(i); |
---|
573 | } |
---|
574 | //step 7 |
---|
575 | mat Pvv=R.to_mat(); |
---|
576 | for (i=0;i<dim2;i++){ |
---|
577 | tmp = Yi.get_col(i)-_yp; |
---|
578 | Pvv+=w(i)*outer_product(tmp,tmp); |
---|
579 | } |
---|
580 | _Ry._Ch() = chol(Pvv); |
---|
581 | |
---|
582 | // step 8 |
---|
583 | mat Pxy=zeros(dim,dimy); |
---|
584 | for (i=0;i<dim2;i++){ |
---|
585 | Pxy+=w(i)*outer_product(Xi.get_col(i)-xp, Yi.get_col(i)-_yp); |
---|
586 | } |
---|
587 | mat iRy=inv(_Ry._Ch()); |
---|
588 | |
---|
589 | //filtering????? -- correction |
---|
590 | mat K=Pxy*iRy*iRy.T(); |
---|
591 | mat K2=Pxy*inv(_Ry.to_mat()); |
---|
592 | |
---|
593 | /////////////// new filtering density |
---|
594 | _mu = xp + K*(yt - _yp); |
---|
595 | |
---|
596 | if ( _mu ( 3 ) >pi ) _mu ( 3 )-=2*pi; |
---|
597 | if ( _mu ( 3 ) <-pi ) _mu ( 3 ) +=2*pi; |
---|
598 | // fill the space in Ppred; |
---|
599 | _P._Ch()=chol(P4-K*_Ry.to_mat()*K.T()); |
---|
600 | } |
---|
601 | void from_setting(const Setting &set){ |
---|
602 | EKFCh::from_setting(set); |
---|
603 | kappa = 1.0; |
---|
604 | UI::get(kappa,set,"kappa"); |
---|
605 | } |
---|
606 | }; |
---|
607 | UIREGISTER(UKFCh); |
---|
608 | |
---|
609 | //////// INstance |
---|
610 | |
---|
611 | /*! \brief (Switching) Multiple Model |
---|
612 | The model runs several models in parallel and evaluates thier weights (fittness). |
---|
613 | |
---|
614 | The statistics of the resulting density are merged using (geometric?) combination. |
---|
615 | |
---|
616 | The next step is performed with the new statistics for all models. |
---|
617 | */ |
---|
618 | class MultiModel: public BM { |
---|
619 | protected: |
---|
620 | //! List of models between which we switch |
---|
621 | Array<EKFCh*> Models; |
---|
622 | //! vector of model weights |
---|
623 | vec w; |
---|
624 | //! cache of model lls |
---|
625 | vec _lls; |
---|
626 | //! type of switching policy [1=maximum,2=...] |
---|
627 | int policy; |
---|
628 | //! internal statistics |
---|
629 | enorm<chmat> est; |
---|
630 | public: |
---|
631 | //! set internal parameters |
---|
632 | void set_parameters ( Array<EKFCh*> A, int pol0 = 1 ) { |
---|
633 | Models = A;//TODO: test if evalll is set |
---|
634 | w.set_length ( A.length() ); |
---|
635 | _lls.set_length ( A.length() ); |
---|
636 | policy = pol0; |
---|
637 | |
---|
638 | est.set_rv ( RV ( "MM", A ( 0 )->posterior().dimension(), 0 ) ); |
---|
639 | est.set_parameters ( A ( 0 )->posterior().mean(), A ( 0 )->posterior()._R() ); |
---|
640 | } |
---|
641 | void bayes ( const vec &yt, const vec &cond = empty_vec ) { |
---|
642 | int n = Models.length(); |
---|
643 | int i; |
---|
644 | for ( i = 0; i < n; i++ ) { |
---|
645 | Models ( i )->bayes ( yt ); |
---|
646 | _lls ( i ) = Models ( i )->_ll(); |
---|
647 | } |
---|
648 | double mlls = max ( _lls ); |
---|
649 | w = exp ( _lls - mlls ); |
---|
650 | w /= sum ( w ); //normalization |
---|
651 | //set statistics |
---|
652 | switch ( policy ) { |
---|
653 | case 1: { |
---|
654 | int mi = max_index ( w ); |
---|
655 | const enorm<chmat> &st = Models ( mi )->posterior() ; |
---|
656 | est.set_parameters ( st.mean(), st._R() ); |
---|
657 | } |
---|
658 | break; |
---|
659 | default: |
---|
660 | bdm_error ( "unknown policy" ); |
---|
661 | } |
---|
662 | // copy result to all models |
---|
663 | for ( i = 0; i < n; i++ ) { |
---|
664 | Models ( i )->set_statistics ( est.mean(), est._R() ); |
---|
665 | } |
---|
666 | } |
---|
667 | //! return correctly typed posterior (covariant return) |
---|
668 | const enorm<chmat>& posterior() const { |
---|
669 | return est; |
---|
670 | } |
---|
671 | |
---|
672 | void from_setting ( const Setting &set ); |
---|
673 | |
---|
674 | }; |
---|
675 | UIREGISTER ( MultiModel ); |
---|
676 | SHAREDPTR ( MultiModel ); |
---|
677 | |
---|
678 | //! conversion of outer ARX model (mlnorm) to state space model |
---|
679 | /*! |
---|
680 | The model is constructed as: |
---|
681 | \f[ x_{t+1} = Ax_t + B u_t + R^{1/2} e_t, y_t=Cx_t+Du_t + R^{1/2}w_t, \f] |
---|
682 | For example, for: |
---|
683 | Using Frobenius form, see []. |
---|
684 | |
---|
685 | For easier use in the future, indices theta_in_A and theta_in_C are set. TODO - explain |
---|
686 | */ |
---|
687 | //template<class sq_T> |
---|
688 | class StateCanonical: public StateSpace<fsqmat> { |
---|
689 | protected: |
---|
690 | //! remember connection from theta ->A |
---|
691 | datalink_part th2A; |
---|
692 | //! remember connection from theta ->C |
---|
693 | datalink_part th2C; |
---|
694 | //! remember connection from theta ->D |
---|
695 | datalink_part th2D; |
---|
696 | //!cached first row of A |
---|
697 | vec A1row; |
---|
698 | //!cached first row of C |
---|
699 | vec C1row; |
---|
700 | //!cached first row of D |
---|
701 | vec D1row; |
---|
702 | |
---|
703 | public: |
---|
704 | //! set up this object to match given mlnorm |
---|
705 | void connect_mlnorm ( const mlnorm<fsqmat> &ml ); |
---|
706 | |
---|
707 | //! fast function to update parameters from ml - not checked for compatibility!! |
---|
708 | void update_from ( const mlnorm<fsqmat> &ml ); |
---|
709 | }; |
---|
710 | /*! |
---|
711 | State-Space representation of multivariate autoregressive model. |
---|
712 | The original model: |
---|
713 | \f[ y_t = heta [\ldots y_{t-k}, \ldots u_{t-l}, \ldots z_{t-m}]' + \Sigma^{-1/2} e_t \f] |
---|
714 | where \f$ k,l,m \f$ are maximum delayes of corresponding variables in the regressor. |
---|
715 | |
---|
716 | The transformed state is: |
---|
717 | \f[ x_t = [y_{t} \ldots y_{t-k-1}, u_{t} \ldots u_{t-l-1}, z_{t} \ldots z_{t-m-1}]\f] |
---|
718 | |
---|
719 | The state accumulates all delayed values starting from time \f$ t \f$ . |
---|
720 | |
---|
721 | |
---|
722 | */ |
---|
723 | class StateFromARX: public StateSpace<chmat> { |
---|
724 | protected: |
---|
725 | //! remember connection from theta ->A |
---|
726 | datalink_part th2A; |
---|
727 | //! remember connection from theta ->B |
---|
728 | datalink_part th2B; |
---|
729 | //!function adds n diagonal elements from given starting point r,c |
---|
730 | void diagonal_part ( mat &A, int r, int c, int n ) { |
---|
731 | for ( int i = 0; i < n; i++ ) { |
---|
732 | A ( r, c ) = 1.0; |
---|
733 | r++; |
---|
734 | c++; |
---|
735 | } |
---|
736 | }; |
---|
737 | //! similar to ARX.have_constant |
---|
738 | bool have_constant; |
---|
739 | public: |
---|
740 | //! set up this object to match given mlnorm |
---|
741 | //! Note that state-space and common mpdf use different meaning of \f$ _t \f$ in \f$ u_t \f$. |
---|
742 | //!While mlnorm typically assumes that \f$ u_t \rightarrow y_t \f$ in state space it is \f$ u_{t-1} \rightarrow y_t \f$ |
---|
743 | //! For consequences in notation of internal variable xt see arx2statespace_notes.lyx. |
---|
744 | void connect_mlnorm ( const mlnorm<chmat> &ml, RV &xrv, RV &urv ); |
---|
745 | |
---|
746 | //! fast function to update parameters from ml - not checked for compatibility!! |
---|
747 | void update_from ( const mlnorm<chmat> &ml ); |
---|
748 | |
---|
749 | //! access function |
---|
750 | bool _have_constant() const { |
---|
751 | return have_constant; |
---|
752 | } |
---|
753 | }; |
---|
754 | |
---|
755 | /////////// INSTANTIATION |
---|
756 | |
---|
757 | template<class sq_T> |
---|
758 | void StateSpace<sq_T>::set_parameters ( const mat &A0, const mat &B0, const mat &C0, const mat &D0, const sq_T &Q0, const sq_T &R0 ) { |
---|
759 | |
---|
760 | A = A0; |
---|
761 | B = B0; |
---|
762 | C = C0; |
---|
763 | D = D0; |
---|
764 | R = R0; |
---|
765 | Q = Q0; |
---|
766 | validate(); |
---|
767 | } |
---|
768 | |
---|
769 | template<class sq_T> |
---|
770 | void StateSpace<sq_T>::validate() { |
---|
771 | bdm_assert ( A.cols() == A.rows(), "KalmanFull: A is not square" ); |
---|
772 | bdm_assert ( B.rows() == A.rows(), "KalmanFull: B is not compatible" ); |
---|
773 | bdm_assert ( C.cols() == A.rows(), "KalmanFull: C is not compatible" ); |
---|
774 | bdm_assert ( ( D.rows() == C.rows() ) && ( D.cols() == B.cols() ), "KalmanFull: D is not compatible" ); |
---|
775 | bdm_assert ( ( Q.cols() == A.rows() ) && ( Q.rows() == A.rows() ), "KalmanFull: Q is not compatible" ); |
---|
776 | bdm_assert ( ( R.cols() == C.rows() ) && ( R.rows() == C.rows() ), "KalmanFull: R is not compatible" ); |
---|
777 | } |
---|
778 | |
---|
779 | } |
---|
780 | #endif // KF_H |
---|
781 | |
---|