ROOT  6.06/09
Reference Guide
stress2D.cxx
Go to the documentation of this file.
1 // @(#)root/test:$Id$
2 // Author: Lorenzo Moneta 06/2005
3 ///////////////////////////////////////////////////////////////////////////////////
4 //
5 // MathCore Benchmark test suite
6 // ==============================
7 //
8 // This program performs tests of ROOT::Math 4D LorentzVectors comparing with TLorentzVector
9 // The time performing various vector operations on a collection of vectors is measured.
10 // The benchmarked operations are:
11 // - vector construction from 4 values
12 // - construction using a setter method
13 // - simple addition of all the vector pairs in the collection
14 // - calculation of deltaR = phi**2 + eta**2 of all vector pairs in the collection
15 // - two simple analysis:
16 // - the first requires some cut (on pt and eta) and on the invariant mass
17 // of the selected pairs
18 // - the second requires just some cut in pt, eta and delta R on all the
19 // vector pair
20 // - conversion between XYZTVectors to PtRhoEtaPhi based vectors
21 //
22 // The two analysis demostrates, especially in the second case, the advantage of using
23 // vector based on cylindrical coordinate, given the fact that the time spent in the conversion is
24 // much less than the time spent in the analysis routine.
25 //
26 // To run the program do:
27 // stressVector : run standard test with collection of 1000 vectors
28 // stressVector 10000 : run with a collection of 10000 vectors
29 //
30 ///////////////////////////////////////////////////////////////////////////////////
31 
32 
33 
34 #include <vector>
35 #include <iostream>
36 #include <algorithm>
37 
38 #include <assert.h>
39 #include <map>
40 
41 #include "TStopwatch.h"
42 
43 #include "TRandom3.h"
44 #include "TVector2.h"
45 
46 #include "Math/Vector2D.h"
47 #include "Math/Point2D.h"
48 #include "Math/SVector.h"
49 
50 #include <cmath>
51 
52 #include "limits"
53 
54 #ifdef HAVE_CLHEP
55 #include "CLHEP/Vector/TwoVector.h"
56 #endif
57 
58 
59 //#define DEBUG
60 
61 using namespace ROOT::Math;
62 
63 
64 
65 
66 class VectorTest {
67 
68 private:
69 
70 // global data variables
71  std::vector<double> dataX;
72  std::vector<double> dataY;
73 
74  size_t nGen;
75  size_t n2Loop ;
76 
77 
78 public:
79 
80  VectorTest(int n1, int n2) :
81  nGen(n1),
82  n2Loop(n2)
83  {}
84 
85 
86 
87 
88  void print(TStopwatch & time, std::string s) {
89  int pr = std::cout.precision(8);
90  std::cout << s << "\t" << " time = " << time.RealTime() << "\t(sec)\t"
91  // << time.CpuTime()
92  << std::endl;
93  std::cout.precision(pr);
94  }
95 
96 
97  int check(std::string name, double s1, double s2, double s3, double scale=1) {
98  double eps = 10*scale*std::numeric_limits<double>::epsilon();
99  if ( std::fabs(s1-s2) < eps*std::fabs(s1) && std::fabs(s1-s3) < eps*std::fabs(s1) ) return 0;
100  int pr = std::cout.precision(16);
101  std::cout << s1 << "\t" << s2 <<"\t" << s3 << "\n";
102  std::cout << "Rel Diff 1-2:\t" << (s1-s2)/std::fabs(s1) << " Diff 1-3:\t" << (s1-s3)/std::fabs(s1) << std::endl;
103  std::cout << "Test " << name << " failed !!\n\n";
104  std::cout.precision(pr);
105  return -1;
106  }
107 
108 
109  void genData() {
110  int n = nGen;
111 
112  // generate 2 d data
113  TRandom3 rdm;
114  for (int i = 0; i < n ; ++ i) {
115 
116  double phi = rdm.Rndm()*3.1415926535897931;
117  double r = rdm.Exp(10.);
118  // fill vectors
119 
120  Polar2DVector q( r, phi);
121  dataX.push_back( q.x() );
122  dataY.push_back( q.y() );
123 
124  }
125  }
126 
127 
128 
129  template <class V>
130  void testCreate( std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
131 
132  int n = dataX.size();
133  dataV.resize(n);
134  tim.Start();
135  for (int i = 0; i < n; ++i) {
136  dataV[i] = new V( dataX[i], dataY[i] );
137  }
138  tim.Stop();
139  t += tim.RealTime();
140  print(tim,s);
141  }
142 
143 
144  template <class V>
145  void testCreate2( std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
146 
147  int n = dataX.size();
148  dataV.resize(n);
149  tim.Start();
150  for (int i = 0; i < n; ++i) {
151  dataV[i] = new V();
152  dataV[i]->SetXY(dataX[i], dataY[i] );
153  }
154  tim.Stop();
155  print(tim,s);
156  t += tim.RealTime();
157  }
158  void testCreate2( std::vector<TVector2 *> & dataV, TStopwatch & tim, double& t, std::string s) {
159 
160  int n = dataX.size();
161  dataV.resize(n);
162  tim.Start();
163  for (int i = 0; i < n; ++i) {
164  dataV[i] = new TVector2();
165  dataV[i]->Set(dataX[i], dataY[i] );
166  }
167  tim.Stop();
168  print(tim,s);
169  t += tim.RealTime();
170  }
171 #ifdef HAVE_CLHEP
172  void testCreate2( std::vector<Hep2Vector *> & dataV, TStopwatch & tim, double& t, std::string s) {
173 
174  int n = dataX.size();
175  dataV.resize(n);
176  tim.Start();
177  for (int i = 0; i < n; ++i) {
178  dataV[i] = new Hep2Vector();
179  dataV[i]->set(dataX[i], dataY[i] );
180  }
181  tim.Stop();
182  print(tim,s);
183  t += tim.RealTime();
184  }
185 #endif
186 
187 
188  template <class V>
189  void clear( std::vector<V *> & dataV ) {
190  for (unsigned int i = 0; i < dataV.size(); ++i) {
191  V * p = dataV[i];
192  delete p;
193  }
194  dataV.clear();
195 
196 }
197 
198 template<class V>
199 inline double addXY(const V & v) {
200  return v.X() + v.Y();
201 }
202 inline double addXY(const SVector<double,2> & v) {
203  return v(0) + v(1);
204 }
205 template<class V>
206 inline double getSum(const V & v1, const V & v2) {
207  return v1.X()+v1.Y() + v2.X() + v2.Y();
208 }
209 
210 inline double getSum(const SVector<double,2> & v1, const SVector<double,2> & v2 ) {
211  return v1(0)+v1(1) + v2(0)+v2(1);
212 }
213 
214 template<class V>
215 inline double dotProd(const V & v1, const V & v2) {
216  return v1 * v2;
217 }
218 
219 inline double dotProd(const XYVector & v1, const XYVector & v2) {
220  return v1.Dot(v2);
221 }
222 
223 inline double dotProd(const SVector<double,2> & v1, const SVector<double,2> & v2 ) {
224  return Dot(v1,v2);
225 }
226 
227 
228 #ifdef HAVE_CLHEP
229 inline double addXY(const Hep2Vector & v) {
230  return v.x() + v.y();
231 }
232 inline double getSum(const Hep2Vector & v1, const Hep2Vector & v2 ) {
233  return v1.x() + v1.y() + v2.x() + v2.y();
234 }
235 #endif
236 
237 template <class V>
238 double testAddition( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
239  unsigned int n = dataV.size();
240  double tot = 0;
241  V v0 = *(dataV[0]);
242  tim.Start();
243  for (unsigned int i = 0; i < n; ++i) {
244  V & v1 = *(dataV[i]);
245  V v3 = v1 + v0;
246  tot += addXY(v3);
247  }
248  tim.Stop();
249  print(tim,s);
250  t += tim.RealTime();
251  return tot;
252 }
253 
254 template <class V>
255 double testAddition2( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
256  unsigned int n = dataV.size();
257  double tot = 0;
258  V v0 = *(dataV[0]);
259  tim.Start();
260  for (unsigned int i = 0; i < n; ++i) {
261  const V & v1 = *(dataV[i]);
262  v0 += v1;
263  tot += addXY(v0);
264  }
265  tim.Stop();
266  print(tim,s);
267  t += tim.RealTime();
268  return tot;
269 }
270 
271 template <class V>
272 double testAddition3( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
273  //unsigned int n = std::min(n2Loop, dataV.size() );
274  unsigned int n = dataV.size();
275  double tot = 0;
276  V v0 = *(dataV[0]);
277  tim.Start();
278  for (unsigned int i = 0; i < n; ++i) {
279  V & v1 = *(dataV[i]);
280 // for (unsigned int j = i +1; j < n; ++j) {
281 // V & v2 = *(dataV[j]);
282 // tot += getSum(v1,v2);
283 // }
284  tot += getSum(v1,v0);
285  }
286  tim.Stop();
287  print(tim,s);
288  t += tim.RealTime();
289  return tot;
290 }
291 
292 
293 template <class V>
294 double testDotProduct( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
295  //unsigned int n = std::min(n2Loop, dataV.size() );
296  double tot = 0;
297  unsigned int n = dataV.size();
298  V v0 = *(dataV[0]);
299  tim.Start();
300  for (unsigned int i = 0; i < n; ++i) {
301  V & v1 = *(dataV[i]);
302 // for (unsigned int j = i +1; j < n; ++j) {
303 // V & v2 = *(dataV[j]);
304 // tot += dotProd(v1,v2);
305 // }
306  tot += dotProd(v1,v0);
307  }
308  tim.Stop();
309  print(tim,s);
310  t += tim.RealTime();
311  return tot;
312 }
313 
314 
315 template <class V>
316 double testScale( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
317  unsigned int n = dataV.size();
318  double tot = 0;
319  tim.Start();
320  for (unsigned int i = 0; i < n; ++i) {
321  V & v1 = *(dataV[i]);
322  // scale
323  V v2 = 2.0*v1;
324  tot += addXY(v2);
325  }
326  tim.Stop();
327  print(tim,s);
328  t += tim.RealTime();
329  return tot;
330 }
331 
332 template <class V>
333 double testScale2( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
334  //unsigned int n = std::min(n2Loop, dataV.size() );
335  unsigned int n = dataV.size();
336  double tot = 0;
337  tim.Start();
338  for (unsigned int i = 0; i < n; ++i) {
339  V & v1 = *(dataV[i]);
340  // scale
341  v1 *= 2.0;
342  tot += addXY(v1);
343  }
344  tim.Stop();
345  print(tim,s);
346  t += tim.RealTime();
347  return tot;
348 }
349 
350 template <class V>
351 double testOperations( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
352  //unsigned int n = std::min(n2Loop, dataV.size() );
353  // test operations like in CMS
354  unsigned int n = dataV.size();
355  double tot = 0;
356  V v0a = *(dataV[0]);
357  V v0b = *(dataV[n-1]);
358  tim.Start();
359  for (unsigned int i = 0; i < n; ++i) {
360  V & v1 = *(dataV[i]);
361  //V v2(v1 - dotProd(v1,v0a)*v0b );
362  double a = dotProd(v1,v0a);
363  V v2(v1 - a*v0b );
364  tot += addXY(v2);
365  }
366  tim.Stop();
367  print(tim,s);
368  t += tim.RealTime();
369  return tot;
370 }
371 
372 
373 template<class V>
374 inline double dPhi(V & v1, V& v2) {
375  return std::abs(v1.Phi() - v2.Phi() );
376 }
377 
378 #ifdef HAVE_CLHEP
379 inline double dPhi(Hep2Vector & v1, Hep2Vector & v2) {
380  return std::abs(v1.phi() - v2.phi() );
381 }
382 #endif
383 
384 
385 template <class V>
386 double testDeltaPhi( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
387  unsigned int n = std::min(n2Loop, dataV.size() );
388  tim.Start();
389  double tot = 0;
390  for (unsigned int i = 0; i < n; ++i) {
391  V & v1 = *(dataV[i]);
392  for (unsigned int j = i +1; j < n; ++j) {
393  V & v2 = *(dataV[j]);
394  double delta = dPhi(v1,v2);
395  tot += delta;
396  }
397  }
398  tim.Stop();
399  print(tim,s);
400  t += tim.RealTime();
401  return tot;
402 }
403 
404 
405 // template <class V>
406 // int testAnalysis( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
407 // int nsel = 0;
408 // int nsel2 = 0;
409 // double deltaMax = 1.;
410 // double ptMin = 1.;
411 // double etaMax = 3.;
412 
413 // unsigned int n = std::min(n2Loop, dataV.size() );
414 // tim.Start();
415 // for (unsigned int i = 0; i < n; ++i) {
416 // V & v1 = *(dataV[i]);
417 // if (cutPtEta(v1,ptMin, etaMax) ) {
418 // double delta;
419 // for (unsigned int j = i +1; j < n; ++j) {
420 // V & v2 = *(dataV[j]);
421 // delta = VectorUtil::DeltaR(v1,v2);
422 // if (delta < deltaMax) {
423 // V v3 = v1 + v2;
424 // nsel++;
425 // if ( cutPtEtaAndMass(v3))
426 // nsel2++;
427 // }
428 // }
429 // }
430 // }
431 // tim.Stop();
432 // print(tim,s);
433 // //std::cout << nsel << "\n";
434 // t += tim.RealTime();
435 // return nsel2;
436 // }
437 
438 
439 
440 // template <class V>
441 // int testAnalysis2( const std::vector<V *> & dataV, TStopwatch & tim, double& t, std::string s) {
442 // int nsel = 0;
443 // double ptMin = 1.;
444 // double etaMax = 3.;
445 // unsigned int n = std::min(n2Loop, dataV.size() );
446 // tim.Start();
447 // //seal::SealTimer t(tim.name(), true, std::cout);
448 // for (unsigned int i = 0; i < n; ++i) {
449 // V & v1 = *(dataV[i]);
450 // if ( cutPtEta(v1, ptMin, etaMax) ) {
451 // for (unsigned int j = i +1; j < n; ++j) {
452 // V & v2 = *(dataV[j]);
453 // if ( VectorUtil::DeltaR(v1,v2) < 0.5) nsel++;
454 // }
455 // }
456 // }
457 // tim.Stop();
458 // print(tim,s);
459 // t += tim.RealTime();
460 // return nsel;
461 // }
462 
463 
464 
465  template <class V1, class V2>
466  void testConversion( std::vector<V1 *> & dataV1, std::vector<V2 *> & dataV2, TStopwatch & tim, double& t, std::string s) {
467 
468  int n = dataX.size();
469  dataV2.resize(n);
470  tim.Start();
471  for (int i = 0; i < n; ++i) {
472  dataV2[i] = new V2( *dataV1[i] );
473  }
474  tim.Stop();
475  print(tim,s);
476  t += tim.RealTime();
477  }
478 
479 
480 
481  // rotation
482  template <class V, class R>
483  double testRotation( std::vector<V *> & dataV, double rotAngle, TStopwatch & tim, double& t, std::string s) {
484 
485  unsigned int n = std::min(n2Loop, dataV.size() );
486  tim.Start();
487  double sum = 0;
488  for (unsigned int i = 0; i < n; ++i) {
489  V & v1 = *(dataV[i]);
490  V v2 = v1;
491  v2.Rotate(rotAngle);
492  sum += addXY(v2);
493  }
494  tim.Stop();
495  print(tim,s);
496  t += tim.RealTime();
497  return sum;
498  }
499 
500 
501 
502 };
503 
504 
505 
506 int main(int argc,const char *argv[]) {
507 
508  int ngen = 1000000;
509  if (argc > 1) ngen = atoi(argv[1]);
510  int nloop2 = int(std::sqrt(2.0*ngen)+0.5);
511  if (argc > 2) nloop2 = atoi(argv[2]);
512 
513  std::cout << "Test with Ngen = " << ngen << " n2loop = " << nloop2 << std::endl;
514 
515 
516  TStopwatch t;
517 
518  VectorTest a(ngen,nloop2);
519 
520  a.genData();
521 
522  int niter = 1;
523  for (int i = 0; i < niter; ++i) {
524 
525 #ifdef DEBUG
526  std::cout << "iteration " << i << std::endl;
527 #endif
528 
529  double t1 = 0;
530  double t2 = 0;
531  double t3 = 0;
532  double t4 = 0;
533  double t5 = 0;
534  double t6 = 0;
535 
536  std::vector<TVector2 *> v1;
537  std::vector<XYVector *> v2;
538  std::vector<Polar2DVector *> v3;
539  std::vector<XYPoint *> v4;
540  std::vector<Polar2DPoint *> v5;
541  std::vector<SVector<double,2> *> v6;
542 
543  a.testCreate (v1, t, t1, "creation TVector2 " );
544  a.testCreate (v2, t, t2, "creation XYVector " );
545  a.testCreate (v3, t, t3, "creation Polar2DVector " );
546  a.testCreate (v4, t, t4, "creation XYPoint " );
547  a.testCreate (v5, t, t5, "creation Polar2DPoint " );
548  a.testCreate (v6, t, t6, "creation SVector<2> " );
549 #ifdef HAVE_CLHEP
550  double t7 = 0;
551  std::vector<Hep2Vector *> v7;
552  a.testCreate (v7, t, t7, "creation Hep2Vector " );
553 #endif
554 
555 
556  a.clear(v3);
557  a.clear(v4);
558  a.clear(v5);
559 
560 #ifdef HAVE_CLHEP
561  a.clear(v7);
562 #endif
563 
564  std::cout << "\n";
565  a.testConversion (v2, v3, t, t3, "Conversion XY->Polar " );
566  a.testConversion (v2, v4, t, t4, "Conversion XYVec->XYPoint " );
567  a.testConversion (v2, v5, t, t5, "Conversion XYVec->PolarP " );
568 
569  a.clear(v1);
570  a.clear(v2);
571  a.clear(v3);
572  a.clear(v4);
573  a.clear(v5);
574  std::cout << "\n";
575 
576  a.testCreate2 (v1, t, t1, "creationSet TVector2 " );
577  a.testCreate2 (v2, t, t2, "creationSet XYVector " );
578  a.testCreate2 (v3, t, t3, "creationSet Polar2DVector " );
579  a.testCreate2 (v4, t, t4, "creationSet XYPoint " );
580  a.testCreate2 (v5, t, t5, "creationSet Polar2DPoint " );
581 // a.testCreate2 (v6, t, t6, "creationSet Polar2DPoint " );
582 #ifdef HAVE_CLHEP
583  a.testCreate2 (v7, t, t7, "creationSet Hep2Vector " );
584 #endif
585 
586  std::cout << "\n";
587 
588  double s1,s2,s3,s4,s5,s6;
589  s1=a.testAddition (v1, t, t1, "Addition TVector2 " );
590  s2=a.testAddition (v2, t, t2, "Addition XYVector " );
591  s3=a.testAddition (v3, t, t3, "Addition Polar2DVector " );
592  s6=a.testAddition (v6, t, t6, "Addition SVector<2> " );
593  a.check("Addition test1",s1,s2,s3);
594  a.check("Addition test2",s1,s2,s6);
595 #ifdef HAVE_CLHEP
596  double s7;
597  s7=a.testAddition (v7, t, t7, "Addition Hep2Vector " );
598  a.check("Addition",s7,s1,s2);
599 #endif
600 
601  double s0 = s2;
602  std::cout << "\n";
603 
604  s1=a.testAddition2 (v1, t, t1, "Addition2 TVector2 " );
605  s2=a.testAddition2 (v2, t, t2, "Addition2 XYVector " );
606  s3=a.testAddition2 (v3, t, t3, "Addition2 Polar2DVector " );
607  s6=a.testAddition2 (v6, t, t6, "Addition2 SVector<2> " );
608  a.check("Addition2 test1",s1,s2,s3,100);
609  a.check("Addition2 test2",s1,s2,s6);
610 #ifdef HAVE_CLHEP
611  s7=a.testAddition2 (v7, t, t7, "Addition2 Hep2Vector " );
612  a.check("Addition2 CLHEP",s7,s1,s2);
613 #endif
614 
615  std::cout << "\n";
616 
617  s1=a.testAddition3 (v1, t, t1, "Addition3 TVector2 " );
618  s2=a.testAddition3 (v2, t, t2, "Addition3 XYVector " );
619  s3=a.testAddition3 (v3, t, t3, "Addition3 Polar2DVector " );
620  s6=a.testAddition3 (v6, t, t6, "Addition3 SVector<2> " );
621  a.check("Addition3 test1",s1,s2,s3);
622  a.check("Addition3 test2",s6,s0,s2);
623 #ifdef HAVE_CLHEP
624  s7=a.testAddition3 (v7, t, t7, "Addition3 Hep2Vector " );
625  a.check("Addition3 CLHEP",s7,s1,s2);
626 #endif
627 
628  std::cout << "\n";
629 
630  s1=a.testDotProduct (v1, t, t1, "DotProduct TVector2 " );
631  s2=a.testDotProduct (v2, t, t2, "DotProduct XYVector " );
632 // s3=a.testDotProduct (v3, t, t3, "DotProduct Polar2DVector " );
633  s6=a.testDotProduct (v6, t, t6, "DotProduct SVector<2> " );
634  a.check("DotProduct test1",s1,s2,s6);
635 // a.check("DotProduct test2",s6,s1,s2);
636 #ifdef HAVE_CLHEP
637  s7=a.testDotProduct (v7, t, t7, "DotProduct Hep2Vector " );
638  a.check("DotProduct CLHEP",s7,s1,s2);
639 #endif
640 
641 
642  std::cout << "\n";
643 
644  s1=a.testDeltaPhi (v1, t, t1, "DeltaPhi TVector2 " );
645  s2=a.testDeltaPhi (v2, t, t2, "DeltaPhi XYVector " );
646  s3=a.testDeltaPhi (v3, t, t3, "DeltaPhi Polar2DVector " );
647  s4=a.testDeltaPhi (v4, t, t4, "DeltaPhi XYPoint " );
648  s5=a.testDeltaPhi (v5, t, t5, "DeltaPhi Polar2DPoint " );
649 #ifdef WIN32
650  //windows is bad here
651  a.check("DeltaPhi",s1,s2,s3,10);
652  a.check("DeltaPhi",s2,s4,s5,10);
653 #else
654  a.check("DeltaPhi",s1,s2,s3);
655  a.check("DeltaPhi",s2,s4,s5);
656 #endif
657 #ifdef HAVE_CLHEP
658  s7=a.testDeltaPhi (v7, t, t7, "DeltaPhi HEP2Vector " );
659  a.check("DeltaPhi",s7,s1,s2);
660 #endif
661 
662  std::cout << "\n";
663  s1=a.testScale (v1, t, t1, "Scale of TVector2 " );
664  s2=a.testScale (v2, t, t2, "Scale of XYVector " );
665  s3=a.testScale (v3, t, t3, "Scale of Polar2DVector " );
666  s4=a.testScale (v4, t, t4, "Scale of XYPoint " );
667  s5=a.testScale (v5, t, t5, "Scale of Polar2DPoint " );
668  a.check("Scaling",s1,s2,s3);
669  a.check("Scaling",s2,s4,s5, 10);
670  s6=a.testScale (v6, t, t6, "Scale of SVector<2> " );
671  a.check("Scaling SV",s6,s1,s2);
672 
673 #ifdef HAVE_CLHEP
674  s7=a.testScale (v7, t, t7, "Scale of HEP2Vector " );
675  a.check("Scaling CLHEP",s7,s2,s3);
676 #endif
677 
678  std::cout << "\n";
679  s1=a.testScale2 (v1, t, t1, "Scale2 of TVector2 " );
680  s2=a.testScale2 (v2, t, t2, "Scale2 of XYVector " );
681  s3=a.testScale2 (v3, t, t3, "Scale2 of Polar2DVector " );
682  s4=a.testScale2 (v4, t, t4, "Scale2 of XYPoint " );
683  s5=a.testScale2 (v5, t, t5, "Scale2 of Polar2DPoint " );
684  a.check("Scaling2",s1,s2,s3);
685  a.check("Scaling2",s2,s4,s5, 10);
686  s6=a.testScale2 (v6, t, t6, "Scale2 of SVector<2> " );
687  a.check("Scaling2 SV",s6,s1,s2);
688 
689 #ifdef HAVE_CLHEP
690  s7=a.testScale2 (v7, t, t7, "Scale2 of HEP2Vector " );
691  a.check("Scaling CLHEP",s7,s2,s3);
692 #endif
693 
694  std::cout << "\n";
695 
696  s1=a.testOperations (v1, t, t1, "Operations of TVector2 " );
697  s2=a.testOperations (v2, t, t2, "Operations of XYVector " );
698  s6=a.testOperations (v6, t, t6, "Operations of SVector<2> " );
699  a.check("Operations testSV",s6,s1,s2);
700 #ifdef HAVE_CLHEP
701  s7=a.testOperations (v7, t, t7, "Operations of HEP2Vector " );
702  a.check("Operations CLHEP",s7,s1,s2);
703 #endif
704 
705 
706 
707 #ifdef LATER
708 
709  int n1, n2, n3,n4,n5;
710  n1 = a.testAnalysis (v1, t, t1, "Analysis1 TVector2 " );
711  n2 = a.testAnalysis (v2, t, t2, "Analysis1 XYVector " );
712  n3 = a.testAnalysis (v3, t, t3, "Analysis1 Polar2DVector " );
713  n4 = a.testAnalysis (v4, t, t4, "Analysis1 XYPoint " );
714  n5 = a.testAnalysis (v5, t, t5, "Analysis1 Polar2DPoint " );
715  a.check("Analysis1",n1,n2,n3);
716  a.check("Analysis1",n2,n4,n5);
717 #ifdef HAVE_CLHEP
718  int n6;
719  n6 = a.testAnalysis (v7, t, t7, "Analysis1 HEP2Vector " );
720  a.check("Analysis2 CLHEP",n6,n1,n2);
721 #endif
722 
723 
724  n1 = a.testAnalysis2 (v1, t, t1, "Analysis2 TVector2 " );
725  n2 = a.testAnalysis2 (v2, t, t2, "Analysis2 XYVector " );
726  n3 = a.testAnalysis2 (v3, t, t3, "Analysis2 Polar2DVector " );
727  n4 = a.testAnalysis2 (v4, t, t4, "Analysis2 XYPoint " );
728  n5 = a.testAnalysis2 (v5, t, t5, "Analysis2 Polar2DPoint " );
729  a.check("Analysis2",n1,n2,n3);
730  a.check("Analysis2",n2,n4,n5);
731 #ifdef HAVE_CLHEP
732  n6 = a.testAnalysis2 (v7, t, t7, "Analysis2 HEP2Vector " );
733  a.check("Analysis2 CLHEP",n6,n1,n2);
734 #endif
735 
736 
737  n1 = a.testAnalysis3 (v1, t, t1, "Analysis3 TVector2 " );
738  n2 = a.testAnalysis3 (v2, t, t2, "Analysis3 XYVector " );
739  n3 = a.testAnalysis3 (v3, t, t3, "Analysis3 Polar2DVector " );
740  n4 = a.testAnalysis3 (v4, t, t4, "Analysis3 XYPoint " );
741  n5 = a.testAnalysis3 (v5, t, t5, "Analysis3 Polar2DPoint " );
742  a.check("Analysis3",n1,n2,n3);
743  a.check("Analysis3",n2,n4,n5);
744 #ifdef HAVE_CLHEP
745  n6 = a.testAnalysis3 (v7, t, t7,"Analysis3 HEP2Vector " );
746  a.check("Analysis3 CLHEP",n6,n1,n2);
747 #endif
748 
749 #endif
750 
751 
752  // clean all at the end
753  a.clear(v1);
754  a.clear(v2);
755  a.clear(v3);
756 
757  std::cout << std::endl;
758  std::cout << "Total Time for TVector2 = " << t1 << "\t(sec)" << std::endl;
759  std::cout << "Total Time for XYVector = " << t2 << "\t(sec)" << std::endl;
760  std::cout << "Total Time for Polar2DVector = " << t3 << "\t(sec)" << std::endl;
761 #ifdef HAVE_CLHEP
762  std::cout << "Total Time for Hep2Vector = " << t7 << "\t(sec)" << std::endl;
763 #endif
764  }
765 
766  //tr.dump();
767 
768 }
769 
770 
771 
T Dot(const SVector< T, D > &lhs, const SVector< T, D > &rhs)
Vector dot product.
Definition: Functions.h:166
Random number generator class based on M.
Definition: TRandom3.h:29
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
Definition: TStopwatch.cxx:108
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:56
const Double_t * v1
Definition: TArcBall.cxx:33
TArc * a
Definition: textangle.C:12
TLatex * t1
Definition: textangle.C:20
double sqrt(double)
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:75
static Vc_ALWAYS_INLINE Vector< T > abs(const Vector< T > &x)
Definition: vector.h:450
virtual Double_t Rndm(Int_t i=0)
Machine independent random number generator.
Definition: TRandom3.cxx:93
int testOperations()
int getSum(const int *x, int n)
Definition: GoFTest.cxx:531
Scalar Dot(const DisplacementVector2D< OtherCoords, Tag > &v) const
Return the scalar (dot) product of two displacement vectors.
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
SVector< double, 2 > v
Definition: Dict.h:5
REAL epsilon
Definition: triangle.c:617
int main(int argc, const char *argv[])
Definition: stress2D.cxx:506
float * q
Definition: THbookFile.cxx:87
const Int_t n
Definition: legend1.C:16
Class describing a generic displacement vector in 2 dimensions.
virtual Double_t Exp(Double_t tau)
Returns an exponential deviate.
Definition: TRandom.cxx:212
SVector: a generic fixed size Vector class.
Stopwatch class.
Definition: TStopwatch.h:30