Logo ROOT   6.08/07
Reference Guide
testKalman.cxx
Go to the documentation of this file.
1 
2 
3 #include "Math/SVector.h"
4 #include "Math/SMatrix.h"
5 
6 #include "TMatrixD.h"
7 #include "TVectorD.h"
8 
9 #include "TRandom3.h"
10 
11 #include "matrix_util.h"
12 
13 
14 #define TEST_SYM
15 
16 //#define HAVE_CLHEP
17 #ifdef HAVE_CLHEP
18 #include "CLHEP/Matrix/SymMatrix.h"
19 #include "CLHEP/Matrix/Matrix.h"
20 #include "CLHEP/Matrix/Vector.h"
21 #endif
22 
23 // #include "SealUtil/SealTimer.h"
24 // #include "SealUtil/SealHRRTChrono.h"
25 // #include "SealUtil/TimingReport.h"
26 
27 #include <iostream>
28 
29 #ifndef NDIM1
30 #define NDIM1 2
31 #endif
32 #ifndef NDIM2
33 #define NDIM2 5
34 #endif
35 
36 #define NITER 1 // number of iterations
37 
38 #define NLOOP 1000000 // number of time the test is repeted
39 
40 using namespace ROOT::Math;
41 
42 #include "TestTimer.h"
43 
45 
46  // need to write explicitly the dimensions
47 
48 
49  typedef SMatrix<double, NDIM1, NDIM1> MnMatrixNN;
50  typedef SMatrix<double, NDIM2, NDIM2> MnMatrixMM;
51  typedef SMatrix<double, NDIM1, NDIM2> MnMatrixNM;
52  typedef SMatrix<double, NDIM2 , NDIM1> MnMatrixMN;
53  typedef SMatrix<double, NDIM1 > MnSymMatrixNN;
54  typedef SMatrix<double, NDIM2 > MnSymMatrixMM;
55  typedef SVector<double, NDIM1> MnVectorN;
56  typedef SVector<double, NDIM2> MnVectorM;
57 
58 
59 
60  int first = NDIM1; //Can change the size of the matrices
61  int second = NDIM2;
62 
63 
64  std::cout << "************************************************\n";
65  std::cout << " SMatrix kalman test " << first << " x " << second << std::endl;
66  std::cout << "************************************************\n";
67 
68 
69 
70 
71  int npass = NITER;
72  TRandom3 r(111);
73  double x2sum = 0, c2sum = 0;
74  for (int k = 0; k < npass; k++) {
75 
76 
77 
78  MnMatrixNM H;
79  MnMatrixMN K0;
80  MnSymMatrixMM Cp;
81  MnSymMatrixNN V;
82  MnVectorN m;
83  MnVectorM xp;
84 
85 
86  {
87  // fill matrices with random data
88  fillRandomMat(r,H,first,second);
89  fillRandomMat(r,K0,second,first);
90  fillRandomSym(r,Cp,second);
91  fillRandomSym(r,V,first);
92  fillRandomVec(r,m,first);
93  fillRandomVec(r,xp,second);
94  }
95 
96 
97  MnSymMatrixMM I;
98  for(int i = 0; i < second; i++)
99  I(i,i) = 1;
100 
101 #ifdef DEBUG
102  std::cout << "pass " << k << std::endl;
103  if (k == 0) {
104  std::cout << " K0 = " << K0 << std::endl;
105  std::cout << " H = " << K0 << std::endl;
106  std::cout << " Cp = " << Cp << std::endl;
107  std::cout << " V = " << V << std::endl;
108  std::cout << " m = " << m << std::endl;
109  std::cout << " xp = " << xp << std::endl;
110  }
111 #endif
112 
113 
114  {
115  double x2 = 0,c2 = 0;
116  test::Timer t("SMatrix Kalman ");
117 
118  MnVectorM x;
119  MnMatrixMN tmp;
120  MnSymMatrixNN Rinv;
121  MnMatrixMN K;
122  MnSymMatrixMM C;
123  MnVectorN vtmp1;
124  MnVectorN vtmp;
125 
126  for (int l = 0; l < NLOOP; l++)
127  {
128 
129 
130 
131  vtmp1 = H*xp -m;
132  //x = xp + K0 * (m- H * xp);
133  x = xp - K0 * vtmp1;
134  tmp = Cp * Transpose(H);
135  Rinv = V; Rinv += H * tmp;
136 
137  bool test = Rinv.Invert();
138  if(!test) {
139  std::cout<<"inversion failed" <<std::endl;
140  std::cout << Rinv << std::endl;
141  }
142 
143  K = tmp * Rinv ;
144  C = Cp; C -= K * Transpose(tmp);
145  //C = ( I - K * H ) * Cp;
146  //x2 = Product(Rinv,m-H*xp); // this does not compile on WIN32
147  vtmp = m-H*xp;
148  x2 = Dot(vtmp, Rinv*vtmp);
149 
150  }
151  //std::cout << k << " chi2 = " << x2 << std::endl;
152  x2sum += x2;
153  c2 = 0;
154  for (int i=0; i<NDIM2; ++i)
155  for (int j=0; j<NDIM2; ++j)
156  c2 += C(i,j);
157  c2sum += c2;
158  }
159  }
160  //tr.dump();
161 
162  std::cerr << "SMatrix: x2sum = " << x2sum << "\tc2sum = " << c2sum << std::endl;
163 
164  return 0;
165 }
166 
167 #ifdef TEST_SYM
169 
170  // need to write explicitly the dimensions
171 
172 
173  typedef SMatrix<double, NDIM1, NDIM1> MnMatrixNN;
174  typedef SMatrix<double, NDIM2, NDIM2> MnMatrixMM;
175  typedef SMatrix<double, NDIM1, NDIM2> MnMatrixNM;
176  typedef SMatrix<double, NDIM2 , NDIM1> MnMatrixMN;
179  typedef SVector<double, NDIM1> MnVectorN;
180  typedef SVector<double, NDIM2> MnVectorM;
181  typedef SVector<double, NDIM1*(NDIM1+1)/2> MnVectorN2;
182  typedef SVector<double, NDIM2*(NDIM2+1)/2> MnVectorM2;
183 
184 
185 
186  int first = NDIM1; //Can change the size of the matrices
187  int second = NDIM2;
188 
189 
190  std::cout << "************************************************\n";
191  std::cout << " SMatrix_SyM kalman test " << first << " x " << second << std::endl;
192  std::cout << "************************************************\n";
193 
194 
195 
196 
197  int npass = NITER;
198  TRandom3 r(111);
199  double x2sum = 0, c2sum = 0;
200  for (int k = 0; k < npass; k++) {
201 
202 
203 
204  MnMatrixNM H;
205  MnMatrixMN K0;
206  MnSymMatrixMM Cp;
207  MnSymMatrixNN V;
208  MnVectorN m;
209  MnVectorM xp;
210 
211 
212  {
213  // fill matrices with random data
214  fillRandomMat(r,H,first,second);
215  fillRandomMat(r,K0,second,first);
216  fillRandomSym(r,Cp,second);
217  fillRandomSym(r,V,first);
218  fillRandomVec(r,m,first);
219  fillRandomVec(r,xp,second);
220  }
221 
222 
223  MnSymMatrixMM I;
224  for(int i = 0; i < second; i++)
225  I(i,i) = 1;
226 
227 #ifdef DEBUG
228  std::cout << "pass " << k << std::endl;
229  if (k == 0) {
230  std::cout << " K0 = " << K0 << std::endl;
231  std::cout << " H = " << K0 << std::endl;
232  std::cout << " Cp = " << Cp << std::endl;
233  std::cout << " V = " << V << std::endl;
234  std::cout << " m = " << m << std::endl;
235  std::cout << " xp = " << xp << std::endl;
236  }
237 #endif
238 
239 
240  {
241  double x2 = 0,c2 = 0;
242  test::Timer t("SMatrix Kalman ");
243 
244  MnVectorM x;
245  MnSymMatrixNN RinvSym;
246  MnMatrixMN K;
247  MnSymMatrixMM C;
248  MnSymMatrixMM Ctmp;
249  MnVectorN vtmp1;
250  MnVectorN vtmp;
251 #define OPTIMIZED_SMATRIX_SYM
252 #ifdef OPTIMIZED_SMATRIX_SYM
253  MnMatrixMN tmp;
254 #endif
255 
256  for (int l = 0; l < NLOOP; l++)
257  {
258 
259 
260 #ifdef OPTIMIZED_SMATRIX_SYM
261  vtmp1 = H*xp -m;
262  //x = xp + K0 * (m- H * xp);
263  x = xp - K0 * vtmp1;
264  tmp = Cp * Transpose(H);
265  // we are sure that H*tmp result is symmetric
266  AssignSym::Evaluate(RinvSym,H*tmp);
267  RinvSym += V;
268 
269  bool test = RinvSym.Invert();
270  if(!test) {
271  std::cout<<"inversion failed" <<std::endl;
272  std::cout << RinvSym << std::endl;
273  }
274 
275  K = tmp * RinvSym ;
276  // we profit from the fact that result of K*tmpT is symmetric
277  AssignSym::Evaluate(Ctmp, K*Transpose(tmp) );
278  C = Cp; C -= Ctmp;
279  //C = ( I - K * H ) * Cp;
280  //x2 = Product(Rinv,m-H*xp); // this does not compile on WIN32
281  vtmp = m-H*xp;
282  x2 = Dot(vtmp, RinvSym*vtmp);
283 #else
284  // use similarity function
285  vtmp1 = H*xp -m;
286  x = xp - K0 * vtmp1;
287  RinvSym = V; RinvSym += Similarity(H,Cp);
288 
289  bool test = RinvSym.Invert();
290  if(!test) {
291  std::cout<<"inversion failed" <<std::endl;
292  std::cout << RinvSym << std::endl;
293  }
294 
295  Ctmp = SimilarityT(H, RinvSym);
296  C = Cp; C -= Similarity(Cp, Ctmp);
297  vtmp = m-H*xp;
298  x2 = Similarity(vtmp, RinvSym);
299 #endif
300 
301  }
302  //std::cout << k << " chi2 = " << x2 << std::endl;
303  x2sum += x2;
304  c2 = 0;
305  for (int i=0; i<NDIM2; ++i)
306  for (int j=0; j<NDIM2; ++j)
307  c2 += C(i,j);
308  c2sum += c2;
309  }
310  }
311  //tr.dump();
312 
313  std::cerr << "SMatrixSym: x2sum = " << x2sum << "\tc2sum = " << c2sum << std::endl;
314 
315  return 0;
316 }
317 
318 #endif
319 
320 
321 // ROOT test
322 
323 
325 
326 
327 
328 
329  typedef TMatrixD MnMatrix;
330  typedef TVectorD MnVector;
331 
332  // typedef boost::numeric::ublas::matrix<double> MnMatrix;
333  //typedef HepSymMatrix MnSymMatrixHep;
334 
335 
336  int first = NDIM1; //Can change the size of the matrices
337  int second = NDIM2;
338 
339 
340  std::cout << "************************************************\n";
341  std::cout << " TMatrix Kalman test " << first << " x " << second << std::endl;
342  std::cout << "************************************************\n";
343 
344 
345 
346  int npass = NITER;
347  TRandom3 r(111);
348  gMatrixCheck = 0;
349  double x2sum = 0,c2sum = 0;
350 
351  for (int k = 0; k < npass; k++)
352  {
353 
354  MnMatrix H(first,second);
355  MnMatrix K0(second,first);
356  MnMatrix Cp(second,second);
357  MnMatrix V(first,first);
358  MnVector m(first);
359  MnVector xp(second);
360 
361  // fill matrices with random data
362  fillRandomMat(r,H,first,second);
363  fillRandomMat(r,K0,second,first);
364  fillRandomSym(r,Cp,second);
365  fillRandomSym(r,V,first);
366  fillRandomVec(r,m,first);
367  fillRandomVec(r,xp,second);
368 
369 
370  MnMatrix I(second,second);//Identity matrix
371  for (int i = 0; i < second; i++)
372  for(int j = 0; j <second; j++) {
373  I(i,j) = 0.0;
374  if (i==j) I(i,i) = 1.0;
375  }
376 
377  // if (k==0) {
378  // std::cout << " Cp " << std::endl;
379  // Cp.Print();
380  // }
381 
382  {
383  double x2 = 0,c2 = 0;
384  TVectorD x(second);
385  TMatrixD Rinv(first,first);
386  TMatrixD Rtmp(first,first);
387  TMatrixDSym RinvSym;
388  TMatrixD K(second,first);
389  TMatrixD C(second,second);
390  TMatrixD Ctmp(second,second);
391  TVectorD tmp1(first);
392  TMatrixD tmp2(second,first);
393 #define OPTIMIZED_TMATRIX
394 #ifndef OPTIMIZED_TMATRIX
395  TMatrixD HT(second,first);
396  TMatrixD tmp2T(first,second);
397 #endif
398 
399  test::Timer t("TMatrix Kalman ");
400  for (Int_t l = 0; l < NLOOP; l++)
401  {
402 #ifdef OPTIMIZED_TMATRIX
403  tmp1 = m; Add(tmp1,-1.0,H,xp);
404  x = xp; Add(x,+1.0,K0,tmp1);
405  tmp2.MultT(Cp,H);
406  Rtmp.Mult(H,tmp2);
407  Rinv.Plus(V,Rtmp);
408  RinvSym.Use(first,Rinv.GetMatrixArray());
409  RinvSym.InvertFast();
410  K.Mult(tmp2,Rinv);
411  Ctmp.MultT(K,tmp2);
412  C.Minus(Cp,Ctmp);
413  x2 = RinvSym.Similarity(tmp1);
414 
415 #else
416  tmp1 = H*xp -m;
417  //x = xp + K0 * (m- H * xp);
418  x = xp - K0 * tmp1;
419  tmp2 = Cp * HT.Transpose(H);
420  Rinv = V; Rinv += H * tmp2;
421  RinvSym.Use(first,Rinv.GetMatrixArray());
422  RinvSym.InvertFast();
423  K= tmp2* Rinv;
424  C = Cp; C -= K*tmp2T.Transpose(tmp2);
425  x2= RinvSym.Similarity(tmp1);
426 #endif
427 
428  }
429  x2sum += x2;
430  c2 = 0;
431  for (int i=0; i<NDIM2; ++i)
432  for (int j=0; j<NDIM2; ++j)
433  c2 += C(i,j);
434  c2sum += c2;
435  }
436 
437  // }
438  }
439  //tr.dump();
440  std::cerr << "TMatrix: x2sum = " << x2sum << "\tc2sum = " << c2sum << std::endl;
441 
442  return 0;
443 }
444 
445 
446 // test CLHEP Kalman
447 
448 #ifdef HAVE_CLHEP
449 int test_clhep_kalman() {
450 
451 
452 
453  typedef HepSymMatrix MnSymMatrix;
454  typedef HepMatrix MnMatrix;
455  typedef HepVector MnVector;
456 
457 
458  // typedef boost::numeric::ublas::matrix<double> MnMatrix;
459  //typedef HepSymMatrix MnSymMatrixHep;
460 
461 
462  int first = NDIM1; //Can change the size of the matrices
463  int second = NDIM2;
464 
465 
466  std::cout << "************************************************\n";
467  std::cout << " CLHEP Kalman test " << first << " x " << second << std::endl;
468  std::cout << "************************************************\n";
469 
470 
471 
472  int npass = NITER;
473  TRandom3 r(111);
474  double x2sum = 0,c2sum = 0;
475 
476  for (int k = 0; k < npass; k++)
477  {
478 
479  // in CLHEP index starts from 1
480  MnMatrix H(first,second);
481  MnMatrix K0(second,first);
482  MnMatrix Cp(second,second);
483  MnMatrix V(first,first);
484  MnVector m(first);
485  MnVector xp(second);
486 
487  // fill matrices with random data
488  fillRandomMat(r,H,first,second,1);
489  fillRandomMat(r,K0,second,first,1);
490  fillRandomSym(r,Cp,second,1);
491  fillRandomSym(r,V,first,1);
492  fillRandomVec(r,m,first);
493  fillRandomVec(r,xp,second);
494 
495  MnSymMatrix I(second,1);//Identity matrix
496 
497  {
498  double x2 = 0,c2 = 0;
499  MnVector x(second);
500  MnMatrix Rinv(first,first);
501  MnSymMatrix RinvSym(first);
502  MnMatrix K(second,first);
503  MnSymMatrix C(second);
504  MnVector vtmp1(first);
505  MnMatrix tmp(second,first);
506 
507  test::Timer t("CLHEP Kalman ");
508  int ifail;
509  for (Int_t l = 0; l < NLOOP; l++)
510  {
511 
512 
513  vtmp1 = H*xp -m;
514  //x = xp + K0 * (m- H * xp);
515  x = xp - K0 * vtmp1;
516  tmp = Cp * H.T();
517  Rinv = V; Rinv += H * tmp;
518  RinvSym.assign(Rinv);
519  RinvSym.invert(ifail);
520  if (ifail !=0) { std::cout << "Error inverting Rinv" << std::endl; break; }
521  K = tmp*RinvSym;
522  //C.assign( (I-K*H)*Cp);
523  //C = (I-K*H)*Cp;
524  C.assign( (I-K*H)*Cp );
525  x2= RinvSym.similarity(vtmp1);
526  if(ifail!=0) { std::cout << "Error inverting Rinv" << std::endl; break; }
527  }
528  // std::cout << k << " chi2 " << x2 << std::endl;
529  x2sum += x2;
530 
531  c2 = 0;
532  for (int i=1; i<=NDIM2; ++i)
533  for (int j=1; j<=NDIM2; ++j)
534  c2 += C(i,j);
535  c2sum += c2;
536  }
537 
538  // }
539  }
540  //tr.dump();
541  std::cerr << "x2sum = " << x2sum << "\tc2sum = " << c2sum << std::endl;
542 
543  return 0;
544 }
545 #endif
546 
547 
548 
549 int testKalman() {
550 
551 #ifdef TEST_SYM
553 #endif
554 
557 #ifdef HAVE_CLHEP
558  test_clhep_kalman();
559 #endif
560 
561  return 0;
562 
563 
564 }
565 
566 int main() {
567  return testKalman();
568 }
T Dot(const SVector< T, D > &lhs, const SVector< T, D > &rhs)
Vector dot product.
Definition: Functions.h:166
T Similarity(const SMatrix< T, D, D, R > &lhs, const SVector< T, D > &rhs)
Similarity Vector - Matrix Product: v^T * A * v returning a scalar value of type T ...
int test_tmatrix_kalman()
Definition: testKalman.cxx:324
const int npass
Definition: testPermute.cxx:21
Random number generator class based on M.
Definition: TRandom3.h:29
void fillRandomSym(TRandom &r, M &m, unsigned int first, unsigned int start=0, double offset=1)
Definition: matrix_util.h:20
int test_smatrix_kalman()
Definition: testKalman.cxx:44
TMatrixT< Element > & Transpose(const TMatrixT< Element > &source)
Transpose matrix source.
Definition: TMatrixT.cxx:1469
TMatrixTSym< Element > & Use(Int_t row_lwb, Int_t row_upb, Element *data)
virtual const Element * GetMatrixArray() const
Definition: TMatrixT.h:223
#define H(x, y, z)
Definition: test.py:1
Expr< TransposeOp< SMatrix< T, D1, D2, R >, T, D1, D2 >, T, D2, D1, typename TranspPolicy< T, D1, D2, R >::RepType > Transpose(const SMatrix< T, D1, D2, R > &rhs)
Matrix Transpose B(i,j) = A(j,i) returning a matrix expression.
int Int_t
Definition: RtypesCore.h:41
#define NDIM2
Definition: testKalman.cxx:33
void Minus(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
General matrix summation. Create a matrix C such that C = A - B.
Definition: TMatrixT.cxx:580
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
SMatrix: a generic fixed size D1 x D2 Matrix class.
void MultT(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
General matrix multiplication. Create a matrix C such that C = A * B^T.
Definition: TMatrixT.cxx:951
void Add(THist< DIMENSIONS, PRECISION_TO, STAT_TO... > &to, THist< DIMENSIONS, PRECISION_FROM, STAT_FROM... > &from)
Add two histograms.
Definition: THist.hxx:327
#define NLOOP
Definition: testKalman.cxx:38
RooCmdArg Timer(Bool_t flag=kTRUE)
static double C[]
TRandom2 r(17)
TMarker * m
Definition: textangle.C:8
int test_smatrix_sym_kalman()
Definition: testKalman.cxx:168
Double_t K()
Definition: TMath.h:95
TLine * l
Definition: textangle.C:4
int testKalman()
Definition: testKalman.cxx:549
R__EXTERN Int_t gMatrixCheck
Definition: TMatrixTBase.h:91
static void Evaluate(SMatrix< T, D, D, MatRepSym< T, D > > &lhs, const Expr< A, T, D, D, R > &rhs)
assign a symmetric matrix from an expression
Definition: HelperOps.h:156
return c2
Definition: legend2.C:14
void fillRandomMat(TRandom &r, M &m, unsigned int first, unsigned int second, unsigned int start=0, double offset=1)
Definition: matrix_util.h:13
TMatrixTSym< Element > & InvertFast(Double_t *det=0)
Invert the matrix and calculate its determinant.
void fillRandomVec(TRandom &r, V &v, unsigned int len, unsigned int start=0, double offset=1)
Definition: matrix_util.h:7
TMatrixTSym< Element > & Similarity(const TMatrixT< Element > &n)
Calculate B * (*this) * B^T , final matrix will be (nrowsb x nrowsb) This is a similarity transform w...
void Plus(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
General matrix summation. Create a matrix C such that C = A + B.
Definition: TMatrixT.cxx:512
SMatrix< T, D2, D2, MatRepSym< T, D2 > > SimilarityT(const SMatrix< T, D1, D2, R > &lhs, const SMatrix< T, D1, D1, MatRepSym< T, D1 > > &rhs)
Transpose Similarity Matrix Product : B = U^T * A * U for A symmetric returning a symmetric matrix ex...
#define NDIM1
Definition: testKalman.cxx:30
Definition: first.py:1
void Mult(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
General matrix multiplication. Create a matrix C such that C = A * B.
Definition: TMatrixT.cxx:648
#define I(x, y, z)
#define NITER
Definition: testKalman.cxx:36
int main()
Definition: testKalman.cxx:566
SVector: a generic fixed size Vector class.