00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "Math/SMatrix.h"
00011 #include "TMatrixD.h"
00012 #include "TMatrixDSym.h"
00013 #include "TRandom3.h"
00014 #include "TFile.h"
00015 #include "TTree.h"
00016 #include "TStopwatch.h"
00017 #include "TSystem.h"
00018
00019 #include <iostream>
00020
00021 #ifdef USE_REFLEX
00022 #include "Cintex/Cintex.h"
00023 #include "Reflex/Reflex.h"
00024 #endif
00025
00026 #include "Track.h"
00027
00028 TRandom3 R;
00029 TStopwatch timer;
00030
00031
00032
00033
00034
00035
00036 #ifdef USE_DOUBLE32
00037 typedef ROOT::Math::SMatrix<Double32_t,5,5,ROOT::Math::MatRepSym<Double32_t,5> > SMatrixSym5;
00038 typedef ROOT::Math::SMatrix<Double32_t,5,5 > SMatrix5;
00039 const std::string sname = "ROOT::Math::SMatrix<Double32_t,5,5,ROOT::Math::MatRepStd<Double32_t,5,5> >";
00040 const std::string sname_sym = "ROOT::Math::SMatrix<Double32_t,5,5,ROOT::Math::MatRepSym<Double32_t,5> >";
00041 double tol = 1.E-6;
00042 #else
00043 typedef ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> > SMatrixSym5;
00044 typedef ROOT::Math::SMatrix<double,5,5 > SMatrix5;
00045 const std::string sname = "ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepStd<double,5,5> >";
00046 const std::string sname_sym = "ROOT::Math::SMatrix<double,5,5,ROOT::Math::MatRepSym<double,5> >";
00047 double tol = 1.E-16;
00048 #endif
00049 double tol32 = 1.E-6;
00050
00051 #ifdef USE_REFLEX
00052 std::string sfile1 = "smatrix_rflx.root";
00053 std::string sfile2 = "smatrix.root";
00054
00055 std::string symfile1 = "smatrixsym_rflx.root";
00056 std::string symfile2 = "smatrixsym.root";
00057 #else
00058 std::string sfile1 = "smatrix.root";
00059 std::string sfile2 = "smatrix_rflx.root";
00060
00061 std::string symfile1 = "smatrixsym.root";
00062 std::string symfile2 = "smatrixsym_rflx.root";
00063 #endif
00064
00065
00066
00067
00068
00069
00070 template<class Matrix>
00071 void FillMatrix(Matrix & m) {
00072 for (int i = 0; i < 5; ++i) {
00073 for (int j = 0; j < 5; ++j) {
00074 m(i,j) = R.Rndm() + 1;
00075 }
00076 }
00077 }
00078
00079 void FillCArray(double * m) {
00080 for (int i = 0; i < 5; ++i) {
00081 for (int j = 0; j < 5; ++j) {
00082 m[i*5+j] = R.Rndm() + 1;
00083 }
00084 }
00085 }
00086
00087 template<class Matrix>
00088 void FillMatrixSym(Matrix & m) {
00089 for (int i = 0; i < 5; ++i) {
00090 for (int j = 0; j < 5; ++j) {
00091 if (j>=i) m(i,j) = R.Rndm() + 1;
00092 else m(i,j) = m(j,i);
00093 }
00094 }
00095 }
00096
00097 template<class R>
00098 double SumSMatrix(ROOT::Math::SMatrix<double,5,5,R> & m) {
00099 double sum = 0;
00100 for (int i = 0; i < 5*5; ++i) {
00101 sum += m.apply(i);
00102 }
00103 return sum;
00104 }
00105
00106 double SumCArray(double * m) {
00107 double sum = 0;
00108 for (int i = 0; i < 5*5; ++i) {
00109 sum += m[i];
00110 }
00111 return sum;
00112 }
00113
00114 template<class TM>
00115 double SumTMatrix(TM & m) {
00116 double sum = 0;
00117 const double * d = m.GetMatrixArray();
00118 for (int i = 0; i < 5*5; ++i) {
00119 sum += d[i];
00120 }
00121 return sum;
00122 }
00123
00124
00125
00126 void initMatrix(int n) {
00127
00128
00129
00130 timer.Start();
00131 SMatrix5 s;
00132 R.SetSeed(1);
00133 for (int i = 0; i < n; ++i) {
00134 FillMatrix(s);
00135 }
00136 timer.Stop();
00137 std::cout << " Time to fill SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00138
00139 timer.Start();
00140 SMatrixSym5 ss;
00141 R.SetSeed(1);
00142 for (int i = 0; i < n; ++i) {
00143 FillMatrixSym(ss);
00144 }
00145 timer.Stop();
00146 std::cout << " Time to fill SMatrix Sym " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00147
00148 timer.Start();
00149 TMatrixD t(5,5);
00150 R.SetSeed(1);
00151 for (int i = 0; i < n; ++i) {
00152 FillMatrix(t);
00153 }
00154 timer.Stop();
00155 std::cout << " Time to fill TMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00156
00157 timer.Start();
00158 TMatrixDSym ts(5);
00159 R.SetSeed(1);
00160 for (int i = 0; i < n; ++i) {
00161 FillMatrixSym(ts);
00162 }
00163 timer.Stop();
00164 std::cout << " Time to fill TMatrix Sym " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00165
00166 }
00167
00168
00169
00170 double writeCArray(int n) {
00171
00172 std::cout << "\n";
00173 std::cout << "**************************************************\n";
00174 std::cout << " Test writing a C Array ........\n";
00175 std::cout << "**************************************************\n";
00176
00177 TFile f1("cmatrix.root","RECREATE");
00178
00179
00180
00181 TTree t1("t1","Tree with C Array");
00182
00183 double m1[25];
00184 t1.Branch("C Array branch",m1,"m1[25]/D");
00185
00186 timer.Start();
00187 double etot = 0;
00188 R.SetSeed(1);
00189 for (int i = 0; i < n; ++i) {
00190 FillCArray(m1);
00191 etot += SumCArray(m1);
00192 t1.Fill();
00193 }
00194
00195 f1.Write();
00196 timer.Stop();
00197
00198 std::cout << " Time to Write CArray " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00199
00200 #ifdef DEBUG
00201 t1.Print();
00202 int pr = std::cout.precision(18);
00203 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00204 std::cout.precision(pr);
00205 #endif
00206
00207 return etot/double(n);
00208 }
00209
00210 double writeSMatrix(int n, const std::string & file) {
00211
00212 std::cout << "\n";
00213 std::cout << "**************************************************\n";
00214 std::cout << " Test writing SMatrix ........\n";
00215 std::cout << "**************************************************\n";
00216
00217 TFile f1(file.c_str(),"RECREATE");
00218
00219
00220 TTree t1("t1","Tree with SMatrix");
00221
00222 SMatrix5 * m1 = new SMatrix5;
00223
00224
00225 t1.Branch("SMatrix branch",sname.c_str(),&m1,16000,0);
00226
00227 timer.Start();
00228 double etot = 0;
00229 R.SetSeed(1);
00230 for (int i = 0; i < n; ++i) {
00231 FillMatrix(*m1);
00232 etot += SumSMatrix(*m1);
00233 t1.Fill();
00234 }
00235
00236 f1.Write();
00237 timer.Stop();
00238
00239 std::cout << " Time to Write SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00240 #ifdef DEBUG
00241 t1.Print();
00242 int pr = std::cout.precision(18);
00243 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00244 std::cout.precision(pr);
00245 #endif
00246
00247 return etot/double(n);
00248 }
00249
00250
00251
00252 double writeSMatrixSym(int n, const std::string & file) {
00253
00254 std::cout << "\n";
00255 std::cout << "**************************************************\n";
00256 std::cout << " Test writing SMatrix Sym.....\n";
00257 std::cout << "**************************************************\n";
00258
00259 TFile f1(file.c_str(),"RECREATE");
00260
00261
00262 TTree t1("t1","Tree with SMatrix");
00263
00264 SMatrixSym5 * m1 = new SMatrixSym5;
00265
00266 t1.Branch("SMatrixSym branch",sname_sym.c_str(),&m1,16000,0);
00267
00268 timer.Start();
00269 double etot = 0;
00270 R.SetSeed(1);
00271 for (int i = 0; i < n; ++i) {
00272 FillMatrixSym(*m1);
00273 etot += SumSMatrix(*m1);
00274 t1.Fill();
00275 }
00276
00277 f1.Write();
00278 timer.Stop();
00279
00280
00281 std::cout << " Time to Write SMatrix Sym " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00282 #ifdef DEBUG
00283 t1.Print();
00284 int pr = std::cout.precision(18);
00285 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00286 std::cout.precision(pr);
00287 #endif
00288
00289 return etot/double(n);
00290 }
00291
00292
00293
00294
00295
00296 double writeTMatrix(int n) {
00297
00298
00299 std::cout << "\n";
00300 std::cout << "**************************************************\n";
00301 std::cout << " Test writing TMatrix........\n";
00302 std::cout << "**************************************************\n";
00303
00304
00305 TFile f2("tmatrix.root","RECREATE");
00306 TTree t2("t2","Tree with TMatrix");
00307
00308 TMatrixD * m2 = new TMatrixD(5,5);
00309 TMatrixD::Class()->IgnoreTObjectStreamer();
00310
00311
00312 t2.Branch("TMatrix branch",&m2,16000,2);
00313
00314 double etot = 0;
00315 timer.Start();
00316 R.SetSeed(1);
00317 for (int i = 0; i < n; ++i) {
00318 FillMatrix(*m2);
00319 etot += SumTMatrix(*m2);
00320 t2.Fill();
00321 }
00322
00323 f2.Write();
00324 timer.Stop();
00325
00326 std::cout << " Time to Write TMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00327 #ifdef DEBUG
00328 t2.Print();
00329 int pr = std::cout.precision(18);
00330 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00331 std::cout.precision(pr);
00332 std::cout << "\n\n\n";
00333 #endif
00334
00335 return etot/double(n);
00336
00337 }
00338
00339 double writeTMatrixSym(int n) {
00340
00341
00342 std::cout << "\n";
00343 std::cout << "**************************************************\n";
00344 std::cout << " Test writing TMatrix.Sym....\n";
00345 std::cout << "**************************************************\n";
00346
00347
00348 TFile f2("tmatrixsym.root","RECREATE");
00349 TTree t2("t2","Tree with TMatrix");
00350
00351 TMatrixDSym * m2 = new TMatrixDSym(5);
00352 TMatrixDSym::Class()->IgnoreTObjectStreamer();
00353
00354
00355 t2.Branch("TMatrixSym branch",&m2,16000,0);
00356
00357 double etot = 0;
00358 timer.Start();
00359 R.SetSeed(1);
00360 for (int i = 0; i < n; ++i) {
00361 FillMatrixSym(*m2);
00362 etot += SumTMatrix(*m2);
00363 t2.Fill();
00364 }
00365
00366 f2.Write();
00367 timer.Stop();
00368
00369 std::cout << " Time to Write TMatrix Sym " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00370 #ifdef DEBUG
00371 t2.Print();
00372 int pr = std::cout.precision(18);
00373 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00374 std::cout.precision(pr);
00375 std::cout << "\n\n\n";
00376 #endif
00377
00378 return etot/double(n);
00379
00380 }
00381
00382
00383
00384
00385 double readTMatrix() {
00386
00387
00388
00389
00390 std::cout << "\n\n";
00391 std::cout << "**************************************************\n";
00392 std::cout << " Test reading TMatrix........\n";
00393 std::cout << "**************************************************\n";
00394
00395
00396 TFile f2("tmatrix.root");
00397 if (f2.IsZombie() ) {
00398 std::cerr << "Error opening the ROOT file" << std::endl;
00399 return -1;
00400 }
00401 TTree *t2 = (TTree*)f2.Get("t2");
00402
00403
00404 TMatrixD * v2 = 0;
00405 t2->SetBranchAddress("TMatrix branch",&v2);
00406
00407 timer.Start();
00408 int n = (int) t2->GetEntries();
00409 double etot = 0;
00410 for (int i = 0; i < n; ++i) {
00411 t2->GetEntry(i);
00412 etot += SumTMatrix(*v2);
00413 }
00414
00415 timer.Stop();
00416 std::cout << " Time for TMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00417 double val = etot/double(n);
00418 #ifdef DEBUG
00419 std::cout << " Tree Entries " << n << std::endl;
00420 int pr = std::cout.precision(18);
00421 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00422 std::cout.precision(pr);
00423 #endif
00424 return val;
00425 }
00426
00427
00428 double readTMatrixSym() {
00429
00430
00431
00432
00433 std::cout << "\n\n";
00434 std::cout << "**************************************************\n";
00435 std::cout << " Test reading TMatrix.Sym....\n";
00436 std::cout << "**************************************************\n";
00437
00438
00439 TFile f2("tmatrixsym.root");
00440 if (f2.IsZombie() ) {
00441 std::cerr << "Error opening the ROOT file" << std::endl;
00442 return -1;
00443 }
00444
00445
00446 TTree *t2 = (TTree*)f2.Get("t2");
00447
00448
00449 TMatrixDSym * v2 = 0;
00450 t2->SetBranchAddress("TMatrixSym branch",&v2);
00451
00452 timer.Start();
00453 int n = (int) t2->GetEntries();
00454 double etot = 0;
00455 for (int i = 0; i < n; ++i) {
00456 t2->GetEntry(i);
00457 etot += SumTMatrix(*v2);
00458 }
00459
00460 timer.Stop();
00461 std::cout << " Time for TMatrix Sym" << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00462 double val = etot/double(n);
00463 #ifdef DEBUG
00464 std::cout << " Tree Entries " << n << std::endl;
00465 int pr = std::cout.precision(18);
00466 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00467 std::cout.precision(pr);
00468 #endif
00469
00470 return val;
00471 }
00472
00473
00474 double readSMatrix(const std::string & file) {
00475
00476
00477 std::cout << "\n";
00478 std::cout << "**************************************************\n";
00479 std::cout << " Test reading SMatrix........\n";
00480 std::cout << "**************************************************\n";
00481
00482
00483 TFile f1(file.c_str());
00484 if (f1.IsZombie() ) {
00485 std::cerr << "Error opening the ROOT file" << file << std::endl;
00486 return -1;
00487 }
00488
00489
00490 TTree *t1 = (TTree*)f1.Get("t1");
00491
00492 SMatrix5 *v1 = 0;
00493 t1->SetBranchAddress("SMatrix branch",&v1);
00494
00495 timer.Start();
00496 int n = (int) t1->GetEntries();
00497 double etot=0;
00498 for (int i = 0; i < n; ++i) {
00499 t1->GetEntry(i);
00500 etot += SumSMatrix(*v1);
00501 }
00502
00503
00504 timer.Stop();
00505 std::cout << " Time for SMatrix : " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00506
00507 #ifdef DEBUG
00508 std::cout << " Tree Entries " << n << std::endl;
00509 int pr = std::cout.precision(18);
00510 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00511 std::cout.precision(pr);
00512 #endif
00513 std::cout << "\n";
00514
00515
00516 return etot/double(n);
00517 }
00518
00519
00520 double readSMatrixSym(const std::string & file) {
00521
00522
00523 std::cout << "\n";
00524 std::cout << "**************************************************\n";
00525 std::cout << " Test reading SMatrix.Sym....\n";
00526 std::cout << "**************************************************\n";
00527
00528
00529 TFile f1(file.c_str());
00530 if (f1.IsZombie() ) {
00531 std::cerr << "Error opening the ROOT file" << file << std::endl;
00532 return -1;
00533 }
00534
00535
00536
00537 TTree *t1 = (TTree*)f1.Get("t1");
00538
00539 SMatrixSym5 *v1 = 0;
00540 t1->SetBranchAddress("SMatrixSym branch",&v1);
00541
00542 timer.Start();
00543 int n = (int) t1->GetEntries();
00544 double etot=0;
00545 for (int i = 0; i < n; ++i) {
00546 t1->GetEntry(i);
00547 etot += SumSMatrix(*v1);
00548 }
00549
00550
00551 timer.Stop();
00552 std::cout << " Time for SMatrix Sym : " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00553
00554 #ifdef DEBUG
00555 std::cout << " Tree Entries " << n << std::endl;
00556 int pr = std::cout.precision(18);
00557 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00558 std::cout.precision(pr);
00559 #endif
00560 std::cout << "\n";
00561
00562 return etot/double(n);
00563 }
00564
00565
00566 double writeTrackD(int n) {
00567
00568 std::cout << "\n";
00569 std::cout << "**************************************************\n";
00570 std::cout << " Test writing Track class........\n";
00571 std::cout << "**************************************************\n";
00572
00573 TFile f1("track.root","RECREATE");
00574
00575
00576 TTree t1("t1","Tree with Track based on SMatrix");
00577
00578 TrackD * m1 = new TrackD();
00579
00580 t1.Branch("Track branch",&m1,16000,0);
00581
00582 timer.Start();
00583 double etot = 0;
00584 R.SetSeed(1);
00585 for (int i = 0; i < n; ++i) {
00586 FillMatrix(m1->CovMatrix());
00587 etot += SumSMatrix(m1->CovMatrix() );
00588 t1.Fill();
00589 }
00590
00591 f1.Write();
00592 timer.Stop();
00593
00594 std::cout << " Time to Write TrackD of SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00595
00596 #ifdef DEBUG
00597 t1.Print();
00598 int pr = std::cout.precision(18);
00599 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00600 std::cout.precision(pr);
00601 #endif
00602
00603 return etot/double(n);
00604 }
00605
00606
00607 double writeTrackD32(int n) {
00608
00609 std::cout << "\n";
00610 std::cout << "**************************************************\n";
00611 std::cout << " Test writing TrackD32 class........\n";
00612 std::cout << "**************************************************\n";
00613
00614 TFile f1("track32.root","RECREATE");
00615
00616
00617 TTree t1("t1","Tree with Track based on SMatrix");
00618
00619 TrackD32 * m1 = new TrackD32();
00620 t1.Branch("Track32 branch",&m1,16000,0);
00621
00622 timer.Start();
00623 double etot = 0;
00624 R.SetSeed(1);
00625 for (int i = 0; i < n; ++i) {
00626 FillMatrix(m1->CovMatrix());
00627 etot += SumSMatrix(m1->CovMatrix() );
00628 t1.Fill();
00629 }
00630
00631 f1.Write();
00632 timer.Stop();
00633
00634 std::cout << " Time to Write TrackD32 of SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00635
00636 #ifdef DEBUG
00637 t1.Print();
00638 int pr = std::cout.precision(18);
00639 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00640 std::cout.precision(pr);
00641 #endif
00642
00643 return etot/double(n);
00644 }
00645
00646
00647 double readTrackD() {
00648
00649 std::cout << "\n";
00650 std::cout << "**************************************************\n";
00651 std::cout << " Test reading Track class........\n";
00652 std::cout << "**************************************************\n";
00653
00654 TFile f1("track.root");
00655 if (f1.IsZombie() ) {
00656 std::cerr << "Error opening the ROOT file" << std::endl;
00657 return -1;
00658 }
00659
00660
00661 TTree *t1 = (TTree*)f1.Get("t1");
00662
00663 TrackD *trk = 0;
00664 t1->SetBranchAddress("Track branch",&trk);
00665
00666 timer.Start();
00667 int n = (int) t1->GetEntries();
00668 double etot=0;
00669 for (int i = 0; i < n; ++i) {
00670 t1->GetEntry(i);
00671 etot += SumSMatrix(trk->CovMatrix());
00672 }
00673
00674 timer.Stop();
00675
00676 std::cout << " Time to Read TrackD of SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00677
00678 #ifdef DEBUG
00679 std::cout << " Tree Entries " << n << std::endl;
00680 int pr = std::cout.precision(18);
00681 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00682 std::cout.precision(pr);
00683 #endif
00684
00685 return etot/double(n);
00686 }
00687
00688 double readTrackD32() {
00689
00690 std::cout << "\n";
00691 std::cout << "**************************************************\n";
00692 std::cout << " Test reading Track32 class........\n";
00693 std::cout << "**************************************************\n";
00694
00695 TFile f1("track32.root");
00696 if (f1.IsZombie() ) {
00697 std::cerr << "Error opening the ROOT file" << std::endl;
00698 return -1;
00699 }
00700
00701
00702 TTree *t1 = (TTree*)f1.Get("t1");
00703
00704 TrackD32 *trk = 0;
00705 t1->SetBranchAddress("Track32 branch",&trk);
00706
00707 timer.Start();
00708 int n = (int) t1->GetEntries();
00709 double etot=0;
00710 for (int i = 0; i < n; ++i) {
00711 t1->GetEntry(i);
00712 etot += SumSMatrix(trk->CovMatrix());
00713 }
00714
00715 timer.Stop();
00716
00717 std::cout << " Time to Read TrackD32 of SMatrix " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
00718
00719 #ifdef DEBUG
00720 std::cout << " Tree Entries " << n << std::endl;
00721 int pr = std::cout.precision(18);
00722 std::cout << " sum " << n<< " " << etot << " " << etot/double(n) << std::endl;
00723 std::cout.precision(pr);
00724 #endif
00725
00726 return etot/double(n);
00727 }
00728
00729
00730
00731
00732
00733 int testWrite(int nEvents, double & w1, double & w2) {
00734
00735 int iret = 0;
00736 double w0 = writeCArray(nEvents);
00737
00738 w1 = writeTMatrix(nEvents);
00739 w2 = writeSMatrix(nEvents,sfile1);
00740 if ( fabs(w1-w2) > tol) {
00741 std::cout << "\nERROR: Differeces SMatrix-TMatrix found when writing" << std::endl;
00742 int pr = std::cout.precision(18); std::cout << w1 << " != " << w2 << std::endl; std::cout.precision(pr);
00743 iret = 1;
00744 }
00745 if ( fabs(w1-w0) > tol) {
00746 std::cout << "\nERROR: Differeces TMatrix-C Array found when writing" << std::endl;
00747 int pr = std::cout.precision(18); std::cout << w1 << " != " << w0 << std::endl; std::cout.precision(pr);
00748 iret = 1;
00749 }
00750
00751 std::cout << "\n\n*************************************************************\n";
00752 if (iret == 0 )
00753 std::cout << " Writing Test:\t" << "OK";
00754 else {
00755 std::cout << " Writing Test:\t" << "FAILED";
00756 }
00757 std::cout << "\n*************************************************************\n\n";
00758
00759
00760 return iret;
00761 }
00762
00763 int testRead(double & r1, double & r2, double & r3) {
00764
00765 int iret = 0;
00766
00767
00768
00769 r1 = readTMatrix();
00770 r2 = readSMatrix(sfile1);
00771 if ( fabs(r1-r2) > tol) {
00772 std::cout << "\nERROR: Differeces SMatrix-TMatrix found when reading " << std::endl;
00773 int pr = std::cout.precision(18); std::cout << r1 << " != " << r2 << std::endl; std::cout.precision(pr);
00774 iret = 2;
00775 }
00776
00777 #ifdef USE_REFLEX
00778 std::cout << "try to read file written with CINT using Reflex Dictionaries " << std::endl;
00779 #else
00780 std::cout << "try to read file written with Reflex using CINT Dictionaries " << std::endl;
00781 #endif
00782 r3 = readSMatrix(sfile2);
00783 if ( r3 != -1. && fabs(r2-r3) > tol) {
00784 std::cout << "\nERROR: Differeces Reflex-CINT found when reading SMatrices" << std::endl;
00785 int pr = std::cout.precision(18); std::cout << r2 << " != " << r3 << std::endl; std::cout.precision(pr);
00786 iret = 3;
00787 }
00788
00789
00790 return iret;
00791 }
00792
00793
00794 int testWriteSym(int nEvents, double & w1, double & w2) {
00795
00796 int iret = 0;
00797
00798
00799
00800 w1 = writeTMatrixSym(nEvents);
00801 w2 = writeSMatrixSym(nEvents,symfile1);
00802 if ( fabs(w1-w2) > tol) {
00803 std::cout << "\nERROR: Differeces found when writing" << std::endl;
00804 int pr = std::cout.precision(18); std::cout << w1 << " != " << w2 << std::endl; std::cout.precision(pr);
00805 iret = 11;
00806 }
00807
00808 std::cout << "\n\n*************************************************************\n";
00809 if (iret == 0 )
00810 std::cout << " Writing Test:\t" << "OK";
00811 else {
00812 std::cout << " Writing Test:\t" << "FAILED";
00813 }
00814 std::cout << "\n*************************************************************\n\n";
00815
00816 return iret;
00817 }
00818
00819 int testReadSym(double & r1, double & r2, double & r3) {
00820
00821 int iret = 0;
00822
00823
00824 r1 = readTMatrixSym();
00825 r2 = readSMatrixSym(symfile1);
00826 if ( fabs(r1-r2) > tol) {
00827 std::cout << "\nERROR: Differeces SMatrixSym-TMAtrixSym found when reading " << std::endl;
00828 int pr = std::cout.precision(18); std::cout << r1 << " != " << r2 << std::endl; std::cout.precision(pr);
00829 iret = 12;
00830 }
00831
00832 #ifdef USE_REFLEX
00833 std::cout << "try to read file written with CINT using Reflex Dictionaries " << std::endl;
00834 #else
00835 std::cout << "try to read file written with Reflex using CINT Dictionaries " << std::endl;
00836 #endif
00837
00838 r3 = readSMatrixSym(symfile2);
00839 if ( r3 != -1. && fabs(r2-r3) > tol) {
00840 std::cout << "\nERROR: Differeces Reflex-CINT found when reading SMatricesSym" << std::endl;
00841 int pr = std::cout.precision(18); std::cout << r2 << " != " << r3 << std::endl; std::cout.precision(pr);
00842 iret = 13;
00843 }
00844
00845
00846 return iret;
00847 }
00848 int testResult(double w1, double r1, double w2, double r2, double r3) {
00849
00850 int iret = 0;
00851
00852 if ( fabs(w1-r1) > tol) {
00853 std::cout << "\nERROR: Differeces found when reading TMatrices" << std::endl;
00854 int pr = std::cout.precision(18); std::cout << w1 << " != " << r1 << std::endl; std::cout.precision(pr);
00855 iret = -1;
00856 }
00857 if ( fabs(w2-r2) > tol) {
00858 std::cout << "\nERROR: Differeces found when reading SMatrices" << std::endl;
00859 int pr = std::cout.precision(18); std::cout << w2 << " != " << r2 << std::endl; std::cout.precision(pr);
00860 iret = -2;
00861 }
00862 if ( r3 != -1. && fabs(w2-r3) > tol) {
00863 std::cout << "\nERROR: Differeces found when reading SMatrices with different Dictionary" << std::endl;
00864 int pr = std::cout.precision(18); std::cout << w2 << " != " << r2 << std::endl; std::cout.precision(pr);
00865 iret = -3;
00866 }
00867 return iret;
00868 }
00869
00870 int testTrack(int nEvents) {
00871
00872 int iret = 0;
00873
00874 double wt1 = writeTrackD(nEvents);
00875
00876 #ifdef USE_REFLEX
00877
00878 gSystem->Load("libSmatrix");
00879 #endif
00880
00881 double wt2 = writeTrackD32(nEvents);
00882
00883 if ( fabs(wt2-wt1) > tol) {
00884 std::cout << "\nERROR: Differeces found when writing Track" << std::endl;
00885 int pr = std::cout.precision(18); std::cout << wt2 << " != " << wt1 << std::endl; std::cout.precision(pr);
00886 iret = 13;
00887 }
00888
00889 double rt1 = readTrackD();
00890 if ( fabs(rt1-wt1) > tol) {
00891 std::cout << "\nERROR: Differeces found when reading Track" << std::endl;
00892 int pr = std::cout.precision(18); std::cout << rt1 << " != " << wt1 << std::endl; std::cout.precision(pr);
00893 iret = 13;
00894 }
00895
00896 double rt2 = readTrackD32();
00897 if ( fabs(rt2-wt2) > tol32) {
00898 std::cout << "\nERROR: Differeces found when reading Track 32" << std::endl;
00899 int pr = std::cout.precision(18); std::cout << rt2 << " != " << wt2 << std::endl; std::cout.precision(pr);
00900 iret = 13;
00901 }
00902
00903 return iret;
00904 }
00905
00906
00907 int testIO() {
00908
00909
00910 int iret = 0;
00911
00912
00913 #ifdef USE_REFLEX
00914
00915
00916 gSystem->Load("libReflex");
00917 gSystem->Load("libCintex");
00918 ROOT::Cintex::Cintex::SetDebug(1);
00919 ROOT::Cintex::Cintex::Enable();
00920
00921 std::cout << "Use Reflex dictionary " << std::endl;
00922
00923 #ifdef USE_REFLEX_SMATRIX
00924 iret |= gSystem->Load("libSmatrixRflx");
00925 #endif
00926 iret |= gSystem->Load("libSmatrix");
00927
00928
00929 #else
00930
00931 iret |= gSystem->Load("libSmatrix");
00932
00933 #endif
00934
00935 iret |= gSystem->Load("libMatrix");
00936
00937
00938 int nEvents = 10000;
00939
00940 initMatrix(nEvents);
00941
00942 double w1, w2 = 0;
00943 iret |= testWrite(nEvents,w1,w2);
00944
00945
00946
00947 double r1, r2, r3 = 0;
00948 int iret2 = 0;
00949 iret2 |= testRead(r1,r2,r3);
00950 iret2 |= testResult(w1,r1,w2,r2,r3);
00951 std::cout << "\n\n*************************************************************\n";
00952 if (iret2 == 0 )
00953 std::cout << " Reading Test:\t" << "OK";
00954 else {
00955 std::cout << " Reading Test:\t" << "FAILED";
00956 }
00957 std::cout << "\n*************************************************************\n\n";
00958
00959 iret |= iret2;
00960
00961
00962 std::cout << "\n*****************************************************\n";
00963 std::cout << " Test Symmetric matrices";
00964 std::cout << "\n*****************************************************\n\n";
00965
00966 iret = testWriteSym(nEvents,w1,w2);
00967 iret2 = testReadSym(r1,r2,r3);
00968 iret2 = testResult(w1,r1,w2,r2,r3);
00969
00970 std::cout << "\n\n*************************************************************\n";
00971 if (iret2 == 0 )
00972 std::cout << " Reading Test:\t" << "OK";
00973 else {
00974 std::cout << " Reading Test:\t" << "FAILED";
00975 }
00976 std::cout << "\n*************************************************************\n\n";
00977
00978 iret |= iret2;
00979
00980
00981 std::cout << "\n*****************************************************\n";
00982 std::cout << " Test Track class";
00983 std::cout << "\n*****************************************************\n\n";
00984
00985
00986 iret |= gSystem->Load("libTrackDict");
00987 if (iret != 0 ) return iret;
00988
00989 iret |= testTrack(nEvents);
00990 std::cout << "\n\n*************************************************************\n";
00991 if (iret2 == 0 )
00992 std::cout << " Track Test:\t" << "OK";
00993 else {
00994 std::cout << " Track Test:\t" << "FAILED";
00995 }
00996 std::cout << "\n*************************************************************\n\n";
00997
00998
00999 return iret;
01000
01001 }
01002
01003
01004
01005 int main() {
01006 int iret = testIO();
01007 std::cout << "\n\n*************************************************************\n";
01008 if (iret != 0) {
01009 std::cerr << "\nERROR !!!!! " << iret << std::endl;
01010 std::cerr << "TESTIO \t FAILED " << std::endl;
01011 }
01012 else
01013 std::cerr << "TESTIO \t OK " << std::endl;
01014 return iret;
01015 std::cout << "*************************************************************\n\n";
01016 }