00001
00002 #define RUN_ALL_POINTS
00003
00004
00005
00006 #include "Math/SVector.h"
00007 #include "Math/SMatrix.h"
00008
00009 #include "TMatrixD.h"
00010 #include "TVectorD.h"
00011
00012 #include "TFile.h"
00013 #include "TSystem.h"
00014 #include "TROOT.h"
00015
00016 #include "TRandom3.h"
00017 #include "TStopwatch.h"
00018
00019 #include <iostream>
00020
00021 #include <map>
00022
00023
00024
00025 #ifdef HAVE_CLHEP
00026 #include "CLHEP/Matrix/SymMatrix.h"
00027 #include "CLHEP/Matrix/Matrix.h"
00028 #include "CLHEP/Matrix/Vector.h"
00029 #endif
00030
00031
00032 #define NITER 1 // number of iterations
00033
00034 #define NLOOP 1000000 // number of time the test is repeted
00035
00036
00037
00038 template <unsigned int NDIM1, unsigned int NDIM2>
00039 class TestRunner {
00040
00041 public:
00042
00043
00044 int test_smatrix_kalman() ;
00045
00046 int test_smatrix_sym_kalman();
00047
00048
00049 int test_tmatrix_kalman();
00050
00051 #ifdef HAVE_CLHEP
00052
00053 int test_clhep_kalman();
00054 #endif
00055
00056
00057 int run() {
00058 int iret = 0;
00059 iret |= test_smatrix_sym_kalman();
00060 iret |= test_smatrix_kalman();
00061 iret |= test_tmatrix_kalman();
00062 #ifdef HAVE_CLHEP
00063 iret |= test_clhep_kalman();
00064 #endif
00065 std::cout << "\n\n";
00066
00067 return iret;
00068 }
00069
00070 private:
00071 double fX2sum;
00072 double fC2sum;
00073
00074 };
00075
00076
00077 #define TR(N1,N2) \
00078 { TestRunner<N1,N2> tr; if (tr.run()) return -1; }
00079
00080
00081 int runTest() {
00082
00083 #ifndef RUN_ALL_POINTS
00084 TR(2,5)
00085
00086 #else
00087 TR(2,2)
00088 TR(2,3)
00089 TR(2,4)
00090 TR(2,5)
00091 TR(2,6)
00092 TR(2,7)
00093 TR(2,8)
00094 TR(3,2)
00095 TR(3,3)
00096 TR(3,4)
00097 TR(3,5)
00098 TR(3,6)
00099 TR(3,7)
00100 TR(3,8)
00101 TR(4,2)
00102 TR(4,3)
00103 TR(4,4)
00104 TR(4,5)
00105 TR(4,6)
00106 TR(4,7)
00107 TR(4,8)
00108 TR(5,2)
00109 TR(5,3)
00110 TR(5,4)
00111 TR(5,5)
00112 TR(5,6)
00113 TR(5,7)
00114 TR(5,8)
00115 TR(6,2)
00116 TR(6,3)
00117 TR(6,4)
00118 TR(6,5)
00119 TR(6,6)
00120 TR(6,7)
00121 TR(6,8)
00122 TR(7,2)
00123 TR(7,3)
00124 TR(7,4)
00125 TR(7,5)
00126 TR(7,6)
00127 TR(7,7)
00128 TR(7,8)
00129 TR(8,2)
00130 TR(8,3)
00131 TR(8,4)
00132 TR(8,5)
00133 TR(8,6)
00134 TR(8,7)
00135 TR(8,8)
00136 TR(9,2)
00137 TR(9,3)
00138 TR(9,4)
00139 TR(9,5)
00140 TR(9,6)
00141 TR(9,7)
00142 TR(9,8)
00143 TR(10,2)
00144 TR(10,3)
00145 TR(10,4)
00146 TR(10,5)
00147 TR(10,6)
00148 TR(10,7)
00149 TR(10,8)
00150 #endif
00151 return 0;
00152 }
00153
00154
00155
00156
00157
00158 template<class V>
00159 void fillRandomVec(TRandom & r, V & v, unsigned int len, unsigned int start = 0, double offset = 1) {
00160 for(unsigned int i = start; i < len+start; ++i)
00161 v[i] = r.Rndm() + offset;
00162 }
00163
00164 template<class M>
00165 void fillRandomMat(TRandom & r, M & m, unsigned int first, unsigned int second, unsigned int start = 0, double offset = 1) {
00166 for(unsigned int i = start; i < first+start; ++i)
00167 for(unsigned int j = start; j < second+start; ++j)
00168 m(i,j) = r.Rndm() + offset;
00169 }
00170
00171 template<class M>
00172 void fillRandomSym(TRandom & r, M & m, unsigned int first, unsigned int start = 0, double offset = 1) {
00173 for(unsigned int i = start; i < first+start; ++i) {
00174 for(unsigned int j = i; j < first+start; ++j) {
00175 if ( i != j ) {
00176 m(i,j) = r.Rndm() + offset;
00177 m(j,i) = m(i,j);
00178 }
00179 else
00180 m(i,i) = r.Rndm() + 3*offset;
00181 }
00182 }
00183 }
00184
00185
00186
00187
00188
00189
00190 void printTime(TStopwatch & time, std::string s) {
00191 int pr = std::cout.precision(4);
00192 std::cout << std::setw(12) << s << "\t" << " Real time = " << time.RealTime() << "\t(sec)\tCPU time = "
00193 << time.CpuTime() << "\t(sec)"
00194 << std::endl;
00195 std::cout.precision(pr);
00196 }
00197
00198
00199 double refTime1[4] = { 40.49, 53.75, 83.21,1000 };
00200 double refTime2[4] = { 393.81, 462.16, 785.50,10000 };
00201
00202 #define NMAX1 9 // matrix storese results from 2 to 10
00203 #define NMAX2 7 // results from 2 to 8
00204
00205 class TimeReport {
00206
00207
00208 typedef std::map<std::string, double> a;
00209 typedef ROOT::Math::SMatrix<double,NMAX1,NMAX2> M;
00210 typedef std::map< std::string, M > ResultTable;
00211
00212 public:
00213
00214 TimeReport() {}
00215
00216 ~TimeReport() { }
00217
00218
00219 void Set(const std::string & name, int dim1, int dim2 ) {
00220 fDim1 = dim1;
00221 fDim2 = dim2;
00222 fName = name;
00223
00224 if (fResult1.find(name) == fResult1.end() )
00225 fResult1.insert(ResultTable::value_type(name, M() ) );
00226 if (fResult2.find(name) == fResult2.end() )
00227 fResult2.insert(ResultTable::value_type(name, M() ) );
00228
00229 }
00230
00231 std::string name() const { return fName; }
00232
00233 void report(double rt, double ct) {
00234 fResult1[fName](fDim1-2,fDim2-2) = rt;
00235 fResult2[fName](fDim1-2,fDim2-2) = ct;
00236 }
00237
00238 double smallSum(const M & m, int cut = 6) {
00239
00240 double sum = 0;
00241 for (int i = 0; i<cut-1; ++i)
00242 for (int j = 0; j<cut-1; ++j)
00243 sum+= m(i,j);
00244
00245 return sum;
00246 }
00247
00248 double largeSum(const M & m, int cut = 6) {
00249
00250 double sum = 0;
00251 for (int i = 0; i<M::kRows; ++i)
00252 for (int j = 0; j<M::kCols; ++j)
00253 if ( i > cut-2 || j > cut-2)
00254 sum+= m(i,j);
00255
00256 return sum;
00257 }
00258
00259 void print(std::ostream & os) {
00260 std::map<std::string,double> r1;
00261 std::map<std::string,double> r2;
00262 os << "Real time results " << std::endl;
00263 for (ResultTable::iterator itr = fResult1.begin(); itr != fResult1.end(); ++itr) {
00264 std::string type = itr->first;
00265 os << " Results for " << type << std::endl;
00266 os << "\n" << itr->second << std::endl << std::endl;
00267 r1[type] = smallSum(itr->second);
00268 r2[type] = largeSum(itr->second);
00269 os << "\nTotal for N1,N2 <= 6 : " << r1[type] << std::endl;
00270 os << "\nTotal for N1,N2 > 6 : " << r2[type] << std::endl;
00271 }
00272 os << "\n\nCPU time results " << std::endl;
00273 for (ResultTable::iterator itr = fResult2.begin(); itr != fResult2.end(); ++itr) {
00274 os << " Results for " << itr->first << std::endl;
00275 os << "\n" << itr->second << std::endl << std::endl;
00276 os << "\nTotal for N1,N2 <= 6 : " << smallSum(itr->second) << std::endl;
00277 os << "\nTotal for N1,N2 > 6 : " << largeSum(itr->second) << std::endl;
00278 }
00279
00280
00281 os << "\n\n****************************************************************************\n";
00282 os << "Root version: " << gROOT->GetVersion() << " "
00283 << gROOT->GetVersionDate() << "/" << gROOT->GetVersionTime() << std::endl;
00284 os << "****************************************************************************\n";
00285 os << "\n\t ROOTMARKS for N1,N2 <= 6 \n\n";
00286 int j = 0;
00287 os.setf(std::ios::right,std::ios::adjustfield);
00288 for (std::map<std::string,double>::iterator i = r1.begin(); i != r1.end(); ++i) {
00289 std::string type = i->first;
00290 os << std::setw(12) << type << "\t=\t" << refTime1[j]*800/r1[type] << std::endl;
00291 j++;
00292 }
00293 os << "\n\t ROOTMARKS for N1,N2 > 6 \n\n";
00294 j = 0;
00295 for (std::map<std::string,double>::iterator i = r1.begin(); i != r1.end(); ++i) {
00296 std::string type = i->first;
00297 os << std::setw(12) << type << "\t=\t" << refTime2[j]*800/r2[type] << std::endl;
00298 j++;
00299 }
00300
00301 }
00302
00303 void save(const std::string & fileName) {
00304 TFile file(fileName.c_str(),"RECREATE");
00305 gSystem->Load("libSmatrix");
00306
00307
00308 for (ResultTable::iterator itr = fResult1.begin(); itr != fResult1.end(); ++itr) {
00309 int ret = file.WriteObject(&(itr->second),(itr->first).c_str() );
00310 if (ret ==0) std::cerr << "==> Error saving results in ROOT file " << fileName << std::endl;
00311 }
00312
00313
00314 for (ResultTable::iterator itr = fResult2.begin(); itr != fResult2.end(); ++itr) {
00315 std::string typeName = itr->first + "_2";
00316 int ret = file.WriteObject(&(itr->second), typeName.c_str() );
00317 if (ret ==0) std::cerr << "==> Error saving results in ROOT file " << fileName << std::endl;
00318 }
00319 file.Close();
00320 }
00321
00322 private:
00323
00324 int fDim1;
00325 int fDim2;
00326 std::string fName;
00327
00328 ResultTable fResult1;
00329 ResultTable fResult2;
00330
00331 };
00332
00333
00334 TimeReport gReporter;
00335
00336 class TestTimer {
00337
00338 public:
00339
00340
00341
00342
00343
00344
00345 TestTimer(TimeReport & r ) :
00346 fTime(0), fRep(&r)
00347 {
00348 fName = fRep->name();
00349 fWatch.Start();
00350 }
00351 TestTimer(double & t, const std::string & s = "") : fName(s), fTime(&t), fRep(0)
00352 {
00353 fWatch.Start();
00354 }
00355
00356 ~TestTimer() {
00357 fWatch.Stop();
00358 printTime(fWatch,fName);
00359 if (fRep) fRep->report( fWatch.RealTime(), fWatch.CpuTime() );
00360 if (fTime) *fTime += fWatch.RealTime();
00361 }
00362
00363
00364 private:
00365
00366 std::string fName;
00367 double * fTime;
00368 TStopwatch fWatch;
00369 TimeReport * fRep;
00370
00371 };
00372
00373
00374
00375 template <unsigned int NDIM1, unsigned int NDIM2>
00376 int TestRunner<NDIM1,NDIM2>::test_smatrix_kalman() {
00377
00378 gReporter.Set("SMatrix",NDIM1,NDIM2);
00379
00380
00381
00382
00383 typedef ROOT::Math::SMatrix<double, NDIM1, NDIM1> MnMatrixNN;
00384 typedef ROOT::Math::SMatrix<double, NDIM2, NDIM2> MnMatrixMM;
00385 typedef ROOT::Math::SMatrix<double, NDIM1, NDIM2> MnMatrixNM;
00386 typedef ROOT::Math::SMatrix<double, NDIM2 , NDIM1> MnMatrixMN;
00387 typedef ROOT::Math::SMatrix<double, NDIM1 > MnSymMatrixNN;
00388 typedef ROOT::Math::SMatrix<double, NDIM2 > MnSymMatrixMM;
00389 typedef ROOT::Math::SVector<double, NDIM1> MnVectorN;
00390 typedef ROOT::Math::SVector<double, NDIM2> MnVectorM;
00391
00392
00393
00394 int first = NDIM1;
00395 int second = NDIM2;
00396
00397
00398 std::cout << "****************************************************************************\n";
00399 std::cout << "\t\tSMatrix kalman test " << first << " x " << second << std::endl;
00400 std::cout << "****************************************************************************\n";
00401
00402
00403
00404
00405 int npass = NITER;
00406 TRandom3 r(111);
00407 double x2sum = 0, c2sum = 0;
00408 for (int k = 0; k < npass; k++) {
00409
00410
00411
00412 MnMatrixNM H;
00413 MnMatrixMN K0;
00414 MnSymMatrixMM Cp;
00415 MnSymMatrixNN V;
00416 MnVectorN m;
00417 MnVectorM xp;
00418
00419
00420 {
00421
00422 fillRandomMat(r,H,first,second);
00423 fillRandomMat(r,K0,second,first);
00424 fillRandomSym(r,Cp,second);
00425 fillRandomSym(r,V,first);
00426 fillRandomVec(r,m,first);
00427 fillRandomVec(r,xp,second);
00428 }
00429
00430
00431 MnSymMatrixMM I;
00432 for(int i = 0; i < second; i++)
00433 I(i,i) = 1;
00434
00435 #ifdef DEBUG
00436 std::cout << "pass " << k << std::endl;
00437 if (k == 0) {
00438 std::cout << " K0 = " << K0 << std::endl;
00439 std::cout << " H = " << K0 << std::endl;
00440 std::cout << " Cp = " << Cp << std::endl;
00441 std::cout << " V = " << V << std::endl;
00442 std::cout << " m = " << m << std::endl;
00443 std::cout << " xp = " << xp << std::endl;
00444 }
00445 #endif
00446
00447
00448 {
00449 double x2 = 0,c2 = 0;
00450 TestTimer t(gReporter);
00451
00452 MnVectorM x;
00453 MnMatrixMN tmp;
00454 MnSymMatrixNN Rinv;
00455 MnMatrixMN K;
00456 MnSymMatrixMM C;
00457 MnVectorN vtmp1;
00458 MnVectorN vtmp;
00459
00460 for (int l = 0; l < NLOOP; l++)
00461 {
00462
00463
00464
00465 vtmp1 = H*xp -m;
00466
00467 x = xp - K0 * vtmp1;
00468 tmp = Cp * Transpose(H);
00469 Rinv = V; Rinv += H * tmp;
00470
00471 bool test = Rinv.Invert();
00472 if(!test) {
00473 std::cout<<"inversion failed" <<std::endl;
00474 std::cout << Rinv << std::endl;
00475 }
00476
00477 K = tmp * Rinv ;
00478 C = Cp; C -= K * Transpose(tmp);
00479
00480
00481 vtmp = m-H*xp;
00482 x2 = Dot(vtmp, Rinv*vtmp);
00483
00484 #ifdef DEBUG
00485 if (l == 0) {
00486 std::cout << " Rinv =\n " << Rinv << std::endl;
00487 std::cout << " C =\n " << C << std::endl;
00488 }
00489 #endif
00490
00491 }
00492
00493 x2sum += x2;
00494 c2 = 0;
00495 for (unsigned int i=0; i<NDIM2; ++i)
00496 for (unsigned int j=0; j<NDIM2; ++j)
00497 c2 += C(i,j);
00498 c2sum += c2;
00499 }
00500 }
00501
00502
00503 int iret = 0;
00504 double d = std::abs(x2sum-fX2sum);
00505 if ( d > 1.E-6 * fX2sum ) {
00506 std::cout << "ERROR: difference found in x2sum = " << x2sum << "\tref = " << fX2sum <<
00507 "\tdiff = " << d << std::endl;
00508 iret = 1;
00509 }
00510 d = std::abs(c2sum-fC2sum);
00511 if ( d > 1.E-6 * fC2sum ) {
00512 std::cout << "ERROR: difference found in c2sum = " << c2sum << "\tref = " << fC2sum <<
00513 "\tdiff = " << d << std::endl;
00514 iret = 1;
00515 }
00516
00517 return iret;
00518 }
00519
00520 template <unsigned int NDIM1, unsigned int NDIM2>
00521 int TestRunner<NDIM1,NDIM2>::test_smatrix_sym_kalman() {
00522
00523 gReporter.Set("SMatrix_sym",NDIM1,NDIM2);
00524
00525
00526
00527
00528 typedef ROOT::Math::SMatrix<double, NDIM1, NDIM1> MnMatrixNN;
00529 typedef ROOT::Math::SMatrix<double, NDIM2, NDIM2> MnMatrixMM;
00530 typedef ROOT::Math::SMatrix<double, NDIM1, NDIM2> MnMatrixNM;
00531 typedef ROOT::Math::SMatrix<double, NDIM2 , NDIM1> MnMatrixMN;
00532 typedef ROOT::Math::SMatrix<double, NDIM1, NDIM1, ROOT::Math::MatRepSym<double, NDIM1> > MnSymMatrixNN;
00533 typedef ROOT::Math::SMatrix<double, NDIM2, NDIM2, ROOT::Math::MatRepSym<double, NDIM2> > MnSymMatrixMM;
00534 typedef ROOT::Math::SVector<double, NDIM1> MnVectorN;
00535 typedef ROOT::Math::SVector<double, NDIM2> MnVectorM;
00536 typedef ROOT::Math::SVector<double, NDIM1*(NDIM1+1)/2> MnVectorN2;
00537 typedef ROOT::Math::SVector<double, NDIM2*(NDIM2+1)/2> MnVectorM2;
00538
00539
00540
00541 int first = NDIM1;
00542 int second = NDIM2;
00543
00544
00545 std::cout << "****************************************************************************\n";
00546 std::cout << "\t\tSMatrix_Sym kalman test " << first << " x " << second << std::endl;
00547 std::cout << "****************************************************************************\n";
00548
00549
00550
00551
00552 int npass = NITER;
00553 TRandom3 r(111);
00554 double x2sum = 0, c2sum = 0;
00555 for (int k = 0; k < npass; k++) {
00556
00557
00558
00559 MnMatrixNM H;
00560 MnMatrixMN K0;
00561 MnSymMatrixMM Cp;
00562 MnSymMatrixNN V;
00563 MnVectorN m;
00564 MnVectorM xp;
00565
00566
00567 {
00568
00569 fillRandomMat(r,H,first,second);
00570 fillRandomMat(r,K0,second,first);
00571 fillRandomSym(r,Cp,second);
00572 fillRandomSym(r,V,first);
00573 fillRandomVec(r,m,first);
00574 fillRandomVec(r,xp,second);
00575 }
00576
00577
00578 MnSymMatrixMM I;
00579 for(int i = 0; i < second; i++)
00580 I(i,i) = 1;
00581
00582 #ifdef DEBUG
00583 std::cout << "pass " << k << std::endl;
00584 if (k == 0) {
00585 std::cout << " K0 = " << K0 << std::endl;
00586 std::cout << " H = " << K0 << std::endl;
00587 std::cout << " Cp = " << Cp << std::endl;
00588 std::cout << " V = " << V << std::endl;
00589 std::cout << " m = " << m << std::endl;
00590 std::cout << " xp = " << xp << std::endl;
00591 }
00592
00593 #endif
00594
00595 {
00596 double x2 = 0,c2 = 0;
00597 TestTimer t(gReporter);
00598
00599 MnVectorM x;
00600 MnSymMatrixNN RinvSym;
00601 MnMatrixMN K;
00602 MnSymMatrixMM C;
00603 MnSymMatrixMM Ctmp;
00604 MnVectorN vtmp1;
00605 MnVectorN vtmp;
00606 #define OPTIMIZED_SMATRIX_SYM
00607 #ifdef OPTIMIZED_SMATRIX_SYM
00608 MnMatrixMN tmp;
00609 #endif
00610
00611 for (int l = 0; l < NLOOP; l++)
00612 {
00613
00614
00615 #ifdef OPTIMIZED_SMATRIX_SYM
00616 vtmp1 = H*xp -m;
00617
00618 x = xp - K0 * vtmp1;
00619 tmp = Cp * Transpose(H);
00620
00621 ROOT::Math::AssignSym::Evaluate(RinvSym,H*tmp);
00622 RinvSym += V;
00623
00624 bool test = RinvSym.Invert();
00625 if(!test) {
00626 std::cout<<"inversion failed" <<std::endl;
00627 std::cout << RinvSym << std::endl;
00628 }
00629
00630 K = tmp * RinvSym ;
00631
00632 ROOT::Math::AssignSym::Evaluate(Ctmp, K*Transpose(tmp) );
00633 C = Cp; C -= Ctmp;
00634
00635
00636 vtmp = m-H*xp;
00637 x2 = ROOT::Math::Dot(vtmp, RinvSym*vtmp);
00638 #else
00639
00640 vtmp1 = H*xp -m;
00641 x = xp - K0 * vtmp1;
00642 RinvSym = V; RinvSym += Similarity(H,Cp);
00643
00644 bool test = RinvSym.Invert();
00645 if(!test) {
00646 std::cout<<"inversion failed" <<std::endl;
00647 std::cout << RinvSym << std::endl;
00648 }
00649
00650 Ctmp = ROOT::Math::SimilarityT(H, RinvSym);
00651 C = Cp; C -= ROOT::Math::Similarity(Cp, Ctmp);
00652 vtmp = m-H*xp;
00653 x2 = ROOT::Math::Similarity(vtmp, RinvSym);
00654 #endif
00655
00656 }
00657
00658 x2sum += x2;
00659 c2 = 0;
00660 for (unsigned int i=0; i<NDIM2; ++i)
00661 for (unsigned int j=0; j<NDIM2; ++j)
00662 c2 += C(i,j);
00663 c2sum += c2;
00664 }
00665 }
00666
00667
00668 fX2sum = x2sum;
00669 fC2sum = c2sum;
00670 if (x2sum == 0 || c2sum == 0) {
00671 std::cout << "WARNING: x2sum = " << x2sum << "\tc2sum = " << c2sum << std::endl;
00672 }
00673
00674 return 0;
00675 }
00676
00677
00678
00679
00680
00681
00682 template <unsigned int NDIM1, unsigned int NDIM2>
00683 int TestRunner<NDIM1,NDIM2>::test_tmatrix_kalman() {
00684
00685 gReporter.Set("TMatrix",NDIM1,NDIM2);
00686
00687
00688 typedef TMatrixD MnMatrix;
00689 typedef TVectorD MnVector;
00690
00691
00692
00693
00694
00695 int first = NDIM1;
00696 int second = NDIM2;
00697
00698
00699 std::cout << "****************************************************************************\n";
00700 std::cout << "\t\tTMatrix Kalman test " << first << " x " << second << std::endl;
00701 std::cout << "****************************************************************************\n";
00702
00703 int npass = NITER;
00704 TRandom3 r(111);
00705 double x2sum = 0,c2sum = 0;
00706
00707 for (int k = 0; k < npass; k++)
00708 {
00709
00710 MnMatrix H(first,second);
00711 MnMatrix K0(second,first);
00712 MnMatrix Cp(second,second);
00713 MnMatrix V(first,first);
00714 MnVector m(first);
00715 MnVector xp(second);
00716
00717
00718 fillRandomMat(r,H,first,second);
00719 fillRandomMat(r,K0,second,first);
00720 fillRandomSym(r,Cp,second);
00721 fillRandomSym(r,V,first);
00722 fillRandomVec(r,m,first);
00723 fillRandomVec(r,xp,second);
00724
00725
00726 MnMatrix I(second,second);
00727 for (int i = 0; i < second; i++)
00728 for(int j = 0; j <second; j++) {
00729 I(i,j) = 0.0;
00730 if (i==j) I(i,i) = 1.0;
00731 }
00732
00733
00734
00735
00736
00737
00738 {
00739 double x2 = 0,c2 = 0;
00740 TVectorD x(second);
00741 TMatrixD Rtmp(first,first);
00742 TMatrixD Rinv(first,first);
00743 TMatrixDSym RinvSym;
00744 TMatrixD K(second,first);
00745 TMatrixD C(second,second);
00746 TMatrixD Ctmp(second,second);
00747 TVectorD tmp1(first);
00748 TMatrixD tmp2(second,first);
00749
00750 TestTimer t(gReporter);
00751 for (Int_t l = 0; l < NLOOP; l++)
00752 {
00753 tmp1 = m; Add(tmp1,-1.0,H,xp);
00754 x = xp; Add(x,+1.0,K0,tmp1);
00755 tmp2.MultT(Cp,H);
00756 Rtmp.Mult(H,tmp2);
00757 Rinv.Plus(V,Rtmp);
00758 RinvSym.Use(first,Rinv.GetMatrixArray());
00759 RinvSym.InvertFast();
00760 K.Mult(tmp2,Rinv);
00761 Ctmp.MultT(K,tmp2);
00762 C.Minus(Cp,Ctmp);
00763 x2 = RinvSym.Similarity(tmp1);
00764
00765 #ifdef DEBUG
00766 if (l == 0) {
00767 std::cout << " Rinv =\n "; Rinv.Print();
00768 std::cout << " RinvSym =\n "; RinvSym.Print();
00769 std::cout << " C =\n "; C.Print();
00770 }
00771 #endif
00772
00773 }
00774 x2sum += x2;
00775 c2 = 0;
00776 for (unsigned int i=0; i<NDIM2; ++i)
00777 for (unsigned int j=0; j<NDIM2; ++j)
00778 c2 += C(i,j);
00779 c2sum += c2;
00780 }
00781
00782
00783 }
00784
00785
00786
00787
00788
00789 int iret = 0;
00790 double d = std::abs(x2sum-fX2sum);
00791 if ( d > 1.E-6 * fX2sum ) {
00792 std::cout << "ERROR: difference found in x2sum = " << x2sum << "\tref = " << fX2sum <<
00793 "\tdiff = " << d << std::endl;
00794 iret = 1;
00795 }
00796 d = std::abs(c2sum-fC2sum);
00797 if ( d > 1.E-6 * fC2sum ) {
00798 std::cout << "ERROR: difference found in c2sum = " << c2sum << "\tref = " << fC2sum <<
00799 "\tdiff = " << d << std::endl;
00800 iret = 1;
00801 }
00802
00803
00804 return iret;
00805 }
00806
00807
00808
00809
00810 #ifdef HAVE_CLHEP
00811 template <unsigned int NDIM1, unsigned int NDIM2>
00812 int TestRunner<NDIM1,NDIM2>::test_clhep_kalman() {
00813
00814
00815 gReporter.Set("HepMatrix",NDIM1,NDIM2);
00816
00817 typedef HepSymMatrix MnSymMatrix;
00818 typedef HepMatrix MnMatrix;
00819 typedef HepVector MnVector;
00820
00821
00822
00823
00824
00825
00826 int first = NDIM1;
00827 int second = NDIM2;
00828
00829
00830 std::cout << "****************************************************************************\n";
00831 std::cout << " CLHEP Kalman test " << first << " x " << second << std::endl;
00832 std::cout << "****************************************************************************\n";
00833
00834 int npass = NITER;
00835 TRandom3 r(111);
00836 double x2sum = 0,c2sum = 0;
00837
00838 for (int k = 0; k < npass; k++)
00839 {
00840
00841
00842 MnMatrix H(first,second);
00843 MnMatrix K0(second,first);
00844 MnMatrix Cp(second,second);
00845 MnMatrix V(first,first);
00846 MnVector m(first);
00847 MnVector xp(second);
00848
00849
00850 fillRandomMat(r,H,first,second,1);
00851 fillRandomMat(r,K0,second,first,1);
00852 fillRandomSym(r,Cp,second,1);
00853 fillRandomSym(r,V,first,1);
00854 fillRandomVec(r,m,first);
00855 fillRandomVec(r,xp,second);
00856
00857 MnSymMatrix I(second,1);
00858
00859 {
00860 double x2 = 0,c2 = 0;
00861 MnVector x(second);
00862 MnMatrix Rinv(first,first);
00863 MnSymMatrix RinvSym(first);
00864 MnMatrix K(second,first);
00865 MnSymMatrix C(second);
00866 MnVector vtmp1(first);
00867 MnMatrix tmp(second,first);
00868
00869 TestTimer t(gReporter);
00870 int ifail;
00871 for (Int_t l = 0; l < NLOOP; l++)
00872 {
00873
00874
00875 vtmp1 = H*xp -m;
00876
00877 x = xp - K0 * vtmp1;
00878 tmp = Cp * H.T();
00879 Rinv = V; Rinv += H * tmp;
00880 RinvSym.assign(Rinv);
00881 RinvSym.invert(ifail);
00882 if (ifail !=0) { std::cout << "Error inverting Rinv" << std::endl; break; }
00883 K = tmp*RinvSym;
00884
00885
00886 C.assign( (I-K*H)*Cp );
00887 x2= RinvSym.similarity(vtmp1);
00888 if(ifail!=0) { std::cout << "Error inverting Rinv" << std::endl; break; }
00889 }
00890 x2sum += x2;
00891
00892 c2 = 0;
00893 for (unsigned int i=1; i<=NDIM2; ++i)
00894 for (unsigned int j=1; j<=NDIM2; ++j)
00895 c2 += C(i,j);
00896 c2sum += c2;
00897 }
00898
00899
00900 }
00901
00902
00903
00904 int iret = 0;
00905 double d = std::abs(x2sum-fX2sum);
00906 if ( d > 1.E-6 * fX2sum ) {
00907 std::cout << "ERROR: difference found in x2sum = " << x2sum << "\tref = " << fX2sum <<
00908 "\tdiff = " << d << std::endl;
00909 iret = 1;
00910 }
00911 d = std::abs(c2sum-fC2sum);
00912 if ( d > 1.E-6 * fC2sum ) {
00913 std::cout << "ERROR: difference found in c2sum = " << c2sum << "\tref = " << fC2sum <<
00914 "\tdiff = " << d << std::endl;
00915 iret = 1;
00916 }
00917
00918 return iret;
00919 }
00920 #endif
00921
00922
00923 int main(int argc, char *argv[]) {
00924
00925
00926 if (runTest() ) {
00927 std::cout << "\nERROR - stressKalman FAILED - exit!" << std::endl ;
00928 return -1;
00929 };
00930
00931 gReporter.print(std::cout);
00932 std::string fname = "kalman";
00933 if (argc > 1) {
00934 std::string platf(argv[1]);
00935 fname = fname + "_" + platf;
00936 }
00937
00938 gReporter.save(fname+".root");
00939
00940 return 0;
00941 }