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