vectorOperation.cxx

Go to the documentation of this file.
00001 // test performance of all vectors operations +,- and * 
00002 
00003 // results on mactelm g++ 4.01 showing ROOT::Math performs best overall
00004 
00005 //v3 = v1+v2 v2 += v1   v3 = v1-v2 v2 -= v1   v2 = a*v1  v1 *= a    v2 = v1/a  v1 /= a   
00006 //    0.59       0.57       0.58       0.56       0.69        0.7       1.65       1.64     2D
00007 //    0.79       0.79       0.78        0.8       0.97       0.95       1.85       1.84     3D 
00008 //    1.07       1.07       1.07       1.07       1.32       1.31       1.72       1.71     4D 
00009 
00010 //  ROOT Physics Vector (TVector's): 
00011 //v3 = v1+v2 v2 += v1   v3 = v1-v2 v2 -= v1   v2 = a*v1  v1 *= a
00012 //    4.4        0.97       4.41       0.96       4.43       1.13          2D
00013 //    5.44       1.25       5.48       1.24       6.12       1.46          3D
00014 //   17.65       7.32      17.65       7.35      10.25       7.79          4D
00015 
00016 //  CLHEP Vector (HepVector's): 
00017 //v3 = v1+v2 v2 += v1   v3 = v1-v2 v2 -= v1   v2 = a*v1  v1 *= a
00018 //    0.57       0.55       0.56       0.55        0.7        0.7                           2D
00019 //     0.8       0.79       0.78       0.77       0.96       0.94        2.7        3.7     3D
00020 //    1.06       1.02       1.06       1.02       1.26       1.26       2.99       3.98     4D
00021 
00022 
00023 
00024 
00025 #include "TRandom2.h"
00026 #include "TStopwatch.h"
00027 
00028 #include <vector>
00029 #include <iostream>
00030 #include <iomanip>
00031 
00032 #ifdef DIM_2
00033 #ifdef USE_POINT
00034 #include "Math/Point2D.h"
00035 typedef ROOT::Math::XYPoint VecType;  
00036 #elif USE_CLHEP
00037 #include "CLHEP/Vector/TwoVector.h"
00038 typedef Hep2Vector VecType;  
00039 #elif USE_ROOT
00040 #include "TVector2.h"
00041 typedef TVector2 VecType;  
00042 #else 
00043 #include "Math/Vector2D.h"
00044 typedef ROOT::Math::XYVector VecType;  
00045 #endif
00046 
00047 #ifndef USE_ROOT
00048 #define VSUM(v) v.x() + v.y() 
00049 #else 
00050 #define VSUM(v) v.X() + v.Y() 
00051 #endif
00052 
00053 
00054 #elif DIM_3 // 3 Dimensions
00055 
00056 #ifdef USE_POINT
00057 #include "Math/Point3D.h"
00058 typedef ROOT::Math::XYZPoint VecType;  
00059 #elif USE_CLHEP
00060 #include "CLHEP/Vector/ThreeVector.h"
00061 typedef Hep3Vector VecType;  
00062 #elif USE_ROOT
00063 #include "TVector3.h"
00064 typedef TVector3 VecType;  
00065 #else 
00066 #include "Math/Vector3D.h"
00067 typedef ROOT::Math::XYZVector VecType;  
00068 #endif
00069 
00070 #ifndef USE_ROOT
00071 #define VSUM(v) v.x() + v.y() + v.z()
00072 #else 
00073 #define VSUM(v) v.X() + v.Y() + v.Z()
00074 #endif
00075 
00076 #else // default is 4D
00077 
00078 #undef USE_POINT
00079 #ifdef USE_CLHEP
00080 #include "CLHEP/Vector/LorentzVector.h"
00081 typedef HepLorentzVector VecType;  
00082 #elif USE_ROOT
00083 #include "TLorentzVector.h"
00084 typedef TLorentzVector VecType;  
00085 #else 
00086 #include "Math/Vector4D.h"
00087 typedef ROOT::Math::XYZTVector VecType;  
00088 #endif
00089 
00090 #ifndef USE_ROOT
00091 #define VSUM(v) v.x() + v.y() + v.z() + v.t()
00092 #else 
00093 #define VSUM(v) v.X() + v.Y() + v.Z() + v.T()
00094 #endif
00095 
00096 #endif
00097 
00098 
00099 const int N = 1000000; 
00100 
00101 template<class Vector>
00102 class TestVector { 
00103 public: 
00104 
00105    TestVector(); 
00106    void Add(); 
00107    void Add2(); 
00108    void Sub(); 
00109    void Sub2(); 
00110    void Scale();
00111    void Scale2();
00112    void Divide();
00113    void Divide2();
00114 
00115    void PrintSummary();
00116 
00117 private: 
00118 
00119    std::vector<Vector> vlist;
00120    std::vector<double> scale;
00121    double fTime[10]; // timing results
00122    int  fTest;
00123 };   
00124 
00125 
00126  
00127 template<class Vector>
00128 TestVector<Vector>::TestVector() : 
00129    vlist(std::vector<Vector>(N) ),
00130    scale(std::vector<double>(N) ), 
00131    fTest(0)
00132 { 
00133 
00134    // create list of vectors and fill them 
00135    
00136    TRandom2 r(111); 
00137    for (int i = 0; i< N; ++i) { 
00138 #ifdef DIM_2
00139       vlist[i] = Vector( r.Uniform(-1,1), r.Uniform(-1,1) ); 
00140 #elif DIM_3
00141       vlist[i] = Vector( r.Uniform(-1,1),r.Uniform(-1,1),r.Uniform(-1,1)  ); 
00142 #else // 4D
00143       vlist[i] = Vector( r.Uniform(-1,1),r.Uniform(-1,1),r.Uniform(-1,1), r.Uniform(2,10) ); 
00144 #endif
00145       scale[i] = r.Uniform(0,1);
00146    }
00147 
00148    std::cout << "test using " << typeid(vlist[0]).name() << std::endl;
00149 }
00150 
00151 template<class Vector>
00152 void TestVector<Vector>::Add()    
00153    // normal addition
00154 { 
00155    TStopwatch w; 
00156    w.Start(); 
00157    double s = 0; 
00158    for (int l = 0; l<100; ++l) { 
00159       for (int i = 1; i< N; ++i) { 
00160          Vector v3 = vlist[i-1] + vlist[i]; 
00161          s += VSUM(v3); 
00162       }
00163    }
00164       
00165    std::cout << "Time for  v3 = v1 + v2 :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00166    std::cout << "value " << s << std::endl << std::endl;
00167    fTime[fTest++] = w.CpuTime();
00168 }
00169 
00170 template<class Vector>
00171 void TestVector<Vector>::Add2()    
00172 { 
00173    // self addition
00174    TStopwatch w; 
00175    Vector v3;
00176    w.Start(); 
00177    double s = 0; 
00178    for (int l = 0; l<100; ++l) { 
00179       for (int i = 0; i< N; ++i) { 
00180          v3 += vlist[i]; 
00181          s += VSUM(v3);
00182       }
00183    }
00184    
00185    std::cout << "Time for  v2 += v1     :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00186    std::cout << "value " << s << std::endl << std::endl;
00187    fTime[fTest++] = w.CpuTime();
00188 }
00189 
00190 template<class Vector>
00191 void TestVector<Vector>::Sub()    
00192 { 
00193    // normal sub
00194    TStopwatch w; 
00195    w.Start(); 
00196    double s = 0; 
00197    for (int l = 0; l<100; ++l) { 
00198       for (int i = 1; i< N; ++i) { 
00199          Vector v3 = vlist[i-1] - vlist[i]; 
00200          s += VSUM(v3);
00201       }
00202    }
00203       
00204    std::cout << "Time for  v3 = v1 - v2 :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00205    std::cout << "value " << s << std::endl << std::endl;
00206    fTime[fTest++] = w.CpuTime();
00207 }
00208 
00209 template<class Vector>
00210 void TestVector<Vector>::Sub2()    
00211 { 
00212    // self subtruction
00213    TStopwatch w; 
00214    Vector v3;
00215    w.Start(); 
00216    double s = 0; 
00217    for (int l = 0; l<100; ++l) { 
00218       for (int i = 0; i< N; ++i) { 
00219          v3 -= vlist[i]; 
00220          s += VSUM(v3);
00221       }
00222    }
00223    
00224    std::cout << "Time for  v2 -= v1     :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00225    std::cout << "value " << s << std::endl << std::endl;
00226    fTime[fTest++] = w.CpuTime();
00227 }
00228 
00229 
00230 template<class Vector>
00231 void TestVector<Vector>::Scale()    
00232 { 
00233 // normal multiply
00234    TStopwatch w; 
00235    w.Start(); 
00236    double s = 0; 
00237    for (int l = 0; l<100; ++l) { 
00238       for (int i = 1; i< N; ++i) { 
00239          Vector v3 = scale[i]*vlist[i]; 
00240          s += VSUM(v3);
00241       }
00242    }
00243       
00244    std::cout << "Time for  v2 = A * v1 :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00245    std::cout << "value " << s << std::endl << std::endl;
00246    fTime[fTest++] = w.CpuTime();
00247 }
00248 
00249 template<class Vector>
00250 void TestVector<Vector>::Scale2()    
00251 { 
00252    // self scale
00253    TStopwatch w; 
00254    Vector v3;
00255    w.Start(); 
00256    double s = 0; 
00257    for (int l = 0; l<100; ++l) { 
00258       for (int i = 0; i< N; ++i) { 
00259          v3 = vlist[i]; 
00260          v3 *= scale[i];
00261          s += VSUM(v3);
00262       }
00263    }
00264       
00265    std::cout << "Time for  v *= a     :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00266    std::cout << "value " << s << std::endl << std::endl;
00267    fTime[fTest++] = w.CpuTime();
00268 }
00269 
00270 template<class Vector>
00271 void TestVector<Vector>::Divide()    
00272 { 
00273 // normal divide
00274    TStopwatch w; 
00275    w.Start(); 
00276    double s = 0; 
00277    for (int l = 0; l<100; ++l) { 
00278       for (int i = 1; i< N; ++i) { 
00279          Vector v3 = vlist[i]/scale[i]; 
00280          s += VSUM(v3);
00281       }
00282    }
00283       
00284    std::cout << "Time for  v2 = v1 / a :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00285    std::cout << "value " << s << std::endl << std::endl;
00286    fTime[fTest++] = w.CpuTime();
00287 }
00288 
00289 template<class Vector> 
00290 void TestVector<Vector>::Divide2()    
00291 { 
00292    // self divide
00293    TStopwatch w; 
00294    Vector v3;
00295    w.Start(); 
00296    double s = 0; 
00297    for (int l = 0; l<100; ++l) { 
00298       for (int i = 0; i< N; ++i) { 
00299          v3 = vlist[i]; 
00300          v3 /= scale[i];
00301          s += VSUM(v3);
00302       }
00303    }
00304    
00305    std::cout << "Time for  v /= a      :\t" << w.RealTime() << "\t" << w.CpuTime() << std::endl;
00306    std::cout << "value " << s << std::endl << std::endl;
00307    fTime[fTest++] = w.CpuTime();
00308 }
00309 
00310 
00311 template<class Vector> 
00312 void TestVector<Vector>::PrintSummary()    
00313 {
00314    std::cout << "\nResults for " << typeid(vlist[0]).name() << std::endl;
00315    std::cout << " v3 = v1+v2"  
00316              << " v2 += v1  "  
00317              << " v3 = v1-v2" 
00318              << " v2 -= v1  " 
00319              << " v2 = a*v1 " 
00320              << " v1 *= a   "
00321              << " v2 = v1/a " 
00322              << " v1 /= a   " << std::endl;
00323 
00324    for (int i = 0; i < fTest; ++i) { 
00325       std::cout << std::setw(8) << fTime[i] << "   "; 
00326    }
00327    std::cout << std::endl << std::endl; 
00328 }
00329 
00330 int main() { 
00331    TestVector<VecType> t; 
00332 #ifndef USE_POINT
00333    t.Add();
00334    t.Add2();
00335    t.Sub();
00336    t.Sub2();
00337 #endif
00338    t.Scale(); 
00339    t.Scale2(); 
00340 #ifndef USE_ROOT
00341    t.Divide();
00342    t.Divide2();
00343 #endif
00344 
00345    // summurize test
00346    t.PrintSummary();
00347 }

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