testOperations.cxx

Go to the documentation of this file.
00001 #define  ENABLE_TEMPORARIES
00002 
00003 #include "Math/SVector.h"
00004 #include "Math/SMatrix.h"
00005 
00006 
00007 #include "TMatrixD.h"
00008 #include "TVectorD.h"
00009 
00010 #include "TRandom3.h"
00011 #include "TH1D.h" 
00012 #include "TProfile.h" 
00013 #include "TFile.h" 
00014 
00015 //#define HAVE_CLHEP
00016 #define TEST_SYM
00017 
00018 #ifdef TEST_ALL_MATRIX_SIZES
00019 #define REPORT_TIME
00020 #endif
00021 #ifndef NITER
00022 #define NITER 1  // number of iterations
00023 #endif
00024 #ifndef NLOOP_MIN
00025 #define NLOOP_MIN 1000;
00026 #endif
00027 
00028 #ifdef HAVE_CLHEP
00029 #include "CLHEP/Matrix/SymMatrix.h"
00030 #include "CLHEP/Matrix/Matrix.h"
00031 #include "CLHEP/Matrix/Vector.h"
00032 #endif
00033 
00034 //#define DEBUG
00035 
00036 #include <iostream>
00037 
00038 // #ifndef NDIM1
00039 // #define NDIM1 5
00040 // #endif
00041 // #ifndef NDIM2
00042 // #define NDIM2 5
00043 // #endif
00044 
00045 
00046 int NLOOP; 
00047 //#define NLOOP 1
00048 
00049 //#define DEBUG
00050 
00051 
00052 #include "matrix_op.h"
00053 #include "matrix_util.h"
00054 #include <map>
00055 
00056 
00057 
00058 
00059 
00060 template<unsigned int NDIM1, unsigned int NDIM2> 
00061 int test_smatrix_op() {
00062     
00063   // need to write explicitly the dimensions
00064    
00065 
00066   typedef SMatrix<double, NDIM1, NDIM1> MnMatrixNN;
00067   typedef SMatrix<double, NDIM2, NDIM2> MnMatrixMM;
00068   typedef SMatrix<double, NDIM1, NDIM2> MnMatrixNM;
00069   typedef SMatrix<double, NDIM2 , NDIM1> MnMatrixMN;
00070   typedef SVector<double, NDIM1> MnVectorN;
00071   typedef SVector<double, NDIM2> MnVectorM;
00072   
00073 
00074 
00075   int first = NDIM1;  //Can change the size of the matrices
00076   int second = NDIM2;
00077   
00078 
00079   std::cout << "************************************************\n";
00080   std::cout << "  SMatrix operations test  "   <<  first << " x " << second  << std::endl;
00081   std::cout << "************************************************\n";
00082 
00083 
00084   double t_veq, t_meq, t_vad, t_mad, t_dot, t_mv, t_gmv, t_mm, t_prd, t_inv, t_vsc, t_msc, t_ama, t_tra = 0;
00085   double totTime1, totTime2; 
00086   
00087   
00088    
00089   double r1,r2;
00090   int npass = NITER; 
00091   TRandom3 r(111);
00092   for (int k = 0; k < npass; k++) {
00093 
00094 
00095     MnMatrixNM A;
00096     MnMatrixMN B;
00097     MnMatrixNN C; 
00098     MnMatrixMM D; 
00099     MnVectorN v;
00100     MnVectorM p;
00101 
00102 
00103     TStopwatch w; 
00104     {       
00105       //seal::SealTimer t(init,false);
00106       // fill matrices with random data
00107       fillRandomMat(r,A,first,second);
00108       fillRandomMat(r,B,second,first);
00109       fillRandomMat(r,C,first,first);
00110       fillRandomMat(r,D,second,second);
00111 
00112       fillRandomVec(r,v,first);
00113       fillRandomVec(r,p,second);
00114 
00115     }
00116 
00117 
00118 //     MnSymMatrixMM I; 
00119 //     for(int i = 0; i < second; i++) 
00120 //       I(i,i) = 1;
00121      
00122 #ifdef DEBUG
00123     std::cout << "pass " << k << std::endl;
00124     if (k == 0) { 
00125       std::cout << " A = " << A << std::endl;
00126       std::cout << " B = " << B << std::endl;
00127       std::cout << " C = " << C << std::endl;
00128       std::cout << " D = " << D << std::endl;
00129       std::cout << " v = " << v << std::endl;
00130       std::cout << " p = " << p << std::endl;
00131     }
00132 #endif
00133 
00134     w.Start(); 
00135 
00136 
00137 
00138                 
00139     MnVectorN   v1;  testMV(A,v,t_mv,v1);
00140     //if (k == 0) v1.Print(std::cout);
00141     MnVectorN   v2;  testGMV(A,v,v1,t_gmv,v2);
00142     //if (k == 0) v2.Print(std::cout);
00143     MnMatrixNN  C0;  testMM(A,B,C,t_mm,C0);
00144     //if (k == 0) C0.Print(std::cout);
00145     MnMatrixNN  C1;  testATBA_S(B,C0,t_ama,C1);
00146     //if (k == 0) C1.Print(std::cout);
00147     MnMatrixNN  C2;  testInv_S(C1,t_inv,C2);
00148     MnVectorN   v3;  testVeq(v,t_veq,v3);
00149     MnVectorN   v4;  testVad(v2,v3,t_vad,v4);
00150     MnVectorN   v5;  testVscale(v4,2.0,t_vsc,v5);
00151     MnMatrixNN  C3;  testMeq(C,t_meq,C3);
00152     MnMatrixNN  C4;  testMad(C2,C3,t_mad,C4);
00153     MnMatrixNN  C5;  testMscale(C4,0.5,t_msc,C5);
00154     MnMatrixNN  C6;  testMT_S(C5,t_tra,C6);
00155 
00156 #ifdef DEBUG
00157     if (k == 0) { 
00158       std::cout << " C6 = " << C5 << std::endl;
00159       std::cout << " v5 = " << v5 << std::endl;
00160     }
00161 #endif
00162 
00163     r1 = testDot_S(v3,v5,t_dot);
00164 
00165     r2 = testInnerProd_S(C6,v5,t_prd);
00166 
00167   
00168     w.Stop();
00169     totTime1 = w.RealTime();
00170     totTime2 = w.CpuTime();
00171   
00172   }
00173   //tr.dump();
00174 
00175   //double totTime = t_veq + t_meq + t_vad + t_mad + t_dot + t_mv + t_gmv + t_mm + t_prd + t_inv + t_vsc + t_msc + t_ama + t_tra; 
00176   std::cout << "Total Time = " << totTime1 << "  (s) " << " cpu " <<  totTime2 << "  (s) " << std::endl; 
00177   std::cerr << "SMatrix:     r1 = " << r1 << " r2 = " << r2 << std::endl; 
00178 
00179   return 0;
00180 }
00181 
00182 
00183 
00184 #ifdef TEST_SYM
00185 template<unsigned int NDIM1, unsigned int NDIM2> 
00186 int test_smatrix_sym_op() {
00187     
00188   // need to write explicitly the dimensions
00189    
00190 
00191   typedef SMatrix<double, NDIM1, NDIM1, MatRepSym<double,NDIM1> > MnSymMatrixNN;
00192   typedef SMatrix<double, NDIM1, NDIM1 > MnMatrixNN;
00193   typedef SVector<double, NDIM1> MnVectorN;
00194   
00195 
00196 
00197   int first = NDIM1;  //Can change the size of the matrices
00198   
00199 
00200   std::cout << "************************************************\n";
00201   std::cout << "  SMatrixSym operations test  "   <<  first << " x " << first << std::endl;
00202   std::cout << "************************************************\n";
00203 
00204 
00205   double t_meq, t_mad, t_mv, t_gmv, t_mm, t_prd, t_inv, t_msc, t_ama = 0;
00206   double totTime1, totTime2; 
00207   
00208   
00209    
00210   double r1;
00211   int npass = NITER; 
00212   TRandom3 r(111);
00213   for (int k = 0; k < npass; k++) {
00214 
00215 
00216     MnSymMatrixNN A;
00217     MnSymMatrixNN B;
00218     MnMatrixNN C;
00219     MnVectorN v;
00220 
00221 
00222     TStopwatch w; 
00223     {       
00224       // fill matrices with random data
00225       fillRandomSym(r,A,first);
00226       fillRandomSym(r,B,first);
00227       fillRandomMat(r,C,first,first);
00228 
00229       fillRandomVec(r,v,first);
00230 
00231     }
00232 
00233      
00234 #ifdef DEBUG
00235     std::cout << "pass " << k << std::endl;
00236     if (k == 0) { 
00237       std::cout << " A = " << A << std::endl;
00238       std::cout << " B = " << B << std::endl;
00239       std::cout << " C = " << C << std::endl;
00240       std::cout << " v = " << v << std::endl;
00241     }
00242 #endif
00243 
00244     w.Start(); 
00245                 
00246     MnVectorN   v1;  testMV(A,v,t_mv,v1);
00247     MnVectorN   v2;  testGMV(A,v,v1,t_gmv,v2);
00248     MnMatrixNN  C0;  testMM(A,B,C,t_mm,C0);
00249     MnSymMatrixNN  C1;  testATBA_S2(C0,B,t_ama,C1);
00250     MnSymMatrixNN  C2;  testInv_S(A,t_inv,C2);
00251     MnSymMatrixNN  C3;  testMeq(C2,t_meq,C3);
00252     MnSymMatrixNN  C4;  testMad(A,C3,t_mad,C4);
00253     MnSymMatrixNN  C5;  testMscale(C4,0.5,t_msc,C5);
00254 
00255     r1 = testInnerProd_S(C5,v2,t_prd);
00256 
00257      
00258 #ifdef DEBUG
00259     std::cout << "output matrices" << std::endl;
00260     if (k == 0) { 
00261       std::cout << " C1 = " << C1 << std::endl;
00262       std::cout << " C3 = " << C3 << std::endl;
00263       std::cout << " C4 = " << C4 << std::endl;
00264       std::cout << " C5 = " << C5 << std::endl;
00265     }
00266 #endif
00267 
00268   
00269     w.Stop();
00270     totTime1 = w.RealTime();
00271     totTime2 = w.CpuTime();
00272 
00273   
00274   }
00275   //tr.dump();
00276 
00277   //double totTime = t_meq + t_mv + t_gmv + t_mm + t_prd + t_inv + t_mad + t_msc + t_ama; 
00278   std::cout << "Total Time = " << totTime1 << "  (s)  -  cpu " <<  totTime2 << "  (s) " << std::endl; 
00279   std::cerr << "SMatrixSym:  r1 = " << r1 << std::endl; 
00280 
00281   return 0;
00282 }
00283 #endif
00284 
00285 
00286 // ROOT test 
00287 
00288 
00289 template<unsigned int NDIM1, unsigned int NDIM2> 
00290 int test_tmatrix_op() {
00291 
00292 
00293     
00294 
00295   typedef TMatrixD MnMatrix;
00296   typedef TVectorD MnVector;
00297   
00298 //   typedef boost::numeric::ublas::matrix<double>  MnMatrix;  
00299   //typedef HepSymMatrix MnSymMatrixHep; 
00300 
00301 
00302   int first = NDIM1;  //Can change the size of the matrices
00303   int second = NDIM2;
00304 
00305 
00306   std::cout << "************************************************\n";
00307   std::cout << "  TMatrix operations test  "   <<  first << " x " << second  << std::endl;
00308   std::cout << "************************************************\n";
00309   
00310   double t_veq, t_meq, t_vad, t_mad, t_dot, t_mv, t_gmv, t_mm, t_prd, t_inv, t_vsc, t_msc, t_ama, t_tra = 0;
00311   double totTime1, totTime2; 
00312    
00313   double r1,r2;
00314   int npass = NITER; 
00315   TRandom3 r(111);
00316   gMatrixCheck = 0;
00317   
00318   for (int k = 0; k < npass; k++) {
00319 
00320 
00321     MnMatrix   A(NDIM1,NDIM2);
00322     MnMatrix   B(NDIM2,NDIM1);
00323     MnMatrix   C(NDIM1,NDIM1); 
00324     MnMatrix   D(NDIM2,NDIM2); 
00325     MnVector   v(NDIM1);
00326     MnVector   p(NDIM2);
00327 
00328 
00329     TStopwatch w; 
00330     { 
00331       // fill matrices with random data
00332       fillRandomMat(r,A,first,second);
00333       fillRandomMat(r,B,second,first);
00334       fillRandomMat(r,C,first,first);
00335       fillRandomMat(r,D,second,second);
00336 
00337       fillRandomVec(r,v,first);
00338       fillRandomVec(r,p,second);
00339 
00340     }
00341      
00342 #ifdef DEBUG
00343     std::cout << "pass " << k << std::endl;
00344     if (k == 0) { 
00345       A.Print(); B.Print(); C.Print(); D.Print(); v.Print(); p.Print();
00346     }
00347 #endif
00348     w.Start(); 
00349         
00350     
00351     MnVector v1(NDIM1);        testMV_T(A,v,t_mv,v1);
00352     //if (k == 0) v1.Print();
00353     MnVector v2(NDIM1);        testGMV_T(A,v,v1,t_gmv,v2);
00354     //if (k == 0) v2.Print();
00355     MnMatrix C0(NDIM1,NDIM1);  testMM_T(A,B,C,t_mm,C0);
00356     //if (k == 0) C0.Print();
00357     MnMatrix C1(NDIM1,NDIM1);  testATBA_T(B,C0,t_ama,C1);
00358     //if (k == 0) C1.Print();
00359     MnMatrix C2(NDIM1,NDIM1);  testInv_T(C1,t_inv,C2);
00360     //if (k == 0) C2.Print();
00361     MnVector v3(NDIM1);        testVeq(v,t_veq,v3);
00362     MnVector v4(NDIM1);        testVad_T(v2,v3,t_vad,v4);
00363     MnVector v5(NDIM1);        testVscale_T(v4,2.0,t_vsc,v5);
00364     MnMatrix C3(NDIM1,NDIM1);  testMeq(C,t_meq,C3);
00365     MnMatrix C4(NDIM1,NDIM1);  testMad_T(C2,C3,t_mad,C4);
00366     //if (k == 0) C4.Print();
00367     MnMatrix C5(NDIM1,NDIM1);  testMscale_T(C4,0.5,t_msc,C5);
00368     //if (k == 0) C5.Print();
00369     MnMatrix C6(NDIM1,NDIM1);  testMT_T(C5,t_tra,C6);
00370 
00371 #ifdef DEBUG
00372     if (k == 0) { 
00373       C6.Print();
00374       v5.Print();
00375     }
00376 #endif
00377 
00378     r1 = testDot_T(v3,v5,t_dot);
00379 
00380     r2 = testInnerProd_T(C6,v5,t_prd);
00381 
00382     //MnMatrix C2b(NDIM1,NDIM1); testInv_T2(C1,t_inv2,C2b);
00383 
00384   
00385     w.Stop();
00386     totTime1 = w.RealTime();
00387     totTime2 = w.CpuTime();
00388   }
00389   //  tr.dump();
00390 
00391   //double totTime = t_veq + t_meq + t_vad + t_mad + t_dot + t_mv + t_gmv + t_mm + t_prd + t_inv + t_inv2 + t_vsc + t_msc + t_ama + t_tra; 
00392   std::cout << "Total Time = " << totTime1 << "  (s)  -  cpu " <<  totTime2 << "  (s) " << std::endl; 
00393   std::cerr << "TMatrix:     r1 = " << r1 << " r2 = " << r2 << std::endl; 
00394 
00395   return 0;
00396 
00397 }
00398 
00399 
00400 
00401 #ifdef TEST_SYM
00402 template<unsigned int NDIM1, unsigned int NDIM2> 
00403 int test_tmatrix_sym_op() {
00404     
00405   // need to write explicitly the dimensions
00406    
00407 
00408   typedef TMatrixDSym MnSymMatrix;
00409   typedef TMatrixD    MnMatrix;
00410   typedef TVectorD MnVector;
00411   
00412 
00413 
00414   int first = NDIM1;  //Can change the size of the matrices
00415   
00416 
00417   std::cout << "************************************************\n";
00418   std::cout << "  TMatrixSym operations test  "   <<  first << " x " << first << std::endl;
00419   std::cout << "************************************************\n";
00420 
00421 
00422   double t_meq, t_mad, t_mv, t_gmv, t_mm, t_prd, t_inv, t_msc, t_ama = 0;
00423   double totTime1, totTime2; 
00424   
00425   
00426    
00427   double r1;
00428   int npass = NITER; 
00429   TRandom3 r(111);
00430   for (int k = 0; k < npass; k++) {
00431 
00432 
00433     MnSymMatrix A(NDIM1);
00434     MnSymMatrix B(NDIM1);
00435     MnMatrix C(NDIM1,NDIM1);
00436     MnVector v(NDIM1);
00437 #define N NDIM1
00438 
00439     TStopwatch w; 
00440 
00441     {       
00442       // fill matrices with random data
00443       fillRandomSym(r,A,first);
00444       fillRandomSym(r,B,first);
00445       fillRandomMat(r,C,first,first);
00446 
00447       fillRandomVec(r,v,first);
00448 
00449     }
00450 
00451      
00452 #ifdef DEBUG
00453     std::cout << "pass " << k << std::endl;
00454     if (k == 0) { 
00455       A.Print(); B.Print(); C.Print();  v.Print(); 
00456     }
00457 #endif
00458 
00459     w.Start(); 
00460                 
00461     MnVector   v1(N);  testMV_T(A,v,t_mv,v1);
00462     MnVector   v2(N);  testGMV_T(A,v,v1,t_gmv,v2);
00463     MnMatrix   C0(N,N);  testMM_T(A,B,C,t_mm,C0);
00464     MnSymMatrix  C1(N);  testATBA_T2(C0,B,t_ama,C1);
00465     MnSymMatrix  C2(N);  testInv_T(A,t_inv,C2);
00466     MnSymMatrix  C3(N);  testMeq(C2,t_meq,C3);
00467     MnSymMatrix  C4(N);  testMad_T(A,C3,t_mad,C4);
00468     MnSymMatrix  C5(N);  testMscale_T(C4,0.5,t_msc,C5);
00469 
00470     r1 = testInnerProd_T(C5,v2,t_prd);
00471 
00472 #ifdef DEBUG
00473     std::cout << "output matrices" << std::endl;
00474     if (k == 0) { 
00475       C1.Print(); C3.Print(); C4.Print(); C5.Print();
00476     }
00477 #endif
00478 
00479     w.Stop();
00480     totTime1 = w.RealTime();
00481     totTime2 = w.CpuTime();
00482   
00483   }
00484   //tr.dump();
00485 
00486   //double totTime = t_meq + t_mv + t_gmv + t_mm + t_prd + t_inv + t_mad + t_msc + t_ama; 
00487   std::cout << "Total Time = " << totTime1 << "  (s)  -  cpu " <<  totTime2 << "  (s) " << std::endl; 
00488   std::cerr << "TMatrixSym:  r1 = " << r1 << std::endl; 
00489 
00490   return 0;
00491 }
00492 #endif  // end TEST_SYM
00493 
00494 #ifdef HAVE_CLHEP
00495 
00496 template<unsigned int NDIM1, unsigned int NDIM2> 
00497 int test_hepmatrix_op() {
00498 
00499 
00500     
00501 
00502   typedef HepMatrix MnMatrix;
00503   typedef HepVector MnVector;
00504   
00505 
00506 
00507   int first = NDIM1;  //Can change the size of the matrices
00508   int second = NDIM2;
00509 
00510 
00511   std::cout << "************************************************\n";
00512   std::cout << "  HepMatrix operations test  "   <<  first << " x " << second  << std::endl;
00513   std::cout << "************************************************\n";
00514   
00515   double t_veq, t_meq, t_vad, t_mad, t_dot, t_mv, t_gmv, t_mm, t_prd, t_inv, t_vsc, t_msc, t_ama, t_tra = 0;
00516   
00517 
00518   double totTime1, totTime2; 
00519    
00520   double r1,r2;
00521   int npass = NITER; 
00522   TRandom3 r(111);
00523 
00524   for (int k = 0; k < npass; k++) {
00525 
00526 
00527     MnMatrix   A(NDIM1,NDIM2);
00528     MnMatrix   B(NDIM2,NDIM1);
00529     MnMatrix   C(NDIM1,NDIM1); 
00530     MnMatrix   D(NDIM2,NDIM2); 
00531     MnVector   v(NDIM1);
00532     MnVector   p(NDIM2);
00533 
00534     TStopwatch w; 
00535 
00536     { 
00537       // fill matrices with random data
00538       fillRandomMat(r,A,first,second,1);
00539       fillRandomMat(r,B,second,first,1);
00540       fillRandomMat(r,C,first,first,1);
00541       fillRandomMat(r,D,second,second,1);
00542 
00543       fillRandomVec(r,v,first);
00544       fillRandomVec(r,p,second);
00545     }
00546 
00547 #ifdef DEBUG
00548     std::cout << "pass " << k << std::endl;
00549     if (k == 0) { 
00550       std::cout << " A = " << A << std::endl;
00551       std::cout << " B = " << B << std::endl;
00552       std::cout << " C = " << C << std::endl;
00553       std::cout << " D = " << D << std::endl;
00554       std::cout << " v = " << v << std::endl;
00555       std::cout << " p = " << p << std::endl;
00556     }
00557 #endif
00558         
00559     w.Start(); 
00560     
00561     MnVector v1(NDIM1);        testMV(A,v,t_mv,v1);
00562     MnVector v2(NDIM1);        testGMV(A,v,v1,t_gmv,v2);
00563     MnMatrix C0(NDIM1,NDIM1);  testMM_C(A,B,C,t_mm,C0);
00564     MnMatrix C1(NDIM1,NDIM1);  testATBA_C(B,C0,t_ama,C1);
00565     //std::cout << " C1 = " << C1 << std::endl;
00566     MnMatrix C2(NDIM1,NDIM1);  testInv_C(C1,t_inv,C2);
00567     //std::cout << " C2 = " << C2 << std::endl;
00568     MnVector v3(NDIM1);        testVeq(v,t_veq,v3);
00569     MnVector v4(NDIM1);        testVad(v2,v3,t_vad,v4);
00570     MnVector v5(NDIM1);        testVscale(v4,2.0,t_vsc,v5);
00571     MnMatrix C3(NDIM1,NDIM1);  testMeq_C(C,t_meq,C3);
00572     MnMatrix C4(NDIM1,NDIM1);  testMad_C(C2,C3,t_mad,C4);
00573     //std::cout << " C4 = " << C4 << std::endl;
00574     MnMatrix C5(NDIM1,NDIM1);  testMscale_C(C4,0.5,t_msc,C5);
00575     //std::cout << " C5 = " << C5 << std::endl;
00576     MnMatrix C6(NDIM1,NDIM1);  testMT_C(C5,t_tra,C6);
00577 
00578 
00579     r1 = testDot_C(v3,v5,t_dot);
00580     r2 = testInnerProd_C(C6,v5,t_prd);
00581 
00582 #ifdef DEBUG
00583     if (k == 0) { 
00584       std::cout << " C6 = " << C6 << std::endl;
00585       std::cout << " v5 = " << v5 << std::endl;
00586     }
00587 #endif
00588 
00589     //    MnMatrix C2b(NDIM1,NDIM1); testInv_T2(C1,t_inv2,C2b);
00590 
00591     w.Stop();
00592     totTime1 = w.RealTime();
00593     totTime2 = w.CpuTime();
00594   
00595   }
00596   //  tr.dump();
00597 
00598   std::cout << "Total Time = " << totTime1 << "  (s)  -  cpu " <<  totTime2 << "  (s) " << std::endl; 
00599   std::cerr << "HepMatrix:   r1 = " << r1 << " r2 = " << r2 << std::endl; 
00600 
00601   return 0;
00602 }
00603 
00604 
00605 #ifdef TEST_SYM
00606 template<unsigned int NDIM1, unsigned int NDIM2> 
00607 int test_hepmatrix_sym_op() {
00608     
00609   // need to write explicitly the dimensions
00610    
00611 
00612   typedef HepSymMatrix MnSymMatrix;
00613   typedef HepMatrix    MnMatrix;
00614   typedef HepVector MnVector;
00615   
00616 
00617 
00618   int first = NDIM1;  //Can change the size of the matrices
00619   
00620 
00621   std::cout << "************************************************\n";
00622   std::cout << "  HepMatrixSym operations test  "   <<  first << " x " << first << std::endl;
00623   std::cout << "************************************************\n";
00624 
00625 
00626   double t_meq, t_mad, t_mv, t_gmv, t_mm, t_prd, t_inv, t_msc, t_ama = 0;
00627   
00628   double totTime1, totTime2; 
00629   
00630    
00631   double r1;
00632   int npass = NITER; 
00633   TRandom3 r(111);
00634   for (int k = 0; k < npass; k++) {
00635 
00636 
00637     MnSymMatrix A(NDIM1);
00638     MnSymMatrix B(NDIM1);
00639     MnMatrix C(NDIM1,NDIM1);
00640     MnVector v(NDIM1);
00641 #define N NDIM1
00642 
00643     TStopwatch w; 
00644 
00645     {       
00646       // fill matrices with random data
00647       fillRandomSym(r,A,first,1);
00648       fillRandomSym(r,B,first,1);
00649       fillRandomMat(r,C,first,first,1);
00650       fillRandomVec(r,v,first);
00651 
00652     }
00653 
00654      
00655 #ifdef DEBUG
00656     std::cout << "pass " << k << std::endl;
00657     if (k == 0) { 
00658     }
00659 #endif
00660     
00661     w.Start(); 
00662                 
00663     MnVector   v1(N);  testMV(A,v,t_mv,v1);
00664     MnVector   v2(N);  testGMV(A,v,v1,t_gmv,v2);
00665     MnMatrix   C0(N,N);  testMM_C(A,B,C,t_mm,C0);
00666     MnSymMatrix  C1(N);  testATBA_C2(C0,B,t_ama,C1);
00667     MnSymMatrix  C2(N);  testInv_C(A,t_inv,C2);
00668     MnSymMatrix  C3(N);  testMeq_C(C2,t_meq,C3);
00669     MnSymMatrix  C4(N);  testMad_C(A,C3,t_mad,C4);
00670     MnSymMatrix  C5(N);  testMscale_C(C4,0.5,t_msc,C5);
00671 
00672     r1 = testInnerProd_C(C5,v2,t_prd);
00673 
00674 #ifdef DEBUG
00675     std::cout << "output matrices" << std::endl;
00676     if (k == 0) { 
00677     }
00678 #endif
00679 
00680     w.Stop();
00681     totTime1 = w.RealTime();
00682     totTime2 = w.CpuTime();
00683   
00684   }
00685   //tr.dump();
00686 
00687   std::cout << "Total Time = " << totTime1 << "  (s)  -  cpu " <<  totTime2 << "  (s) " << std::endl; 
00688   std::cerr << "HepMatrixSym: r1 = " << r1 << std::endl; 
00689 
00690   return 0;
00691 }
00692 
00693 #endif  // TEST_SYM
00694 #endif  // HAVE_CLHEP
00695 
00696 
00697 #if defined(HAVE_CLHEP) && defined (TEST_SYM)
00698 #define NTYPES 6
00699 #define TEST(N) \
00700    MATRIX_SIZE=N;  \
00701    TEST_TYPE=0; test_smatrix_op<N,N>(); \
00702    TEST_TYPE=1; test_tmatrix_op<N,N>();     \
00703    TEST_TYPE=2; test_hepmatrix_op<N,N>();   \
00704    TEST_TYPE=3; test_smatrix_sym_op<N,N>(); \
00705    TEST_TYPE=4; test_tmatrix_sym_op<N,N>();     \
00706    TEST_TYPE=5; test_hepmatrix_sym_op<N,N>();
00707 #elif !defined(HAVE_CLHEP) && defined (TEST_SYM)
00708 #define NTYPES 4
00709 #define TEST(N) \
00710    MATRIX_SIZE=N;  \
00711    TEST_TYPE=0; test_smatrix_op<N,N>(); \
00712    TEST_TYPE=1; test_tmatrix_op<N,N>();     \
00713    TEST_TYPE=2; test_smatrix_sym_op<N,N>(); \
00714    TEST_TYPE=3; test_tmatrix_sym_op<N,N>();     
00715 #elif defined(HAVE_CLHEP) && !defined (TEST_SYM)
00716 #define NTYPES 3
00717 #define TEST(N) \
00718    MATRIX_SIZE=N;  \
00719    TEST_TYPE=0; test_smatrix_op<N,N>(); \
00720    TEST_TYPE=1; test_tmatrix_op<N,N>();     \
00721    TEST_TYPE=2; test_hepmatrix_op<N,N>();
00722 #else 
00723 #define NTYPES 2
00724 #define TEST(N) \
00725    TEST_TYPE=0; test_smatrix_op<N,N>(); \
00726    TEST_TYPE=1; test_tmatrix_op<N,N>();     
00727 #endif
00728 
00729 
00730 
00731 int TEST_TYPE; 
00732 int MATRIX_SIZE; 
00733 #ifdef REPORT_TIME
00734 std::vector< std::map<std::string, TH1D *> > testTimeResults(NTYPES); 
00735 std::vector< std::string > typeNames(NTYPES);
00736 
00737 void ROOT::Math::test::reportTime(std::string s, double time) { 
00738   assert( TEST_TYPE >= 0 && TEST_TYPE < NTYPES );
00739   std::map<std::string, TH1D * > & result = testTimeResults[TEST_TYPE];
00740   
00741   std::map<std::string, TH1D * >::iterator pos = result.find(s);   
00742   TH1D * h = 0; 
00743   if (  pos != result.end() ) { 
00744     h = pos->second; 
00745   }
00746   else { 
00747     // add new elements in map
00748     //std::cerr << "insert element in map" << s << typeNames[TEST_TYPE] << std::endl;
00749     std::string name = typeNames[TEST_TYPE] + "_" + s; 
00750     h = new TProfile(name.c_str(), name.c_str(),100,0.5,100.5);
00751     //result.insert(std::map<std::string, TH1D * >::value_type(s,h) ); 
00752     result[s] = h;
00753   }
00754   double scale=1; 
00755   if (s.find("dot") != std::string::npos || 
00756       s.find("V=V") != std::string::npos || 
00757       s.find("V+V") != std::string::npos ) scale = 10;  
00758   h->Fill(double(MATRIX_SIZE),time/double(NLOOP*NITER*scale) ); 
00759 }
00760 #endif
00761 
00762 int testOperations() { 
00763 
00764   NLOOP = 1000*NLOOP_MIN 
00765   initValues();
00766 
00767   TEST(5)
00768 
00769    return 0;   
00770 }
00771 
00772 
00773 int main(int argc , char *argv[] ) { 
00774 
00775 
00776   std::string fname = "testOperations"; 
00777   if (argc > 1) { 
00778     std::string platf(argv[1]); 
00779     fname = fname + "_" + platf; 
00780   }
00781   fname = fname + ".root"; 
00782 
00783 
00784 #ifdef REPORT_TIME
00785   TFile * f = new TFile(fname.c_str(),"recreate");
00786 
00787   typeNames[0] = "SMatrix";
00788   typeNames[1] = "TMatrix";
00789 #if !defined(HAVE_CLHEP) && defined (TEST_SYM)
00790   typeNames[2] = "SMatrix_sym";
00791   typeNames[3] = "TMatrix_sym";
00792 #elif defined(HAVE_CLHEP) && defined (TEST_SYM)
00793   typeNames[2] = "HepMatrix";
00794   typeNames[3] = "SMatrix_sym";
00795   typeNames[4] = "TMatrix_sym";
00796   typeNames[5] = "HepMatrix_sym";
00797 #elif defined(HAVE_CLHEP) && !defined (TEST_SYM)
00798   typeNames[2] = "HepMatrix";
00799 #endif
00800 
00801 #endif
00802 
00803 #ifndef TEST_ALL_MATRIX_SIZES
00804 //   NLOOP = 1000*NLOOP_MIN 
00805 //   initValues();
00806 
00807 //   TEST(5)
00808 //   NLOOP = 50*NLOOP_MIN;
00809 //   TEST(30);
00810 
00811   return testOperations();
00812 
00813 #else
00814   NLOOP = 5000*NLOOP_MIN; 
00815   initValues();
00816 
00817 
00818 
00819   TEST(2);
00820   TEST(3);
00821   TEST(4);
00822   NLOOP = 1000*NLOOP_MIN 
00823   TEST(5);
00824   TEST(6);
00825   TEST(7);
00826   TEST(10); 
00827   NLOOP = 100*NLOOP_MIN; 
00828   TEST(15); 
00829   TEST(20);
00830   NLOOP = 50*NLOOP_MIN; 
00831   TEST(30);
00832 //   NLOOP = NLOOP_MIN;
00833 //   TEST(50);
00834 //   TEST(75);
00835 //   TEST(100);
00836 #endif
00837 
00838 #ifdef REPORT_TIME
00839   f->Write();
00840   f->Close();
00841 #endif
00842 
00843 }
00844 

Generated on Tue Jul 5 14:43:12 2011 for ROOT_528-00b_version by  doxygen 1.5.1