Logo ROOT   6.07/09
Reference Guide
mathcoreGenVector.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// \notebook -nodraw
4 /// Example macro testing available methods and operation of the GenVector classes.
5 /// The results are compared and check at the
6 /// numerical precision levels.
7 /// Some small discrepancy can appear when the macro
8 /// is executed on different architectures where it has been calibrated (Power PC G5)
9 /// The macro is divided in 4 parts:
10 /// - testVector3D : tests of the 3D Vector classes
11 /// - testPoint3D : tests of the 3D Point classes
12 /// - testLorentzVector : tests of the 4D LorentzVector classes
13 /// - testVectorUtil : tests of the utility functions of all the vector classes
14 ///
15 /// To execute the macro type in:
16 ///
17 /// ~~~{.cpp}
18 /// root[0] .x mathcoreGenVector.C
19 /// ~~~
20 ///
21 /// \macro_output
22 /// \macro_code
23 ///
24 /// \author Lorenzo Moneta
25 
26 #include "TMatrixD.h"
27 #include "TVectorD.h"
28 #include "TMath.h"
29 
30 #include "Math/Point3D.h"
31 #include "Math/Vector3D.h"
32 #include "Math/Vector4D.h"
42 #include "Math/GenVector/Boost.h"
43 #include "Math/GenVector/BoostX.h"
44 #include "Math/GenVector/BoostY.h"
45 #include "Math/GenVector/BoostZ.h"
47 #include "Math/GenVector/Plane3D.h"
49 
50 using namespace ROOT::Math;
51 
52 int ntest = 0;
53 int nfail = 0;
54 int ok = 0;
55 
56 int compare( double v1, double v2, const char* name, double Scale = 1.0) {
57  ntest = ntest + 1;
58 
59  // numerical double limit for epsilon
60  double eps = Scale* 2.22044604925031308e-16;
61  int iret = 0;
62  double delta = v2 - v1;
63  double d = 0;
64  if (delta < 0 ) delta = - delta;
65  if (v1 == 0 || v2 == 0) {
66  if (delta > eps ) {
67  iret = 1;
68  }
69  }
70  // skip case v1 or v2 is infinity
71  else {
72  d = v1;
73 
74  if ( v1 < 0) d = -d;
75  // add also case when delta is small by default
76  if ( delta/d > eps && delta > eps )
77  iret = 1;
78  }
79 
80  if (iret == 0)
81  std::cout << ".";
82  else {
83  int pr = std::cout.precision (18);
84  int discr;
85  if (d != 0)
86  discr = int(delta/d/eps);
87  else
88  discr = int(delta/eps);
89 
90  std::cout << "\nDiscrepancy in " << name << "() : " << v1 << " != " << v2 << " discr = " << discr
91  << " (Allowed discrepancy is " << eps << ")\n";
92  std::cout.precision (pr);
93  nfail = nfail + 1;
94  }
95  return iret;
96 }
97 
98 int testVector3D() {
99  std::cout << "\n************************************************************************\n "
100  << " Vector 3D Test"
101  << "\n************************************************************************\n";
102 
103  gSystem->Load("libGenVector");
104 
105  XYZVector v1(0.01, 0.02, 16);
106 
107  std::cout << "Test Cartesian-Polar : " ;
108 
109  Polar3DVector v2(v1.R(), v1.Theta(), v1.Phi() );
110 
111  ok = 0;
112  ok+= compare(v1.X(), v2.X(), "x");
113  ok+= compare(v1.Y(), v2.Y(), "y");
114  ok+= compare(v1.Z(), v2.Z(), "z");
115  ok+= compare(v1.Phi(), v2.Phi(), "phi");
116  ok+= compare(v1.Theta(), v2.Theta(), "theta");
117  ok+= compare(v1.R(), v2.R(), "r");
118  ok+= compare(v1.Eta(), v2.Eta(), "eta");
119  ok+= compare(v1.Rho(), v2.Rho(), "rho");
120 
121  if (ok == 0) std::cout << "\t OK " << std::endl;
122 
123  std::cout << "Test Cartesian-CylindricalEta : ";
124 
125  RhoEtaPhiVector v3( v1.Rho(), v1.Eta(), v1.Phi() );
126 
127  ok = 0;
128  ok+= compare(v1.X(), v3.X(), "x");
129  ok+= compare(v1.Y(), v3.Y(), "y");
130  ok+= compare(v1.Z(), v3.Z(), "z");
131  ok+= compare(v1.Phi(), v3.Phi(), "phi");
132  ok+= compare(v1.Theta(), v3.Theta(), "theta");
133  ok+= compare(v1.R(), v3.R(), "r");
134  ok+= compare(v1.Eta(), v3.Eta(), "eta");
135  ok+= compare(v1.Rho(), v3.Rho(), "rho");
136 
137  if (ok == 0) std::cout << "\t OK " << std::endl;
138 
139  std::cout << "Test Cartesian-Cylindrical : ";
140 
141  RhoZPhiVector v4( v1.Rho(), v1.Z(), v1.Phi() );
142 
143  ok = 0;
144  ok+= compare(v1.X(), v4.X(), "x");
145  ok+= compare(v1.Y(), v4.Y(), "y");
146  ok+= compare(v1.Z(), v4.Z(), "z");
147  ok+= compare(v1.Phi(), v4.Phi(), "phi");
148  ok+= compare(v1.Theta(), v4.Theta(), "theta");
149  ok+= compare(v1.R(), v4.R(), "r");
150  ok+= compare(v1.Eta(), v4.Eta(), "eta");
151  ok+= compare(v1.Rho(), v4.Rho(), "rho");
152 
153  if (ok == 0) std::cout << "\t OK " << std::endl;
154 
155  std::cout << "Test Operations : " ;
156 
157  ok = 0;
158  double Dot = v1.Dot(v2);
159  ok+= compare( Dot, v1.Mag2(),"dot" );
160  XYZVector vcross = v1.Cross(v2);
161  ok+= compare( vcross.R(), 0,"cross" );
162 
163  XYZVector vscale1 = v1*10;
164  XYZVector vscale2 = vscale1/10;
165  ok+= compare( v1.R(), vscale2.R(), "scale");
166 
167  XYZVector vu = v1.Unit();
168  ok+= compare(v2.Phi(),vu.Phi(),"unit Phi");
169  ok+= compare(v2.Theta(),vu.Theta(),"unit Theta");
170  ok+= compare(1.0,vu.R(),"unit ");
171 
172  XYZVector q1 = v1;
173  RhoEtaPhiVector q2(1.0,1.0,1.0);
174 
175  XYZVector q3 = q1 + q2;
176  XYZVector q4 = q3 - q2;
177 
178  ok+= compare( q4.X(), q1.X(), "op X" );
179  ok+= compare( q4.Y(), q1.Y(), "op Y" );
180  ok+= compare( q4.Z(), q1.Z(), "op Z" );
181 
182  // test operator ==
183  XYZVector w1 = v1;
184  Polar3DVector w2 = v2;
185  RhoEtaPhiVector w3 = v3;
186  RhoZPhiVector w4 = v4;
187  ok+= compare( w1 == v1, static_cast<double>(true), "== XYZ");
188  ok+= compare( w2 == v2, static_cast<double>(true), "== Polar");
189  ok+= compare( w3 == v3, static_cast<double>(true), "== RhoEtaPhi");
190  ok+= compare( w4 == v4, static_cast<double>(true), "== RhoZPhi");
191 
192  if (ok == 0) std::cout << "\t OK " << std::endl;
193 
194  //test setters
195  std::cout << "Test Setters : " ;
196 
197  q2.SetXYZ(q1.X(), q1.Y(), q1.Z() );
198 
199  ok+= compare( q2.X(), q1.X(), "setXYZ X" );
200  ok+= compare( q2.Y(), q1.Y(), "setXYZ Y" );
201  ok+= compare( q2.Z(), q1.Z(), "setXYZ Z" );
202 
203  q2.SetCoordinates( 2.0*q1.Rho(), q1.Eta(), q1.Phi() );
204  XYZVector q1s = 2.0*q1;
205  ok+= compare( q2.X(), q1s.X(), "set X" );
206  ok+= compare( q2.Y(), q1s.Y(), "set Y" );
207  ok+= compare( q2.Z(), q1s.Z(), "set Z" );
208 
209  if (ok == 0) std::cout << "\t\t OK " << std::endl;
210 
211  std::cout << "Test Linear Algebra conversion: " ;
212 
213  XYZVector vxyz1(1.,2.,3.);
214 
215  TVectorD vla1(3);
216  vxyz1.Coordinates().GetCoordinates(vla1.GetMatrixArray() );
217 
218  TVectorD vla2(3);
219  vla2[0] = 1.; vla2[1] = -2.; vla2[2] = 1.;
220 
221  XYZVector vxyz2;
222  vxyz2.SetCoordinates(&vla2[0]);
223 
224  ok = 0;
225  double prod1 = vxyz1.Dot(vxyz2);
226  double prod2 = vla1*vla2;
227  ok+= compare( prod1, prod2, "la test" );
228 
229  if (ok == 0) std::cout << "\t\t OK " << std::endl;
230 
231  return ok;
232 }
233 
234 int testPoint3D() {
235  std::cout << "\n************************************************************************\n "
236  << " Point 3D Tests"
237  << "\n************************************************************************\n";
238 
239  XYZPoint p1(1.0, 2.0, 3.0);
240 
241  std::cout << "Test Cartesian-Polar : ";
242 
243  Polar3DPoint p2(p1.R(), p1.Theta(), p1.Phi() );
244 
245  ok = 0;
246  ok+= compare(p1.x(), p2.X(), "x");
247  ok+= compare(p1.y(), p2.Y(), "y");
248  ok+= compare(p1.z(), p2.Z(), "z");
249  ok+= compare(p1.phi(), p2.Phi(), "phi");
250  ok+= compare(p1.theta(), p2.Theta(), "theta");
251  ok+= compare(p1.r(), p2.R(), "r");
252  ok+= compare(p1.eta(), p2.Eta(), "eta");
253  ok+= compare(p1.rho(), p2.Rho(), "rho");
254 
255  if (ok == 0) std::cout << "\t OK " << std::endl;
256 
257  std::cout << "Test Polar-CylindricalEta : ";
258 
259  RhoEtaPhiPoint p3( p2.Rho(), p2.Eta(), p2.Phi() );
260 
261  ok = 0;
262  ok+= compare(p2.X(), p3.X(), "x");
263  ok+= compare(p2.Y(), p3.Y(), "y");
264  ok+= compare(p2.Z(), p3.Z(), "z",3);
265  ok+= compare(p2.Phi(), p3.Phi(), "phi");
266  ok+= compare(p2.Theta(), p3.Theta(), "theta");
267  ok+= compare(p2.R(), p3.R(), "r");
268  ok+= compare(p2.Eta(), p3.Eta(), "eta");
269  ok+= compare(p2.Rho(), p3.Rho(), "rho");
270 
271  if (ok == 0) std::cout << "\t OK " << std::endl;
272 
273  std::cout << "Test operations : ";
274 
275  //std::cout << "\nTest Dot and Cross products with Vectors : ";
276  Polar3DVector vperp(1.,p1.Theta() + TMath::PiOver2(),p1.Phi() );
277  double Dot = p1.Dot(vperp);
278  ok+= compare( Dot, 0.0,"dot", 10 );
279 
280  XYZPoint vcross = p1.Cross(vperp);
281  ok+= compare( vcross.R(), p1.R(),"cross mag" );
282  ok+= compare( vcross.Dot(vperp), 0.0,"cross dir" );
283 
284  XYZPoint pscale1 = 10*p1;
285  XYZPoint pscale2 = pscale1/10;
286  ok+= compare( p1.R(), pscale2.R(), "scale");
287 
288  // test operator ==
289  ok+= compare( p1 == pscale2, static_cast<double>(true), "== Point");
290 
291  //RhoEtaPhiPoint q1 = p1; ! constructor yet not working in CINT
292  RhoEtaPhiPoint q1; q1 = p1;
293  q1.SetCoordinates(p1.Rho(),2.0, p1.Phi() );
294 
295  Polar3DVector v2(p1.R(), p1.Theta(),p1.Phi());
296 
297  if (ok == 0) std::cout << "\t OK " << std::endl;
298 
299  return ok;
300 }
301 
302 int testLorentzVector() {
303  std::cout << "\n************************************************************************\n "
304  << " Lorentz Vector Tests"
305  << "\n************************************************************************\n";
306 
307  XYZTVector v1(1.0, 2.0, 3.0, 4.0);
308 
309  std::cout << "Test XYZT - PtEtaPhiE Vectors: ";
310 
311  PtEtaPhiEVector v2( v1.Rho(), v1.Eta(), v1.Phi(), v1.E() );
312 
313  ok = 0;
314  ok+= compare(v1.Px(), v2.X(), "x");
315  ok+= compare(v1.Py(), v2.Y(), "y");
316  ok+= compare(v1.Pz(), v2.Z(), "z", 2);
317  ok+= compare(v1.E(), v2.T(), "e");
318  ok+= compare(v1.Phi(), v2.Phi(), "phi");
319  ok+= compare(v1.Theta(), v2.Theta(), "theta");
320  ok+= compare(v1.Pt(), v2.Pt(), "pt");
321  ok+= compare(v1.M(), v2.M(), "mass", 5);
322  ok+= compare(v1.Et(), v2.Et(), "et");
323  ok+= compare(v1.Mt(), v2.Mt(), "mt", 3);
324 
325  if (ok == 0) std::cout << "\t OK " << std::endl;
326 
327  std::cout << "Test XYZT - PtEtaPhiM Vectors: ";
328 
329  PtEtaPhiMVector v3( v1.Rho(), v1.Eta(), v1.Phi(), v1.M() );
330 
331  ok = 0;
332  ok+= compare(v1.Px(), v3.X(), "x");
333  ok+= compare(v1.Py(), v3.Y(), "y");
334  ok+= compare(v1.Pz(), v3.Z(), "z", 2);
335  ok+= compare(v1.E(), v3.T(), "e");
336  ok+= compare(v1.Phi(), v3.Phi(), "phi");
337  ok+= compare(v1.Theta(), v3.Theta(), "theta");
338  ok+= compare(v1.Pt(), v3.Pt(), "pt");
339  ok+= compare(v1.M(), v3.M(), "mass", 5);
340  ok+= compare(v1.Et(), v3.Et(), "et");
341  ok+= compare(v1.Mt(), v3.Mt(), "mt", 3);
342 
343  if (ok == 0) std::cout << "\t OK " << std::endl;
344 
345  std::cout << "Test PtEtaPhiE - PxPyPzM Vect.: ";
346 
347  PxPyPzMVector v4( v3.X(), v3.Y(), v3.Z(), v3.M() );
348 
349  ok = 0;
350  ok+= compare(v4.Px(), v3.X(), "x");
351  ok+= compare(v4.Py(), v3.Y(), "y");
352  ok+= compare(v4.Pz(), v3.Z(), "z",2);
353  ok+= compare(v4.E(), v3.T(), "e");
354  ok+= compare(v4.Phi(), v3.Phi(), "phi");
355  ok+= compare(v4.Theta(), v3.Theta(), "theta");
356  ok+= compare(v4.Pt(), v3.Pt(), "pt");
357  ok+= compare(v4.M(), v3.M(), "mass",5);
358  ok+= compare(v4.Et(), v3.Et(), "et");
359  ok+= compare(v4.Mt(), v3.Mt(), "mt",3);
360 
361  if (ok == 0) std::cout << "\t OK " << std::endl;
362 
363  std::cout << "Test operations : ";
364 
365  ok = 0;
366  double Dot = v1.Dot(v2);
367  ok+= compare( Dot, v1.M2(),"dot" , 10 );
368 
369  XYZTVector vscale1 = v1*10;
370  XYZTVector vscale2 = vscale1/10;
371  ok+= compare( v1.M(), vscale2.M(), "scale");
372 
373  XYZTVector q1 = v1;
374  PtEtaPhiEVector q2(1.0,1.0,1.0,5.0);
375 
376  XYZTVector q3 = q1 + q2;
377  XYZTVector q4 = q3 - q2;
378 
379  ok+= compare( q4.x(), q1.X(), "op X" );
380  ok+= compare( q4.y(), q1.Y(), "op Y" );
381  ok+= compare( q4.z(), q1.Z(), "op Z" );
382  ok+= compare( q4.t(), q1.E(), "op E" );
383 
384  // test operator ==
385  XYZTVector w1 = v1;
386  PtEtaPhiEVector w2 = v2;
387  PtEtaPhiMVector w3 = v3;
388  PxPyPzMVector w4 = v4;
389  ok+= compare( w1 == v1, static_cast<double>(true), "== PxPyPzE");
390  ok+= compare( w2 == v2, static_cast<double>(true), "== PtEtaPhiE");
391  ok+= compare( w3 == v3, static_cast<double>(true), "== PtEtaPhiM");
392  ok+= compare( w4 == v4, static_cast<double>(true), "== PxPyPzM");
393 
394  // test gamma beta and boost
395  XYZVector b = q1.BoostToCM();
396  double beta = q1.Beta();
397  double gamma = q1.Gamma();
398 
399  ok += compare( b.R(), beta, "beta" );
400  ok += compare( gamma, 1./sqrt( 1 - beta*beta ), "gamma");
401 
402  if (ok == 0) std::cout << "\t OK " << std::endl;
403 
404  //test setters
405  std::cout << "Test Setters : " ;
406 
407  q2.SetXYZT(q1.Px(), q1.Py(), q1.Pz(), q1.E() );
408 
409  ok+= compare( q2.X(), q1.X(), "setXYZT X" );
410  ok+= compare( q2.Y(), q1.Y(), "setXYZT Y" );
411  ok+= compare( q2.Z(), q1.Z(), "setXYZT Z" ,2);
412  ok+= compare( q2.T(), q1.E(), "setXYZT E" );
413 
414  q2.SetCoordinates( 2.0*q1.Rho(), q1.Eta(), q1.Phi(), 2.0*q1.E() );
415  XYZTVector q1s = q1*2.0;
416  ok+= compare( q2.X(), q1s.X(), "set X" );
417  ok+= compare( q2.Y(), q1s.Y(), "set Y" );
418  ok+= compare( q2.Z(), q1s.Z(), "set Z" ,2);
419  ok+= compare( q2.T(), q1s.T(), "set E" );
420 
421  if (ok == 0) std::cout << "\t OK " << std::endl;
422 
423  return ok;
424 }
425 
426 int testVectorUtil() {
427  std::cout << "\n************************************************************************\n "
428  << " Utility Function Tests"
429  << "\n************************************************************************\n";
430 
431  std::cout << "Test Vector utility functions : ";
432 
433  XYZVector v1(1.0, 2.0, 3.0);
434  Polar3DVector v2pol(v1.R(), v1.Theta()+TMath::PiOver2(), v1.Phi() + 1.0);
435  // mixedmethods not yet impl.
436  XYZVector v2; v2 = v2pol;
437 
438  ok = 0;
439  ok += compare( VectorUtil::DeltaPhi(v1,v2), 1.0, "deltaPhi Vec");
440 
441  RhoEtaPhiVector v2cyl(v1.Rho(), v1.Eta() + 1.0, v1.Phi() + 1.0);
442  v2 = v2cyl;
443 
444  ok += compare( VectorUtil::DeltaR(v1,v2), sqrt(2.0), "DeltaR Vec");
445 
446  XYZVector vperp = v1.Cross(v2);
447  ok += compare( VectorUtil::CosTheta(v1,vperp), 0.0, "costheta Vec");
448  ok += compare( VectorUtil::Angle(v1,vperp), TMath::PiOver2(), "angle Vec");
449 
450  if (ok == 0) std::cout << "\t\t OK " << std::endl;
451 
452 
453  std::cout << "Test Point utility functions : ";
454 
455  XYZPoint p1(1.0, 2.0, 3.0);
456  Polar3DPoint p2pol(p1.R(), p1.Theta()+TMath::PiOver2(), p1.Phi() + 1.0);
457  // mixedmethods not yet impl.
458  XYZPoint p2; p2 = p2pol;
459 
460  ok = 0;
461  ok += compare( VectorUtil::DeltaPhi(p1,p2), 1.0, "deltaPhi Point");
462 
463  RhoEtaPhiPoint p2cyl(p1.Rho(), p1.Eta() + 1.0, p1.Phi() + 1.0);
464  p2 = p2cyl;
465  ok += compare( VectorUtil::DeltaR(p1,p2), sqrt(2.0), "DeltaR Point");
466 
467  XYZPoint pperp(vperp.X(), vperp.Y(), vperp.Z());
468  ok += compare( VectorUtil::CosTheta(p1,pperp), 0.0, "costheta Point");
469  ok += compare( VectorUtil::Angle(p1,pperp), TMath::PiOver2(), "angle Point");
470 
471  if (ok == 0) std::cout << "\t\t OK " << std::endl;
472 
473  std::cout << "LorentzVector utility funct.: ";
474 
475  XYZTVector q1(1.0, 2.0, 3.0,4.0);
476  PtEtaPhiEVector q2cyl(q1.Pt(), q1.Eta()+1.0, q1.Phi() + 1.0, q1.E() );
477  XYZTVector q2; q2 = q2cyl;
478 
479  ok = 0;
480  ok += compare( VectorUtil::DeltaPhi(q1,q2), 1.0, "deltaPhi LVec");
481  ok += compare( VectorUtil::DeltaR(q1,q2), sqrt(2.0), "DeltaR LVec");
482 
483  XYZTVector qsum = q1+q2;
484  ok += compare( VectorUtil::InvariantMass(q1,q2), qsum.M(), "InvMass");
485 
486  if (ok == 0) std::cout << "\t\t OK " << std::endl;
487 
488  return ok;
489 }
490 
491 int testRotation() {
492  std::cout << "\n************************************************************************\n "
493  << " Rotation and Transformation Tests"
494  << "\n************************************************************************\n";
495 
496  std::cout << "Test Vector Rotations : ";
497  ok = 0;
498 
499  XYZPoint v(1.,2,3.);
500 
501  double pi = TMath::Pi();
502  // initiate rotation with some non -trivial angles to test all matrix
503  EulerAngles r1( pi/2.,pi/4., pi/3 );
504  Rotation3D r2(r1);
505  // only operator= is in CINT for the other rotations
506  Quaternion r3; r3 = r2;
507  AxisAngle r4; r4 = r3;
508  RotationZYX r5; r5 = r2;
509 
510  XYZPoint v1 = r1 * v;
511  XYZPoint v2 = r2 * v;
512  XYZPoint v3 = r3 * v;
513  XYZPoint v4 = r4 * v;
514  XYZPoint v5 = r5 * v;
515 
516  ok+= compare(v1.X(), v2.X(), "x",2);
517  ok+= compare(v1.Y(), v2.Y(), "y",2);
518  ok+= compare(v1.Z(), v2.Z(), "z",2);
519 
520  ok+= compare(v1.X(), v3.X(), "x",2);
521  ok+= compare(v1.Y(), v3.Y(), "y",2);
522  ok+= compare(v1.Z(), v3.Z(), "z",2);
523 
524  ok+= compare(v1.X(), v4.X(), "x",5);
525  ok+= compare(v1.Y(), v4.Y(), "y",5);
526  ok+= compare(v1.Z(), v4.Z(), "z",5);
527 
528  ok+= compare(v1.X(), v5.X(), "x",2);
529  ok+= compare(v1.Y(), v5.Y(), "y",2);
530  ok+= compare(v1.Z(), v5.Z(), "z",2);
531 
532  // test with matrix
533  double rdata[9];
534  r2.GetComponents(rdata, rdata+9);
535  TMatrixD m(3,3,rdata);
536  double vdata[3];
537  v.GetCoordinates(vdata);
538  TVectorD q(3,vdata);
539  TVectorD q2 = m*q;
540 
541  XYZPoint v6;
542  v6.SetCoordinates( q2.GetMatrixArray() );
543 
544  ok+= compare(v1.X(), v6.X(), "x");
545  ok+= compare(v1.Y(), v6.Y(), "y");
546  ok+= compare(v1.Z(), v6.Z(), "z");
547 
548  if (ok == 0) std::cout << "\t OK " << std::endl;
549  else std::cout << std::endl;
550 
551  std::cout << "Test Axial Rotations : ";
552  ok = 0;
553 
554  RotationX rx( pi/3);
555  RotationY ry( pi/4);
556  RotationZ rz( 4*pi/5);
557 
558  Rotation3D r3x(rx);
559  Rotation3D r3y(ry);
560  Rotation3D r3z(rz);
561 
562  Quaternion qx; qx = rx;
563  Quaternion qy; qy = ry;
564  Quaternion qz; qz = rz;
565 
566  RotationZYX rzyx( rz.Angle(), ry.Angle(), rx.Angle() );
567 
568  XYZPoint vrot1 = rx * ry * rz * v;
569  XYZPoint vrot2 = r3x * r3y * r3z * v;
570 
571  ok+= compare(vrot1.X(), vrot2.X(), "x");
572  ok+= compare(vrot1.Y(), vrot2.Y(), "y");
573  ok+= compare(vrot1.Z(), vrot2.Z(), "z");
574 
575  vrot2 = qx * qy * qz * v;
576 
577  ok+= compare(vrot1.X(), vrot2.X(), "x",2);
578  ok+= compare(vrot1.Y(), vrot2.Y(), "y",2);
579  ok+= compare(vrot1.Z(), vrot2.Z(), "z",2);
580 
581  vrot2 = rzyx * v;
582 
583  ok+= compare(vrot1.X(), vrot2.X(), "x");
584  ok+= compare(vrot1.Y(), vrot2.Y(), "y");
585  ok+= compare(vrot1.Z(), vrot2.Z(), "z");
586 
587  // now inverse (first x then y then z)
588  vrot1 = rz * ry * rx * v;
589  vrot2 = r3z * r3y * r3x * v;
590 
591  ok+= compare(vrot1.X(), vrot2.X(), "x");
592  ok+= compare(vrot1.Y(), vrot2.Y(), "y");
593  ok+= compare(vrot1.Z(), vrot2.Z(), "z");
594 
595  XYZPoint vinv1 = rx.Inverse()*ry.Inverse()*rz.Inverse()*vrot1;
596 
597  ok+= compare(vinv1.X(), v.X(), "x",2);
598  ok+= compare(vinv1.Y(), v.Y(), "y");
599  ok+= compare(vinv1.Z(), v.Z(), "z");
600 
601  if (ok == 0) std::cout << "\t OK " << std::endl;
602  else std::cout << std::endl;
603 
604 
605  std::cout << "Test Rotations by a PI angle : ";
606  ok = 0;
607 
608  double b[4] = { 6,8,10,3.14159265358979323 };
609  AxisAngle arPi(b,b+4 );
610  Rotation3D rPi(arPi);
611  AxisAngle a1; a1 = rPi;
612  ok+= compare(arPi.Axis().X(), a1.Axis().X(),"x");
613  ok+= compare(arPi.Axis().Y(), a1.Axis().Y(),"y");
614  ok+= compare(arPi.Axis().Z(), a1.Axis().Z(),"z");
615  ok+= compare(arPi.Angle(), a1.Angle(),"angle");
616 
617  EulerAngles ePi; ePi=rPi;
618  EulerAngles e1; e1=Rotation3D(a1);
619  ok+= compare(ePi.Phi(), e1.Phi(),"phi");
620  ok+= compare(ePi.Theta(), e1.Theta(),"theta");
621  ok+= compare(ePi.Psi(), e1.Psi(),"ps1");
622 
623  if (ok == 0) std::cout << "\t\t OK " << std::endl;
624  else std::cout << std::endl;
625 
626  std::cout << "Test Inversions : ";
627  ok = 0;
628 
629  EulerAngles s1 = r1.Inverse();
630  Rotation3D s2 = r2.Inverse();
631  Quaternion s3 = r3.Inverse();
632  AxisAngle s4 = r4.Inverse();
633  RotationZYX s5 = r5.Inverse();
634 
635  // Euler angles not yet impl.
636  XYZPoint p = s2 * r2 * v;
637 
638  ok+= compare(p.X(), v.X(), "x",10);
639  ok+= compare(p.Y(), v.Y(), "y",10);
640  ok+= compare(p.Z(), v.Z(), "z",10);
641 
642 
643  p = s3 * r3 * v;
644 
645  ok+= compare(p.X(), v.X(), "x",10);
646  ok+= compare(p.Y(), v.Y(), "y",10);
647  ok+= compare(p.Z(), v.Z(), "z",10);
648 
649  p = s4 * r4 * v;
650  // axis angle inversion not very precise
651  ok+= compare(p.X(), v.X(), "x",1E9);
652  ok+= compare(p.Y(), v.Y(), "y",1E9);
653  ok+= compare(p.Z(), v.Z(), "z",1E9);
654 
655  p = s5 * r5 * v;
656 
657  ok+= compare(p.X(), v.X(), "x",10);
658  ok+= compare(p.Y(), v.Y(), "y",10);
659  ok+= compare(p.Z(), v.Z(), "z",10);
660 
661 
662  Rotation3D r6(r5);
663  Rotation3D s6 = r6.Inverse();
664 
665  p = s6 * r6 * v;
666 
667  ok+= compare(p.X(), v.X(), "x",10);
668  ok+= compare(p.Y(), v.Y(), "y",10);
669  ok+= compare(p.Z(), v.Z(), "z",10);
670 
671  if (ok == 0) std::cout << "\t OK " << std::endl;
672  else std::cout << std::endl;
673 
674  // test Rectify
675  std::cout << "Test rectify : ";
676  ok = 0;
677 
678  XYZVector u1(0.999498,-0.00118212,-0.0316611);
679  XYZVector u2(0,0.999304,-0.0373108);
680  XYZVector u3(0.0316832,0.0372921,0.998802);
681  Rotation3D rr(u1,u2,u3);
682  // check orto-normality
683  XYZPoint vrr = rr* v;
684  ok+= compare(v.R(), vrr.R(), "R",1.E9);
685 
686  if (ok == 0) std::cout << "\t\t OK " << std::endl;
687  else std::cout << std::endl;
688 
689  std::cout << "Test Transform3D : ";
690  ok = 0;
691 
692  XYZVector d(1.,-2.,3.);
693  Transform3D t(r2,d);
694 
695  XYZPoint pd = t * v;
696  // apply directly rotation
697  XYZPoint vd = r2 * v + d;
698 
699  ok+= compare(pd.X(), vd.X(), "x");
700  ok+= compare(pd.Y(), vd.Y(), "y");
701  ok+= compare(pd.Z(), vd.Z(), "z");
702 
703  // test with matrix
704  double tdata[12];
705  t.GetComponents(tdata);
706  TMatrixD mt(3,4,tdata);
707  double vData[4]; // needs a vector of dim 4
708  v.GetCoordinates(vData);
709  vData[3] = 1;
710  TVectorD q0(4,vData);
711 
712  TVectorD qt = mt*q0;
713 
714  ok+= compare(pd.X(), qt(0), "x");
715  ok+= compare(pd.Y(), qt(1), "y");
716  ok+= compare(pd.Z(), qt(2), "z");
717 
718  // test inverse
719 
720  Transform3D tinv = t.Inverse();
721 
722  p = tinv * t * v;
723 
724  ok+= compare(p.X(), v.X(), "x",10);
725  ok+= compare(p.Y(), v.Y(), "y",10);
726  ok+= compare(p.Z(), v.Z(), "z",10);
727 
728  // test construct inverse from translation first
729 
730  Transform3D tinv2 ( r2.Inverse(), r2.Inverse() *( -d) ) ;
731  p = tinv2 * t * v;
732 
733  ok+= compare(p.X(), v.X(), "x",10);
734  ok+= compare(p.Y(), v.Y(), "y",10);
735  ok+= compare(p.Z(), v.Z(), "z",10);
736 
737  // test from only rotation and only translation
738  Transform3D ta( EulerAngles(1.,2.,3.) );
739  Transform3D tb( XYZVector(1,2,3) );
740  Transform3D tc( Rotation3D(EulerAngles(1.,2.,3.)) , XYZVector(1,2,3) );
741  Transform3D td( ta.Rotation(), ta.Rotation() * XYZVector(1,2,3) ) ;
742 
743  ok+= compare( tc == tb*ta, static_cast<double>(true), "== Rot*Tra");
744  ok+= compare( td == ta*tb, static_cast<double>(true), "== Rot*Tra");
745 
746 
747  if (ok == 0) std::cout << "\t OK " << std::endl;
748  else std::cout << std::endl;
749 
750  std::cout << "Test Plane3D : ";
751  ok = 0;
752 
753  // test transfrom a 3D plane
754 
755  XYZPoint p1(1,2,3);
756  XYZPoint p2(-2,-1,4);
757  XYZPoint p3(-1,3,2);
758  Plane3D plane(p1,p2,p3);
759 
760  XYZVector n = plane.Normal();
761  // normal is perpendicular to vectors on the planes obtained from subtracting the points
762  ok+= compare(n.Dot(p2-p1), 0.0, "n.v12",10);
763  ok+= compare(n.Dot(p3-p1), 0.0, "n.v13",10);
764  ok+= compare(n.Dot(p3-p2), 0.0, "n.v23",10);
765 
766  Plane3D plane1 = t(plane);
767 
768  // transform the points
769  XYZPoint pt1 = t(p1);
770  XYZPoint pt2 = t(p2);
771  XYZPoint pt3 = t(p3);
772  Plane3D plane2(pt1,pt2,pt3);
773 
774  XYZVector n1 = plane1.Normal();
775  XYZVector n2 = plane2.Normal();
776 
777  ok+= compare(n1.X(), n2.X(), "a",10);
778  ok+= compare(n1.Y(), n2.Y(), "b",10);
779  ok+= compare(n1.Z(), n2.Z(), "c",10);
780  ok+= compare(plane1.HesseDistance(), plane2.HesseDistance(), "d",10);
781 
782  // check distances
783  ok += compare(plane1.Distance(pt1), 0.0, "distance",10);
784 
785  if (ok == 0) std::cout << "\t OK " << std::endl;
786  else std::cout << std::endl;
787 
788  std::cout << "Test LorentzRotation : ";
789  ok = 0;
790 
791  XYZTVector lv(1.,2.,3.,4.);
792 
793  // test from rotx (using boosts and 3D rotations not yet impl.)
794  // rx,ry and rz already defined
795  Rotation3D r3d = rx*ry*rz;
796 
797  LorentzRotation rlx(rx);
798  LorentzRotation rly(ry);
799  LorentzRotation rlz(rz);
800 
801  LorentzRotation rl0 = rlx*rly*rlz;
802  LorentzRotation rl1( r3d);
803 
804  XYZTVector lv0 = rl0 * lv;
805 
806  XYZTVector lv1 = rl1 * lv;
807 
808  XYZTVector lv2 = r3d * lv;
809 
810  ok+= compare(lv1== lv2,true,"V0==V2");
811  ok+= compare(lv1== lv2,true,"V1==V2");
812 
813  double rlData[16];
814  rl0.GetComponents(rlData);
815  TMatrixD ml(4,4,rlData);
816  double lvData[4];
817  lv.GetCoordinates(lvData);
818  TVectorD ql(4,lvData);
819 
820  TVectorD qlr = ml*ql;
821 
822  ok+= compare(lv1.X(), qlr(0), "x");
823  ok+= compare(lv1.Y(), qlr(1), "y");
824  ok+= compare(lv1.Z(), qlr(2), "z");
825  ok+= compare(lv1.E(), qlr(3), "t");
826 
827  // test inverse
828  lv0 = rl0 * rl0.Inverse() * lv;
829 
830  ok+= compare(lv0.X(), lv.X(), "x");
831  ok+= compare(lv0.Y(), lv.Y(), "y");
832  ok+= compare(lv0.Z(), lv.Z(), "z");
833  ok+= compare(lv0.E(), lv.E(), "t");
834 
835  if (ok == 0) std::cout << "\t OK " << std::endl;
836  else std::cout << std::endl;
837 
838  // test Boosts
839  std::cout << "Test Boost : ";
840  ok = 0;
841 
842  Boost bst( 0.3,0.4,0.5); // boost (must be <= 1)
843 
844  XYZTVector lvb = bst ( lv );
845 
846  LorentzRotation rl2 (bst);
847 
848  XYZTVector lvb2 = rl2 (lv);
849 
850  // test with lorentz rotation
851  ok+= compare(lvb.X(), lvb2.X(), "x");
852  ok+= compare(lvb.Y(), lvb2.Y(), "y");
853  ok+= compare(lvb.Z(), lvb2.Z(), "z");
854  ok+= compare(lvb.E(), lvb2.E(), "t");
855  ok+= compare(lvb.M(), lv.M(), "m",50); // m must stay constant
856 
857  // test inverse
858  lv0 = bst.Inverse() * lvb;
859 
860  ok+= compare(lv0.X(), lv.X(), "x",5);
861  ok+= compare(lv0.Y(), lv.Y(), "y",5);
862  ok+= compare(lv0.Z(), lv.Z(), "z",3);
863  ok+= compare(lv0.E(), lv.E(), "t",3);
864 
865  XYZVector brest = lv.BoostToCM();
866  bst.SetComponents( brest.X(), brest.Y(), brest.Z() );
867 
868  XYZTVector lvr = bst * lv;
869 
870  ok+= compare(lvr.X(), 0.0, "x",10);
871  ok+= compare(lvr.Y(), 0.0, "y",10);
872  ok+= compare(lvr.Z(), 0.0, "z",10);
873  ok+= compare(lvr.M(), lv.M(), "m",10);
874 
875  if (ok == 0) std::cout << "\t OK " << std::endl;
876  else std::cout << std::endl;
877  return ok;
878 }
879 
880 void mathcoreGenVector() {
881 
882 
883  testVector3D();
884  testPoint3D();
885  testLorentzVector();
886  testVectorUtil();
887  testRotation();
888 
889  std::cout << "\n\nNumber of tests " << ntest << " failed = " << nfail << std::endl;
890 }
T Dot(const SVector< T, D > &lhs, const SVector< T, D > &rhs)
Vector dot product.
Definition: Functions.h:166
Quaternion Inverse() const
Return inverse of a rotation.
Definition: Quaternion.h:259
Class describing a geometrical plane in 3 dimensions.
Definition: Plane3D.h:47
Scalar E() const
return 4-th component (time, or energy for a 4-momentum vector)
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
Definition: LorentzVector.h:54
int testVector3D()
DisplacementVector3D< CoordSystem, Tag > & SetCoordinates(const Scalar src[])
Set internal data based on a C-style array of 3 Scalar numbers.
Scalar Phi() const
Return Phi Euler angle.
Definition: EulerAngles.h:212
static double p3(double t, double a, double b, double c, double d)
RotationZYX Inverse() const
Return inverse of a rotation.
Definition: RotationZYX.h:283
const double pi
Rotation class representing a 3D rotation about the Z axis by the angle of rotation.
Definition: RotationZ.h:43
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
Scalar R() const
Polar R, converting if necessary from internal coordinate system.
int compare(double v1, double v2, const std::string &name="", double scale=1.0)
Class describing a generic position vector (point) in 3 dimensions.
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1818
Scalar M() const
return magnitude (mass) using the (-,-,-,+) metric.
int testPoint3D()
Rotation class with the (3D) rotation represented by angles describing first a rotation of an angle p...
Definition: RotationZYX.h:71
AxisVector Axis() const
accesss to rotation axis
Definition: AxisAngle.h:183
Rotation class with the (3D) rotation represented by a unit quaternion (u, i, j, k).
Definition: Quaternion.h:47
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
Scalar Dot(const DisplacementVector3D< OtherCoords, Tag > &v) const
Return the scalar (dot) product of two displacement vectors.
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
double beta(double x, double y)
Calculates the beta function.
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
AxisAngle class describing rotation represented with direction axis (3D Vector) and an angle of rotat...
Definition: AxisAngle.h:41
Rotation3D Inverse() const
Return inverse of a rotation.
Definition: Rotation3D.h:420
Rotation class representing a 3D rotation about the Y axis by the angle of rotation.
Definition: RotationY.h:43
unsigned int r3[N_CITIES]
Definition: simanTSP.cxx:323
Vector1::Scalar DeltaR(const Vector1 &v1, const Vector2 &v2)
Find difference in pseudorapidity (Eta) and Phi betwen two generic vectors The only requirements on t...
Definition: VectorUtil.h:97
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
Lorentz boost class with the (4D) transformation represented internally by a 4x4 orthosymplectic matr...
Definition: Boost.h:46
static double p2(double t, double a, double b, double c)
Vector1::Scalar InvariantMass(const Vector1 &v1, const Vector2 &v2)
return the invariant mass of two LorentzVector The only requirement on the LorentzVector is that they...
Definition: VectorUtil.h:217
PositionVector3D< CoordSystem, Tag > & SetCoordinates(const Scalar src[])
Set internal data based on a C-style array of 3 Scalar numbers.
Scalar R() const
Polar R, converting if necessary from internal coordinate system.
Scalar Angle() const
access to rotation angle
Definition: AxisAngle.h:188
AxisAngle Inverse() const
Return inverse of an AxisAngle rotation.
Definition: AxisAngle.h:266
Element * GetMatrixArray()
Definition: TVectorT.h:84
Lorentz transformation class with the (4D) transformation represented by a 4x4 orthosymplectic matrix...
double gamma(double x)
Vector1::Scalar DeltaPhi(const Vector1 &v1, const Vector2 &v2)
Find aximutal Angle difference between two generic vectors ( v2.Phi() - v1.Phi() ) The only requireme...
Definition: VectorUtil.h:61
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
SVector< double, 2 > v
Definition: Dict.h:5
double Angle(const Vector1 &v1, const Vector2 &v2)
Find Angle between two vectors.
Definition: VectorUtil.h:140
LorentzRotation Inverse() const
Return inverse of a rotation.
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
Rotation class representing a 3D rotation about the X axis by the angle of rotation.
Definition: RotationX.h:43
Rotation class with the (3D) rotation represented by a 3x3 orthogonal matrix.
Definition: Rotation3D.h:65
TMarker * m
Definition: textangle.C:8
void GetComponents(Foreign4Vector &v1, Foreign4Vector &v2, Foreign4Vector &v3, Foreign4Vector &v4) const
Get components into four 4-vectors which will be the (orthosymplectic) columns of the rotation matrix...
Basic 3D Transformation class describing a rotation and then a translation The internal data are a 3D...
Definition: Transform3D.h:85
static double p1(double t, double a, double b)
Scalar Theta() const
Return Theta Euler angle.
Definition: EulerAngles.h:222
Double_t Pi()
Definition: TMath.h:44
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
EulerAngles class describing rotation as three angles (Euler Angles).
Definition: EulerAngles.h:43
Scalar Dot(const DisplacementVector3D< OtherCoords, Tag > &v) const
Return the scalar (Dot) product of this with a displacement vector in any coordinate system...
DisplacementVector3D< Cartesian3D< double >, DefaultCoordinateSystemTag > XYZVector
3D Vector based on the cartesian coordinates x,y,z in double precision
Definition: Vector3Dfwd.h:34
DisplacementVector3D Cross(const DisplacementVector3D< OtherCoords, Tag > &v) const
Return vector (cross) product of two displacement vectors, as a vector in the coordinate system of th...
Double_t PiOver2()
Definition: TMath.h:46
Transform3D Inverse() const
Return the inverse of the transformation.
Definition: Transform3D.h:576
int testVectorUtil()
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
TLine * lv
Definition: textalign.C:5
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
DisplacementVector3D Unit() const
return unit vector parallel to this
float * q
Definition: THbookFile.cxx:87
const Int_t n
Definition: legend1.C:16
double CosTheta(const Vector1 &v1, const Vector2 &v2)
Find CosTheta Angle between two generic 3D vectors pre-requisite: vectors implement the X()...
Definition: VectorUtil.h:114
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
char name[80]
Definition: TGX11.cxx:109
Scalar Psi() const
Return Psi Euler angle.
Definition: EulerAngles.h:232