Logo ROOT  
Reference Guide
VectorUtil.h
Go to the documentation of this file.
1 // @(#)root/mathcore:$Id: 9ef2a4a7bd1b62c1293920c2af2f64791c75bdd8 $
2 // Authors: W. Brown, M. Fischler, L. Moneta 2005
3 
4 
5 /**********************************************************************
6  * *
7  * Copyright (c) 2005 , LCG ROOT MathLib Team *
8  * *
9  * *
10  **********************************************************************/
11 
12 // Header file for Vector Utility functions
13 //
14 // Created by: moneta at Tue May 31 21:10:29 2005
15 //
16 // Last update: Tue May 31 21:10:29 2005
17 //
18 #ifndef ROOT_Math_GenVector_VectorUtil
19 #define ROOT_Math_GenVector_VectorUtil 1
20 
21 #include "Math/Math.h"
22 
23 
24 #include "Math/GenVector/Boost.h"
25 
26 namespace ROOT {
27 
28  namespace Math {
29 
30 
31  // utility functions for vector classes
32 
33 
34 
35  /**
36  Global Helper functions for generic Vector classes. Any Vector classes implementing some defined member functions,
37  like Phi() or Eta() or mag() can use these functions.
38  The functions returning a scalar value, returns always double precision number even if the vector are
39  based on another precision type
40 
41  @ingroup GenVector
42  */
43 
44 
45  namespace VectorUtil {
46 
47 
48  // methods for 3D vectors
49 
50  /**
51  Find aximutal Angle difference between two generic vectors ( v2.Phi() - v1.Phi() )
52  The only requirements on the Vector classes is that they implement the Phi() method
53  \param v1 Vector of any type implementing the Phi() operator
54  \param v2 Vector of any type implementing the Phi() operator
55  \return Phi difference
56  \f[ \Delta \phi = \phi_2 - \phi_1 \f]
57  */
58  template <class Vector1, class Vector2>
59  inline typename Vector1::Scalar DeltaPhi( const Vector1 & v1, const Vector2 & v2) {
60  typename Vector1::Scalar dphi = v2.Phi() - v1.Phi();
61  if ( dphi > M_PI ) {
62  dphi -= 2.0*M_PI;
63  } else if ( dphi <= -M_PI ) {
64  dphi += 2.0*M_PI;
65  }
66  return dphi;
67  }
68 
69 
70 
71  /**
72  Find square of the difference in pseudorapidity (Eta) and Phi betwen two generic vectors
73  The only requirements on the Vector classes is that they implement the Phi() and Eta() method
74  \param v1 Vector 1
75  \param v2 Vector 2
76  \return Angle between the two vectors
77  \f[ \Delta R2 = ( \Delta \phi )^2 + ( \Delta \eta )^2 \f]
78  */
79  template <class Vector1, class Vector2>
80  inline typename Vector1::Scalar DeltaR2( const Vector1 & v1, const Vector2 & v2) {
81  typename Vector1::Scalar dphi = DeltaPhi(v1,v2);
82  typename Vector1::Scalar deta = v2.Eta() - v1.Eta();
83  return dphi*dphi + deta*deta;
84  }
85 
86  /**
87  Find square of the difference in true rapidity (y) and Phi betwen two generic vectors
88  The only requirements on the Vector classes is that they implement the Phi() and Rapidity() method
89  \param v1 Vector 1
90  \param v2 Vector 2
91  \return Angle between the two vectors
92  \f[ \Delta R2 = ( \Delta \phi )^2 + ( \Delta \y )^2 \f]
93  */
94  template <class Vector1, class Vector2>
95  inline typename Vector1::Scalar DeltaR2RapidityPhi( const Vector1 & v1, const Vector2 & v2) {
96  typename Vector1::Scalar dphi = DeltaPhi(v1,v2);
97  typename Vector1::Scalar drap = v2.Rapidity() - v1.Rapidity();
98  return dphi*dphi + drap*drap;
99  }
100 
101  /**
102  Find difference in pseudorapidity (Eta) and Phi betwen two generic vectors
103  The only requirements on the Vector classes is that they implement the Phi() and Eta() method
104  \param v1 Vector 1
105  \param v2 Vector 2
106  \return Angle between the two vectors
107  \f[ \Delta R = \sqrt{ ( \Delta \phi )^2 + ( \Delta \eta )^2 } \f]
108  */
109  template <class Vector1, class Vector2>
110  inline typename Vector1::Scalar DeltaR( const Vector1 & v1, const Vector2 & v2) {
111  return std::sqrt( DeltaR2(v1,v2) );
112  }
113 
114  /**
115  Find difference in Rapidity (y) and Phi betwen two generic vectors
116  The only requirements on the Vector classes is that they implement the Phi() and Rapidity() method
117  \param v1 Vector 1
118  \param v2 Vector 2
119  \return Angle between the two vectors
120  \f[ \Delta R = \sqrt{ ( \Delta \phi )^2 + ( \Delta y )^2 } \f],
121  */
122  template <class Vector1, class Vector2>
123  inline typename Vector1::Scalar DeltaRapidityPhi( const Vector1 & v1, const Vector2 & v2) {
124  return std::sqrt( DeltaR2RapidityPhi(v1,v2) );
125  }
126 
127  /**
128  Find CosTheta Angle between two generic 3D vectors
129  pre-requisite: vectors implement the X(), Y() and Z()
130  \param v1 Vector v1
131  \param v2 Vector v2
132  \return cosine of Angle between the two vectors
133  \f[ \cos \theta = \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
134  */
135  // this cannot be made all generic since Mag2() for 2, 3 or 4 D is different
136  // need to have a specialization for polar Coordinates ??
137  template <class Vector1, class Vector2>
138  double CosTheta( const Vector1 & v1, const Vector2 & v2) {
139  double arg;
140  double v1_r2 = v1.X()*v1.X() + v1.Y()*v1.Y() + v1.Z()*v1.Z();
141  double v2_r2 = v2.X()*v2.X() + v2.Y()*v2.Y() + v2.Z()*v2.Z();
142  double ptot2 = v1_r2*v2_r2;
143  if(ptot2 <= 0) {
144  arg = 0.0;
145  }else{
146  double pdot = v1.X()*v2.X() + v1.Y()*v2.Y() + v1.Z()*v2.Z();
147  arg = pdot/std::sqrt(ptot2);
148  if(arg > 1.0) arg = 1.0;
149  if(arg < -1.0) arg = -1.0;
150  }
151  return arg;
152  }
153 
154 
155  /**
156  Find Angle between two vectors.
157  Use the CosTheta() function
158  \param v1 Vector v1
159  \param v2 Vector v2
160  \return Angle between the two vectors
161  \f[ \theta = \cos ^{-1} \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
162  */
163  template <class Vector1, class Vector2>
164  inline double Angle( const Vector1 & v1, const Vector2 & v2) {
165  return std::acos( CosTheta(v1, v2) );
166  }
167 
168  /**
169  Find the projection of v along the given direction u.
170  \param v Vector v for which the propjection is to be found
171  \param u Vector specifying the direction
172  \return Vector projection (same type of v)
173  \f[ \vec{proj} = \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
174  Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
175  */
176  template <class Vector1, class Vector2>
177  Vector1 ProjVector( const Vector1 & v, const Vector2 & u) {
178  double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
179  if (magU2 == 0) return Vector1(0,0,0);
180  double d = v.Dot(u)/magU2;
181  return Vector1( u.X() * d, u.Y() * d, u.Z() * d);
182  }
183 
184  /**
185  Find the vector component of v perpendicular to the given direction of u
186  \param v Vector v for which the perpendicular component is to be found
187  \param u Vector specifying the direction
188  \return Vector component of v which is perpendicular to u
189  \f[ \vec{perp} = \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
190  Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
191  */
192  template <class Vector1, class Vector2>
193  inline Vector1 PerpVector( const Vector1 & v, const Vector2 & u) {
194  return v - ProjVector(v,u);
195  }
196 
197  /**
198  Find the magnitude square of the vector component of v perpendicular to the given direction of u
199  \param v Vector v for which the perpendicular component is to be found
200  \param u Vector specifying the direction
201  \return square value of the component of v which is perpendicular to u
202  \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} |^2 \f]
203  Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
204  */
205  template <class Vector1, class Vector2>
206  inline double Perp2( const Vector1 & v, const Vector2 & u) {
207  double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
208  double prjvu = v.Dot(u);
209  double magV2 = v.Dot(v);
210  return magU2 > 0.0 ? magV2-prjvu*prjvu/magU2 : magV2;
211  }
212 
213  /**
214  Find the magnitude of the vector component of v perpendicular to the given direction of u
215  \param v Vector v for which the perpendicular component is to be found
216  \param u Vector specifying the direction
217  \return value of the component of v which is perpendicular to u
218  \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} | \f]
219  Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
220  */
221  template <class Vector1, class Vector2>
222  inline double Perp( const Vector1 & v, const Vector2 & u) {
223  return std::sqrt(Perp2(v,u) );
224  }
225 
226 
227 
228  // Lorentz Vector functions
229 
230 
231  /**
232  return the invariant mass of two LorentzVector
233  The only requirement on the LorentzVector is that they need to implement the
234  X() , Y(), Z() and E() methods.
235  \param v1 LorenzVector 1
236  \param v2 LorenzVector 2
237  \return invariant mass M
238  \f[ M_{12} = \sqrt{ (\vec{v1} + \vec{v2} ) \cdot (\vec{v1} + \vec{v2} ) } \f]
239  */
240  template <class Vector1, class Vector2>
241  inline typename Vector1::Scalar InvariantMass( const Vector1 & v1, const Vector2 & v2) {
242  typedef typename Vector1::Scalar Scalar;
243  Scalar ee = (v1.E() + v2.E() );
244  Scalar xx = (v1.X() + v2.X() );
245  Scalar yy = (v1.Y() + v2.Y() );
246  Scalar zz = (v1.Z() + v2.Z() );
247  Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
248  return mm2 < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
249  // PxPyPzE4D<double> q(xx,yy,zz,ee);
250  // return q.M();
251  //return ( v1 + v2).mag();
252  }
253 
254  template <class Vector1, class Vector2>
255  inline typename Vector1::Scalar InvariantMass2( const Vector1 & v1, const Vector2 & v2) {
256  typedef typename Vector1::Scalar Scalar;
257  Scalar ee = (v1.E() + v2.E() );
258  Scalar xx = (v1.X() + v2.X() );
259  Scalar yy = (v1.Y() + v2.Y() );
260  Scalar zz = (v1.Z() + v2.Z() );
261  Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
262  return mm2 ; // < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
263  // PxPyPzE4D<double> q(xx,yy,zz,ee);
264  // return q.M();
265  //return ( v1 + v2).mag();
266  }
267 
268  // rotation and transformations
269 
270 
271 #ifndef __CINT__
272  /**
273  rotation along X axis for a generic vector by an Angle alpha
274  returning a new vector.
275  The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
276  and SetXYZ methods.
277  */
278  template <class Vector>
279  Vector RotateX(const Vector & v, double alpha) {
280  double sina = sin(alpha);
281  double cosa = cos(alpha);
282  double y2 = v.Y() * cosa - v.Z()*sina;
283  double z2 = v.Z() * cosa + v.Y() * sina;
284  Vector vrot;
285  vrot.SetXYZ(v.X(), y2, z2);
286  return vrot;
287  }
288 
289  /**
290  rotation along Y axis for a generic vector by an Angle alpha
291  returning a new vector.
292  The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
293  and SetXYZ methods.
294  */
295  template <class Vector>
296  Vector RotateY(const Vector & v, double alpha) {
297  double sina = sin(alpha);
298  double cosa = cos(alpha);
299  double x2 = v.X() * cosa + v.Z() * sina;
300  double z2 = v.Z() * cosa - v.X() * sina;
301  Vector vrot;
302  vrot.SetXYZ(x2, v.Y(), z2);
303  return vrot;
304  }
305 
306  /**
307  rotation along Z axis for a generic vector by an Angle alpha
308  returning a new vector.
309  The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
310  and SetXYZ methods.
311  */
312  template <class Vector>
313  Vector RotateZ(const Vector & v, double alpha) {
314  double sina = sin(alpha);
315  double cosa = cos(alpha);
316  double x2 = v.X() * cosa - v.Y() * sina;
317  double y2 = v.Y() * cosa + v.X() * sina;
318  Vector vrot;
319  vrot.SetXYZ(x2, y2, v.Z());
320  return vrot;
321  }
322 
323 
324  /**
325  rotation on a generic vector using a generic rotation class.
326  The only requirement on the vector is that implements the
327  X(), Y(), Z() and SetXYZ methods.
328  The requirement on the rotation matrix is that need to implement the
329  (i,j) operator returning the matrix element with R(0,0) = xx element
330  */
331  template<class Vector, class RotationMatrix>
332  Vector Rotate(const Vector &v, const RotationMatrix & rot) {
333  double xX = v.X();
334  double yY = v.Y();
335  double zZ = v.Z();
336  double x2 = rot(0,0)*xX + rot(0,1)*yY + rot(0,2)*zZ;
337  double y2 = rot(1,0)*xX + rot(1,1)*yY + rot(1,2)*zZ;
338  double z2 = rot(2,0)*xX + rot(2,1)*yY + rot(2,2)*zZ;
339  Vector vrot;
340  vrot.SetXYZ(x2,y2,z2);
341  return vrot;
342  }
343 
344  /**
345  Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost
346  The only requirement on the vector is that implements the
347  X(), Y(), Z(), T() and SetXYZT methods.
348  The requirement on the boost vector is that needs to implement the
349  X(), Y() , Z() retorning the vector elements describing the boost
350  The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
351  */
352  template <class LVector, class BoostVector>
353  LVector boost(const LVector & v, const BoostVector & b) {
354  double bx = b.X();
355  double by = b.Y();
356  double bz = b.Z();
357  double b2 = bx*bx + by*by + bz*bz;
358  if (b2 >= 1) {
359  GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
360  return LVector();
361  }
362  double gamma = 1.0 / std::sqrt(1.0 - b2);
363  double bp = bx*v.X() + by*v.Y() + bz*v.Z();
364  double gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
365  double x2 = v.X() + gamma2*bp*bx + gamma*bx*v.T();
366  double y2 = v.Y() + gamma2*bp*by + gamma*by*v.T();
367  double z2 = v.Z() + gamma2*bp*bz + gamma*bz*v.T();
368  double t2 = gamma*(v.T() + bp);
369  LVector lv;
370  lv.SetXYZT(x2,y2,z2,t2);
371  return lv;
372  }
373 
374 
375  /**
376  Boost a generic Lorentz Vector class along the X direction with a factor beta
377  The only requirement on the vector is that implements the
378  X(), Y(), Z(), T() and SetXYZT methods.
379  The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
380  */
381  template <class LVector, class T>
382  LVector boostX(const LVector & v, T beta) {
383  if (beta >= 1) {
384  GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
385  return LVector();
386  }
387  T gamma = 1.0/ std::sqrt(1.0 - beta*beta);
388  typename LVector::Scalar x2 = gamma * v.X() + gamma * beta * v.T();
389  typename LVector::Scalar t2 = gamma * beta * v.X() + gamma * v.T();
390 
391  LVector lv;
392  lv.SetXYZT(x2,v.Y(),v.Z(),t2);
393  return lv;
394  }
395 
396  /**
397  Boost a generic Lorentz Vector class along the Y direction with a factor beta
398  The only requirement on the vector is that implements the
399  X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
400  The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
401  */
402  template <class LVector>
403  LVector boostY(const LVector & v, double beta) {
404  if (beta >= 1) {
405  GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
406  return LVector();
407  }
408  double gamma = 1.0/ std::sqrt(1.0 - beta*beta);
409  double y2 = gamma * v.Y() + gamma * beta * v.T();
410  double t2 = gamma * beta * v.Y() + gamma * v.T();
411  LVector lv;
412  lv.SetXYZT(v.X(),y2,v.Z(),t2);
413  return lv;
414  }
415 
416  /**
417  Boost a generic Lorentz Vector class along the Z direction with a factor beta
418  The only requirement on the vector is that implements the
419  X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
420  The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
421  */
422  template <class LVector>
423  LVector boostZ(const LVector & v, double beta) {
424  if (beta >= 1) {
425  GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
426  return LVector();
427  }
428  double gamma = 1.0/ std::sqrt(1.0 - beta*beta);
429  double z2 = gamma * v.Z() + gamma * beta * v.T();
430  double t2 = gamma * beta * v.Z() + gamma * v.T();
431  LVector lv;
432  lv.SetXYZT(v.X(),v.Y(),z2,t2);
433  return lv;
434  }
435 
436 #endif
437 
438 
439 
440 
441  // MATRIX VECTOR MULTIPLICATION
442  // cannot define an operator * otherwise conflicts with rotations
443  // operations like Rotation3D * vector use Mult
444 
445  /**
446  Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
447  Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no check on the matrix size !!
448  */
449  template<class Matrix, class CoordSystem, class U>
450  inline
453  vret.SetXYZ( m(0,0) * v.x() + m(0,1) * v.y() + m(0,2) * v.z() ,
454  m(1,0) * v.x() + m(1,1) * v.y() + m(1,2) * v.z() ,
455  m(2,0) * v.x() + m(2,1) * v.y() + m(2,2) * v.z() );
456  return vret;
457  }
458 
459 
460  /**
461  Multiplications of a generic matrices with a generic PositionVector
462  Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no check on the matrix size !!
463  */
464  template<class Matrix, class CoordSystem, class U>
465  inline
468  pret.SetXYZ( m(0,0) * p.x() + m(0,1) * p.y() + m(0,2) * p.z() ,
469  m(1,0) * p.x() + m(1,1) * p.y() + m(1,2) * p.z() ,
470  m(2,0) * p.x() + m(2,1) * p.y() + m(2,2) * p.z() );
471  return pret;
472  }
473 
474 
475  /**
476  Multiplications of a generic matrices with a LorentzVector described
477  in any coordinate system.
478  Assume that the matrix implements the operator( i,j) and that it has at least 4 columns and 4 rows. There is no check on the matrix size !!
479  */
480  // this will not be ambigous with operator*(Scalar, LorentzVector) since that one // Scalar is passed by value
481  template<class CoordSystem, class Matrix>
482  inline
485  vret.SetXYZT( m(0,0)*v.x() + m(0,1)*v.y() + m(0,2)*v.z() + m(0,3)* v.t() ,
486  m(1,0)*v.x() + m(1,1)*v.y() + m(1,2)*v.z() + m(1,3)* v.t() ,
487  m(2,0)*v.x() + m(2,1)*v.y() + m(2,2)*v.z() + m(2,3)* v.t() ,
488  m(3,0)*v.x() + m(3,1)*v.y() + m(3,2)*v.z() + m(3,3)* v.t() );
489  return vret;
490  }
491 
492 
493 
494  // non-template utility functions for all objects
495 
496 
497  /**
498  Return a phi angle in the interval (0,2*PI]
499  */
500  double Phi_0_2pi(double phi);
501  /**
502  Returns phi angle in the interval (-PI,PI]
503  */
504  double Phi_mpi_pi(double phi);
505 
506 
507 
508  } // end namespace Vector Util
509 
510 
511 
512  } // end namespace Math
513 
514 } // end namespace ROOT
515 
516 
517 #endif /* ROOT_Math_GenVector_VectorUtil */
ROOT::Math::Cephes::gamma
double gamma(double x)
Definition: SpecFuncCephes.cxx:339
m
auto * m
Definition: textangle.C:8
ROOT::Math::VectorUtil::DeltaR2RapidityPhi
Vector1::Scalar DeltaR2RapidityPhi(const Vector1 &v1, const Vector2 &v2)
Find square of the difference in true rapidity (y) and Phi betwen two generic vectors The only requir...
Definition: VectorUtil.h:95
ROOT::Math::VectorUtil::RotateY
Vector RotateY(const Vector &v, double alpha)
rotation along Y axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:296
ROOT::Math::LorentzVector::SetXYZT
LorentzVector< CoordSystem > & SetXYZT(Scalar xx, Scalar yy, Scalar zz, Scalar tt)
set the values of the vector from the cartesian components (x,y,z,t) (if the vector is held in anothe...
Definition: LorentzVector.h:235
ROOT::Math::VectorUtil::PerpVector
Vector1 PerpVector(const Vector1 &v, const Vector2 &u)
Find the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:193
ROOT::Math::PositionVector3D
Class describing a generic position vector (point) in 3 dimensions.
Definition: PositionVector3D.h:53
ROOT::Math::VectorUtil::InvariantMass2
Vector1::Scalar InvariantMass2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:255
ROOT::Math::VectorUtil::boostZ
LVector boostZ(const LVector &v, double beta)
Boost a generic Lorentz Vector class along the Z direction with a factor beta The only requirement on...
Definition: VectorUtil.h:423
ROOT::Math::VectorUtil::Perp2
double Perp2(const Vector1 &v, const Vector2 &u)
Find the magnitude square of the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:206
ROOT::Math::VectorUtil::boostX
LVector boostX(const LVector &v, T beta)
Boost a generic Lorentz Vector class along the X direction with a factor beta The only requirement on...
Definition: VectorUtil.h:382
sin
double sin(double)
cos
double cos(double)
ROOT::Math::VectorUtil::DeltaR2
Vector1::Scalar DeltaR2(const Vector1 &v1, const Vector2 &v2)
Find square of the difference in pseudorapidity (Eta) and Phi betwen two generic vectors The only req...
Definition: VectorUtil.h:80
ROOT::Math::VectorUtil::Angle
double Angle(const Vector1 &v1, const Vector2 &v2)
Find Angle between two vectors.
Definition: VectorUtil.h:164
ROOT::Math::VectorUtil::DeltaR
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:110
ROOT::Math::VectorUtil::RotateX
Vector RotateX(const Vector &v, double alpha)
rotation along X axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:279
ROOT::Math::beta
double beta(double x, double y)
Calculates the beta function.
Definition: SpecFuncMathCore.cxx:111
v
@ v
Definition: rootcling_impl.cxx:3643
b
#define b(i)
Definition: RSha256.hxx:100
lv
auto * lv
Definition: textalign.C:5
ROOT::Math::VectorUtil::boostY
LVector boostY(const LVector &v, double beta)
Boost a generic Lorentz Vector class along the Y direction with a factor beta The only requirement on...
Definition: VectorUtil.h:403
ROOT::Math::GenVector::Throw
void Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
Definition: GenVector_exception.h:80
Boost.h
Math.h
ROOT::Math::VectorUtil::DeltaPhi
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:59
ROOT::Math::PositionVector3D::x
Scalar x() const
Definition: PositionVector3D.h:445
ROOT::Math::PositionVector3D::z
Scalar z() const
Definition: PositionVector3D.h:447
ROOT::Math::VectorUtil::CosTheta
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:138
TGeant4Unit::mm2
static constexpr double mm2
Definition: TGeant4SystemOfUnits.h:109
ROOT::Math::VectorUtil::Phi_0_2pi
double Phi_0_2pi(double phi)
Return a phi angle in the interval (0,2*PI].
Definition: VectorUtil.cxx:22
sqrt
double sqrt(double)
ROOT::Math::VectorUtil::DeltaRapidityPhi
Vector1::Scalar DeltaRapidityPhi(const Vector1 &v1, const Vector2 &v2)
Find difference in Rapidity (y) and Phi betwen two generic vectors The only requirements on the Vecto...
Definition: VectorUtil.h:123
acos
double acos(double)
ROOT::Math::VectorUtil::Mult
DisplacementVector3D< CoordSystem, U > Mult(const Matrix &m, const DisplacementVector3D< CoordSystem, U > &v)
Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
Definition: VectorUtil.h:451
ROOT::Math::VectorUtil::boost
LVector boost(const LVector &v, const BoostVector &b)
Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost The only re...
Definition: VectorUtil.h:353
v1
@ v1
Definition: rootcling_impl.cxx:3645
ROOT::Math::Scalar
Rotation3D::Scalar Scalar
Definition: Rotation3DxAxial.cxx:69
v2
@ v2
Definition: rootcling_impl.cxx:3646
M_PI
#define M_PI
Definition: Rotated.cxx:105
ROOT::Math::VectorUtil::ProjVector
Vector1 ProjVector(const Vector1 &v, const Vector2 &u)
Find the projection of v along the given direction u.
Definition: VectorUtil.h:177
ROOT::Math::DisplacementVector3D::SetXYZ
DisplacementVector3D< CoordSystem, Tag > & SetXYZ(Scalar a, Scalar b, Scalar c)
set the values of the vector from the cartesian components (x,y,z) (if the vector is held in polar or...
Definition: DisplacementVector3D.h:249
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:34
ROOT::Math::DisplacementVector3D
Class describing a generic displacement vector in 3 dimensions.
Definition: DisplacementVector3D.h:56
d
#define d(i)
Definition: RSha256.hxx:102
ROOT::Math::VectorUtil::Phi_mpi_pi
double Phi_mpi_pi(double phi)
Returns phi angle in the interval (-PI,PI].
Definition: VectorUtil.cxx:36
x2
static const double x2[5]
Definition: RooGaussKronrodIntegrator1D.cxx:364
ROOT::Math::PositionVector3D::y
Scalar y() const
Definition: PositionVector3D.h:446
ROOT::Math::VectorUtil::Rotate
Vector Rotate(const Vector &v, const RotationMatrix &rot)
rotation on a generic vector using a generic rotation class.
Definition: VectorUtil.h:332
ROOT::Math::VectorUtil::RotateZ
Vector RotateZ(const Vector &v, double alpha)
rotation along Z axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:313
ROOT::Math::VectorUtil::Perp
double Perp(const Vector1 &v, const Vector2 &u)
Find the magnitude of the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:222
ROOT
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Definition: EExecutionPolicy.hxx:4
ROOT::Math::VectorUtil::InvariantMass
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:241
Math
Namespace for new Math classes and functions.
ROOT::Math::LorentzVector
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
Definition: LorentzVector.h:59