stressVector.cxx

Go to the documentation of this file.
00001 // @(#)root/test:$Id: stressVector.cxx 20882 2007-11-19 11:31:26Z rdm $
00002 // Author: Lorenzo Moneta   06/2005 
00003 ///////////////////////////////////////////////////////////////////////////////////
00004 //
00005 //  MathCore Benchmark test suite
00006 //  ==============================
00007 //
00008 //  This program performs tests of ROOT::Math 4D LorentzVectors comparing with TLorentzVector
00009 //  The time performing various vector operations on a collection of vectors is measured. 
00010 //  The benchmarked operations are: 
00011 //      - vector construction from 4 values
00012 //      - construction using a setter method
00013 //      - simple addition of all the vector pairs in the collection 
00014 //      - calculation of deltaR = phi**2 + eta**2 of all vector pairs in the collection            
00015 //      -  two simple analysis: 
00016 //         - the first requires some cut (on pt and eta) and on the invariant mass 
00017 //           of the selected pairs
00018 //         - the second requires just some cut in pt, eta and delta R on all the 
00019 //           vector pair
00020 //      - conversion between XYZTVectors to PtRhoEtaPhi based vectors
00021 //
00022 //   The two analysis demostrates, especially in the second case, the advantage of using 
00023 //   vector based on cylindrical coordinate, given the fact that the time spent in the conversion is 
00024 //   much less than the time spent in the analysis routine. 
00025 //
00026 //   To run the program do: 
00027 //   stressVector          : run standard test with collection of 1000 vectors
00028 //   stressVector  10000   : run with a collection of 10000 vectors
00029 //
00030 ///////////////////////////////////////////////////////////////////////////////////
00031 
00032 
00033 
00034 #include <vector>
00035 #include <iostream>
00036 #include <algorithm>
00037 
00038 #include <assert.h>
00039 #include <map>
00040 
00041 #include "TStopwatch.h"
00042 
00043 #include "TRandom3.h"
00044 #include "TLorentzVector.h"
00045 #include "TRotation.h"
00046 #include "TLorentzRotation.h"
00047 
00048 #include "TMatrixD.h"
00049 
00050 #include "Math/Vector3D.h"
00051 #include "Math/Vector4D.h"
00052 #include "Math/VectorUtil.h"
00053 #include "Math/LorentzRotation.h"
00054 #include "Math/Rotation3D.h"
00055 #include "Math/RotationX.h"
00056 #include "Math/RotationY.h"
00057 #include "Math/RotationZ.h"
00058 
00059 #include "Math/SMatrix.h"
00060 
00061 
00062 #include "limits"
00063 
00064 //#define DEBUG
00065 
00066 using namespace ROOT::Math;
00067 
00068 
00069 
00070 
00071 class VectorTest { 
00072 
00073 private: 
00074 
00075 // global data variables 
00076   std::vector<double> dataX; 
00077   std::vector<double> dataY;  
00078   std::vector<double> dataZ;  
00079   std::vector<double> dataE; 
00080 
00081   size_t n2Loop ;
00082   size_t nGen;
00083 
00084 
00085 public: 
00086   
00087   VectorTest(int n1, int n2) : 
00088     n2Loop(n1),
00089     nGen(n2)
00090   {}
00091     
00092 
00093   
00094 
00095   void print(TStopwatch & time, std::string s) { 
00096     int pr = std::cout.precision(8);
00097     std::cout << s << "\t" << " time = " << time.RealTime() << "\t(sec)\t" 
00098       //    << time.CpuTime() 
00099               << std::endl;
00100     std::cout.precision(pr);
00101   }
00102 
00103     
00104   int check(std::string name, double s1, double s2, double s3, double scale=1) {
00105     double eps = 10*scale*std::numeric_limits<double>::epsilon();
00106     if (  std::fabs(s1-s2) < eps*std::fabs(s1) && std::fabs(s1-s3)  < eps*std::fabs(s1) ) return 0; 
00107     std::cout.precision(16);
00108     std::cout << s1 << "\t" << s2 <<"\t" << s3 << "\n";
00109     std::cout << "Test " << name << " failed !!\n\n"; 
00110     return -1; 
00111   }
00112 
00113 
00114   void genData() { 
00115     int n = nGen;
00116 
00117   // generate n -4 momentum quantities 
00118     TRandom3 r;
00119     for (int i = 0; i < n ; ++ i) { 
00120  
00121       double phi = r.Rndm()*3.1415926535897931; 
00122       double eta = r.Uniform(-5.,5.); 
00123       double pt = r.Exp(10.);
00124       double m = r.Uniform(0,10.); 
00125       if ( i%50 == 0 ) 
00126         m = r.BreitWigner(1.,0.01); 
00127 
00128       double E = sqrt( m*m + pt*pt*cosh(eta)*cosh(eta) );
00129     
00130     // fill vectors 
00131     
00132       PtEtaPhiEVector q( pt, eta, phi, E); 
00133       dataX.push_back( q.x() ); 
00134       dataY.push_back( q.y() ); 
00135       dataZ.push_back( q.z() ); 
00136       dataE.push_back( q.t() ); 
00137 
00138     }
00139   }
00140 
00141 
00142 
00143   template <class V> 
00144   void testCreate( std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00145     
00146     int n = dataX.size(); 
00147     dataV.resize(n); 
00148     tim.Start();
00149     for (int i = 0; i < n; ++i) { 
00150       dataV[i] = new V( dataX[i], dataY[i], dataZ[i], dataE[i] ); 
00151     }
00152     tim.Stop();
00153     t += tim.RealTime();
00154     print(tim,s);
00155   }
00156 
00157 
00158   template <class V> 
00159   void testCreate2( std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00160     
00161     int n = dataX.size(); 
00162     dataV.resize(n); 
00163     tim.Start();
00164     for (int i = 0; i < n; ++i) { 
00165       dataV[i] = new V(); 
00166       dataV[i]->SetXYZT(dataX[i], dataY[i], dataZ[i], dataE[i] );
00167     }
00168     tim.Stop();
00169     print(tim,s);
00170     t += tim.RealTime();
00171   }
00172 
00173 
00174 
00175   template <class V> 
00176   void clear(  std::vector<V *> & dataV ) { 
00177     for (unsigned int i = 0; i < dataV.size(); ++i) { 
00178       V * p = dataV[i]; 
00179       delete p; 
00180     }
00181     dataV.clear(); 
00182   
00183 }
00184 
00185 template <class V>
00186 double testAddition( const std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00187   unsigned int n = std::min(n2Loop, dataV.size() );
00188   double tot = 0;
00189   tim.Start(); 
00190   for (unsigned int i = 0; i < n; ++i) { 
00191     V  & v1 = *(dataV[i]); 
00192     for (unsigned int j = i +1; j < n; ++j) {
00193       V & v2 = *(dataV[j]); 
00194       V v3 = v1 + v2;
00195       tot += v3.E();
00196     }
00197   }
00198   tim.Stop();
00199   print(tim,s);
00200   t += tim.RealTime();
00201   return tot; 
00202 }  
00203 
00204 template <class V>
00205 double testScale( const std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00206   unsigned int n = std::min(n2Loop, dataV.size() );
00207   double tot = 0;
00208   tim.Start(); 
00209   for (unsigned int i = 0; i < n; ++i) { 
00210     V  & v1 = *(dataV[i]); 
00211     // scale
00212     v1 = 2.0*v1;
00213     tot += v1.E();
00214   }
00215   tim.Stop();
00216   print(tim,s);
00217   t += tim.RealTime();
00218   return tot; 
00219 }  
00220 
00221 
00222 template< class V> 
00223 bool cutPtEtaAndMass(const V & v) { 
00224  double pt = v.Pt();
00225  double mass = v.M();
00226  double eta = v.Eta();
00227  return ( pt > 3. && fabs(mass - 1.) < 0.5 && fabs(eta) < 3. );
00228 }
00229 
00230 
00231 template< class V> 
00232 bool cutPtEta(const V & v,double ptMin, double etaMax) { 
00233  double pt = v.Pt();
00234  double eta = v.Eta();
00235  return ( pt > ptMin && fabs(eta) < etaMax );
00236 }
00237 
00238 
00239 
00240 template <class V>
00241 double testDeltaR( const std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00242   unsigned int n = std::min(n2Loop, dataV.size() );
00243   tim.Start(); 
00244   double tot = 0;
00245   for (unsigned int i = 0; i < n; ++i) { 
00246     V  & v1 = *(dataV[i]); 
00247     for (unsigned int j = i +1; j < n; ++j) {
00248       V & v2 = *(dataV[j]); 
00249       double delta = VectorUtil::DeltaR(v1,v2);
00250       tot += delta;
00251     }
00252   }
00253   tim.Stop();
00254   print(tim,s);
00255   t += tim.RealTime();
00256   return tot;
00257 }  
00258 
00259 
00260 template <class V>
00261 int testAnalysis( const std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00262   int nsel = 0;  
00263   int nsel2 = 0; 
00264   double deltaMax = 1.;
00265   double ptMin = 1.;
00266   double etaMax = 3.;
00267   
00268   unsigned int n = std::min(n2Loop, dataV.size() );
00269   tim.Start(); 
00270   for (unsigned int i = 0; i < n; ++i) { 
00271     V  & v1 = *(dataV[i]); 
00272     if (cutPtEta(v1,ptMin, etaMax) ) { 
00273       double delta; 
00274       for (unsigned int j = i +1; j < n; ++j) {
00275         V & v2 = *(dataV[j]); 
00276         delta = VectorUtil::DeltaR(v1,v2);
00277         if (delta < deltaMax) { 
00278           V v3 = v1 + v2; 
00279           nsel++;
00280           if ( cutPtEtaAndMass(v3)) 
00281             nsel2++; 
00282         }
00283         
00284       }
00285     }
00286   }
00287   tim.Stop();
00288   print(tim,s);
00289   //std::cout << nsel << "\n"; 
00290   t += tim.RealTime();
00291   return nsel2; 
00292 }  
00293 
00294 
00295 
00296 template <class V>
00297 int testAnalysis2( const std::vector<V *> & dataV, TStopwatch & tim, double& t,  std::string s) { 
00298   int nsel = 0; 
00299   double ptMin = 1.;
00300   double etaMax = 3.;
00301   unsigned int n = std::min(n2Loop, dataV.size() );
00302   tim.Start();
00303   //seal::SealTimer t(tim.name(), true, std::cout); 
00304   for (unsigned int i = 0; i < n; ++i) { 
00305     V  & v1 = *(dataV[i]); 
00306     if ( cutPtEta(v1, ptMin, etaMax) ) { 
00307       for (unsigned int j = i +1; j < n; ++j) {
00308         V & v2 = *(dataV[j]); 
00309         if ( VectorUtil::DeltaR(v1,v2) < 0.5) nsel++;
00310       }
00311     }
00312   }
00313   tim.Stop();
00314   print(tim,s);
00315   t += tim.RealTime();
00316   return nsel; 
00317 }  
00318 
00319 
00320 
00321   template <class V1, class V2> 
00322   void testConversion( std::vector<V1 *> & dataV1, std::vector<V2 *> & dataV2, TStopwatch & tim, double& t,  std::string s) { 
00323     
00324     int n = dataX.size(); 
00325     dataV2.resize(n); 
00326     tim.Start();
00327     for (int i = 0; i < n; ++i) { 
00328       dataV2[i] = new V2( *dataV1[i] ); 
00329     }
00330     tim.Stop();
00331     print(tim,s);
00332     t += tim.RealTime();
00333   }
00334 
00335 
00336   // rotation using rotation classes 
00337   template <class V, class R> 
00338   double testRotation( std::vector<V *> & dataV, const R & rot, TStopwatch & tim, double& t,  std::string s) { 
00339     
00340     unsigned int n = std::min(n2Loop, dataV.size() );
00341     tim.Start();
00342     double sum = 0;
00343     for (unsigned int i = 0; i < n; ++i) { 
00344       V  & v1 = *(dataV[i]);
00345       V v2 = rot * v1; 
00346       sum += v2.X() + v2.Y() + v2.Z();
00347     }
00348     tim.Stop();
00349     print(tim,s);
00350     t += tim.RealTime();
00351     return sum;
00352   }
00353 
00354 
00355   // test matrix vector multiplication
00356   template <class V, class M> 
00357   double testMatVec( std::vector<V *> & dataV, const M & mat, TStopwatch & tim, double& t,  std::string s) { 
00358     
00359     unsigned int n = std::min(n2Loop, dataV.size() );
00360     tim.Start();
00361     double sum = 0;
00362     for (unsigned int i = 0; i < n; ++i) { 
00363       V  & v1 = *(dataV[i]);
00364       V v2 = VectorUtil::Mult(mat, v1 ); 
00365       sum += v2.X() + v2.Y() + v2.Z();
00366     }
00367     tim.Stop();
00368     print(tim,s);
00369     t += tim.RealTime();
00370     return sum;
00371   }
00372 
00373 
00374   // Boost using boost  classes 
00375   template <class V, class B> 
00376   double testBoost1( std::vector<V *> & dataV, const B & bv, TStopwatch & tim, double& t,  std::string s) { 
00377     
00378     unsigned int n = std::min(n2Loop, dataV.size() );
00379     tim.Start();
00380     double sum = 0;
00381     for (unsigned int i = 0; i < n; ++i) { 
00382       V  & v1 = *(dataV[i]);
00383       Boost b(bv); 
00384       V v2 = b(v1); 
00385       sum += v2.X() + v2.Y() + v2.Z() + v2.T();
00386     }
00387     tim.Stop();
00388     print(tim,s);
00389     t += tim.RealTime();
00390     return sum;
00391   }
00392 
00393 
00394   // Boost using vector util function
00395   template <class V, class B> 
00396   double testBoost2( std::vector<V *> & dataV, const B & bv, TStopwatch & tim, double& t,  std::string s) { 
00397     
00398     unsigned int n = std::min(n2Loop, dataV.size() );
00399     tim.Start();
00400     double sum = 0;
00401     for (unsigned int i = 0; i < n; ++i) { 
00402       V  & v1 = *(dataV[i]);
00403       V v2 = VectorUtil::boost(v1,bv); 
00404       sum += v2.X() + v2.Y() + v2.Z() + v2.T();
00405     }
00406     tim.Stop();
00407     print(tim,s);
00408     t += tim.RealTime();
00409     return sum;
00410   }
00411 
00412   // Boost using TLorentzVector
00413   double testBoost_TL( std::vector<TLorentzVector *> & dataV, const TVector3 & bv, TStopwatch & tim, double& t,  std::string s) { 
00414     
00415     unsigned int n = std::min(n2Loop, dataV.size() );
00416     tim.Start();
00417     double sum = 0;
00418     for (unsigned int i = 0; i < n; ++i) { 
00419       TLorentzVector  v2 = *(dataV[i]);
00420       v2.Boost(bv);  
00421       sum += v2.X() + v2.Y() + v2.Z() + v2.T();
00422     }
00423     tim.Stop();
00424     print(tim,s);
00425     t += tim.RealTime();
00426     return sum;
00427   }
00428 
00429 
00430 
00431   // Boost using boost  classes 
00432   template <class V> 
00433   double testBoostX1( std::vector<V *> & dataV, double beta, TStopwatch & tim, double& t,  std::string s) { 
00434     
00435     unsigned int n = std::min(n2Loop, dataV.size() );
00436     tim.Start();
00437     double sum = 0;
00438     for (unsigned int i = 0; i < n; ++i) { 
00439       V  & v1 = *(dataV[i]);
00440       BoostX b(beta); 
00441       V v2 = b(v1); 
00442       sum += v2.X() + v2.Y() + v2.Z() + v2.T();
00443     }
00444     tim.Stop();
00445     print(tim,s);
00446     t += tim.RealTime();
00447     return sum;
00448   }
00449 
00450   // Boost using vector util function
00451   template <class V> 
00452   double testBoostX2( std::vector<V *> & dataV, double beta, TStopwatch & tim, double& t,  std::string s) { 
00453     
00454     unsigned int n = std::min(n2Loop, dataV.size() );
00455     tim.Start();
00456     double sum = 0;
00457     for (unsigned int i = 0; i < n; ++i) { 
00458       V  & v1 = *(dataV[i]);
00459       V v2 = VectorUtil::boostX(v1,beta); 
00460       sum += v2.X() + v2.Y() + v2.Z() + v2.T();
00461     }
00462     tim.Stop();
00463     print(tim,s);
00464     t += tim.RealTime();
00465     return sum;
00466   }
00467 
00468 
00469 
00470 };
00471 
00472 int main(int argc,const char *argv[]) { 
00473 
00474   int ngen = 1000;
00475   if (argc > 1)  ngen = atoi(argv[1]);
00476   int nloop2 = ngen;
00477   if (argc > 2)  nloop2 = atoi(argv[1]);
00478 
00479 
00480   TStopwatch t;
00481 
00482   VectorTest a(ngen,nloop2);
00483 
00484   a.genData(); 
00485 
00486   int niter = 1;
00487   for (int i = 0; i < niter; ++i) { 
00488 
00489 #ifdef DEBUG
00490       std::cout << "iteration " << i << std::endl;
00491 #endif
00492     
00493       double t1 = 0;
00494       double t2 = 0;
00495       double t3 = 0;
00496 
00497       std::vector<TLorentzVector *> v1;
00498       std::vector<XYZTVector *> v2;
00499       std::vector<PtEtaPhiEVector *> v3;
00500 
00501       a.testCreate     (v1, t, t1,    "creation TLorentzVector      " ); 
00502       a.testCreate     (v2, t, t2,    "creation XYZTVector          " ); 
00503       a.testCreate     (v3, t, t3,     "creation PtEtaPhiEVector     " ); 
00504 
00505       
00506       a.clear(v3);
00507       std::cout << "\n";
00508       a.testConversion  (v2, v3, t, t3,   "Conversion XYZT->PtEtaPhiE " ); 
00509 
00510       a.clear(v1);
00511       a.clear(v2);
00512       a.clear(v3); 
00513 
00514       a.testCreate2     (v1, t, t1,   "creationSet TLorentzVector   " ); 
00515       a.testCreate2     (v2, t, t2,   "creationSet  XYZTVector      " ); 
00516       a.testCreate2     (v3, t, t3,   "creationSet  PtEtaPhiEVector " ); 
00517 
00518 
00519       double s1,s2,s3;
00520       s1=a.testAddition   (v1, t, t1, "Addition TLorentzVector      " );  
00521       s2=a.testAddition   (v2, t, t2, "Addition XYZTVector          "  ); 
00522       s3=a.testAddition   (v3, t, t3, "Addition PtEtaPhiEVector     " );       
00523       a.check("Addition",s1,s2,s3);
00524 
00525 
00526       s1=a.testScale   (v1, t, t1, "Scale of TLorentzVector      " );  
00527       s2=a.testScale   (v2, t, t2, "Scale of XYZTVector          "  ); 
00528       s3=a.testScale   (v3, t, t3, "Scale of PtEtaPhiEVector     " ); 
00529       a.check("Scaling",s1,s2,s3);
00530 
00531       s1=a.testDeltaR   (v1, t, t1,      "DeltaR   TLorentzVector      " );  
00532       s2=a.testDeltaR   (v2, t, t2,      "DeltaR   XYZTVector          " ); 
00533       s3=a.testDeltaR   (v3, t, t3,      "DeltaR   PtEtaPhiEVector     " ); 
00534 #ifdef WIN32
00535       //windows is bad here 
00536       a.check("DeltaR",s1,s2,s3,10);      
00537 #else
00538       a.check("DeltaR",s1,s2,s3);      
00539 #endif
00540 
00541 
00542       int n1, n2, n3; 
00543       n1 = a.testAnalysis (v1, t, t1, "Analysis1 TLorentzVector     " ); 
00544       n2 = a.testAnalysis (v2, t, t2, "Analysis1 XYZTVector         " ); 
00545       n3 = a.testAnalysis (v3, t, t3, "Analysis1 PtEtaPhiEVector    " ); 
00546       a.check("Analysis1",n1,n2,n3);      
00547 
00548 
00549 
00550       n1 = a.testAnalysis2 (v1, t, t1, "Analysis2 TLorentzVector    " ); 
00551       n2 = a.testAnalysis2 (v2, t, t2, "Analysis2 XYZTVector        " ); 
00552       n3 = a.testAnalysis2 (v3, t, t3, "Analysis2 PtEtaPhiEVector   " ); 
00553       a.check("Analysis2",n1,n2,n3);      
00554 
00555 
00556       // test Rotations on Vectors 
00557       TRotation r1; r1.RotateX(1.); r1.RotateY(2.); r1.RotateZ(3);
00558       Rotation3D  r2 = RotationZ(3)*RotationY(2)*RotationX(1.);
00559       // need to go through LorentzRotation since ROOT does not have ROtation*4D Vector
00560       TLorentzRotation lr1(r1);
00561       LorentzRotation lr2(r2);
00562       // apply also a boost
00563       XYZVector bVec(0.4,0.5,0.6); // bVec.R() must be < 1
00564       assert( bVec.R() <= 1);
00565       lr1.Boost(bVec.x(), bVec.y(), bVec.z());
00566       Boost b(bVec);
00567       LorentzRotation lrb (b);
00568 
00569       lr2 =  lrb * lr2; 
00570 #ifdef DEBUG
00571       // for TLR need to loop 
00572       std::cout << " TLorentzRotation: " << std::endl;
00573       for (int i = 0; i < 4; ++i) {
00574         for (int j = 0; j < 4; ++j)
00575           std::cout << lr1(i,j) << "  ";
00576         std::cout << "\n";
00577       }
00578       std::cout << "\n";
00579       std::cout << "LorentzRotation :\n"  << lr2 << std::endl; 
00580 #endif
00581 
00582       s1=a.testRotation   (v1, lr1, t, t1,      "TRotation  TLorentzVector      " );  
00583       s2=a.testRotation   (v2, lr2, t, t2,      "Rotation3D XYZTVector          " ); 
00584       s3=a.testRotation   (v3, lr2, t, t3,      "Rotation3D PtEtaPhiEVector     " ); 
00585       a.check("Rotation",s1,s2,s3,10);
00586 
00587 
00588       double s0 = s1;
00589       // test rotations using the matrix for multiplications
00590       double rotData[16];
00591       lr2.GetComponents(rotData, rotData+16);
00592       TMatrixD m1(4,4,rotData);
00593       SMatrix<double,4,4>  m2(rotData, rotData+16);
00594 #ifdef DEBUG
00595       m1.Print();
00596       std::cout << m2 << std::endl; 
00597 #endif
00598 
00599       s1=a.testMatVec    (v2, m1, t, t1,        "TMatrix * XYZTVector           " );  
00600       s2=a.testMatVec    (v2, m2, t, t2,        "SMatrix * XYZTVector           " ); 
00601       s3=a.testMatVec    (v3, m2, t, t3,        "SMatrix * PtEtaPhiEVector      " ); 
00602       a.check("Matrix mult",s1,s2,s3,10);
00603 
00604 
00605       // test Boost
00606       TVector3 tVec(bVec.X(), bVec.Y(), bVec.Z() );
00607       s1 = a.testBoost_TL (v1, tVec, t, t1,     "Boost TLorentzVector           ");
00608       s2 = a.testBoost1   (v2, bVec, t, t2,     "Boost XYZTVector               ");
00609       s3 = a.testBoost1   (v3, bVec, t, t3,     "Boost PtEtaPhiEVector          ");
00610       a.check("Boost1",s1,s2,s3,10);
00611 
00612       // test Boost (2)
00613       s0 = s1;
00614       s1 = a.testBoost2  (v1, tVec, t, t1,      "Boost2 TLorentzVector          ");
00615       s2 = a.testBoost2  (v2, bVec, t, t2,      "Boost2 XYZTVector              ");
00616       s3 = a.testBoost2  (v3, bVec, t, t3,      "Boost2 PtEtaPhiEVector         ");
00617       a.check("Boost1-2",s0,s1,s2);
00618       a.check("Boost2",s1,s2,s3,10);
00619 
00620       // test BoostX
00621       double beta = 0.8;
00622       s1 = a.testBoostX2   (v1, beta, t, t1,    "BoostX TLorentzVector          ");
00623       s2 = a.testBoostX1   (v2, beta, t, t2,    "BoostX XYZTVector              ");
00624       s3 = a.testBoostX1   (v3, beta, t, t3,    "BoostX PtEtaPhiEVector         ");
00625       a.check("BoostX1",s1,s2,s3,10);
00626 
00627       // test Boost (2)
00628       s0 = s2;
00629       s1 = a.testBoostX2  (v1, beta, t, t1,     "BoostX2 TLorentzVector         ");
00630       s2 = a.testBoostX2  (v2, beta, t, t2,     "BoostX2 XYZTVector             ");
00631       s3 = a.testBoostX2  (v3, beta, t, t3,     "BoostX2 PtEtaPhiEVector        ");
00632       a.check("BoostX1-2",s0,s1,s2);
00633       a.check("BoostX2",s1,s2,s3,10);
00634      
00635 
00636 
00637       // clean all at the end
00638       a.clear(v1); 
00639       a.clear(v2);
00640       a.clear(v3);
00641 
00642       std::cout << std::endl;
00643       std::cout << "Total Time for  TLorentzVector        = " << t1 << "\t(sec)" << std::endl;
00644       std::cout << "Total Time for  XYZTVector            = " << t2 << "\t(sec)" << std::endl;
00645       std::cout << "Total Time for  PtRhoPhiEtaVector     = " << t3 << "\t(sec)" << std::endl;
00646    }
00647 
00648   //tr.dump(); 
00649 
00650 }
00651 

Generated on Tue Jul 5 15:16:01 2011 for ROOT_528-00b_version by  doxygen 1.5.1