Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
25
26namespace 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 using std::sqrt;
112 return sqrt( DeltaR2(v1,v2) );
113 }
114
115 /**
116 Find difference in Rapidity (y) and Phi betwen two generic vectors
117 The only requirements on the Vector classes is that they implement the Phi() and Rapidity() method
118 \param v1 Vector 1
119 \param v2 Vector 2
120 \return Angle between the two vectors
121 \f[ \Delta R = \sqrt{ ( \Delta \phi )^2 + ( \Delta y )^2 } \f],
122 */
123 template <class Vector1, class Vector2>
124 inline typename Vector1::Scalar DeltaRapidityPhi( const Vector1 & v1, const Vector2 & v2) {
125 using std::sqrt;
126 return sqrt( DeltaR2RapidityPhi(v1,v2) );
127 }
128
129 /**
130 Find CosTheta Angle between two generic 3D vectors
131 pre-requisite: vectors implement the X(), Y() and Z()
132 \param v1 Vector v1
133 \param v2 Vector v2
134 \return cosine of Angle between the two vectors
135 \f[ \cos \theta = \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
136 */
137 // this cannot be made all generic since Mag2() for 2, 3 or 4 D is different
138 // need to have a specialization for polar Coordinates ??
139 template <class Vector1, class Vector2>
140 double CosTheta( const Vector1 & v1, const Vector2 & v2) {
141 double arg;
142 double v1_r2 = v1.X()*v1.X() + v1.Y()*v1.Y() + v1.Z()*v1.Z();
143 double v2_r2 = v2.X()*v2.X() + v2.Y()*v2.Y() + v2.Z()*v2.Z();
144 double ptot2 = v1_r2*v2_r2;
145 if(ptot2 <= 0) {
146 arg = 0.0;
147 }else{
148 double pdot = v1.X()*v2.X() + v1.Y()*v2.Y() + v1.Z()*v2.Z();
149 using std::sqrt;
150 arg = pdot/sqrt(ptot2);
151 if(arg > 1.0) arg = 1.0;
152 if(arg < -1.0) arg = -1.0;
153 }
154 return arg;
155 }
156
157
158 /**
159 Find Angle between two vectors.
160 Use the CosTheta() function
161 \param v1 Vector v1
162 \param v2 Vector v2
163 \return Angle between the two vectors
164 \f[ \theta = \cos ^{-1} \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
165 */
166 template <class Vector1, class Vector2>
167 inline double Angle( const Vector1 & v1, const Vector2 & v2) {
168 using std::acos;
169 return acos( CosTheta(v1, v2) );
170 }
171
172 /**
173 Find the projection of v along the given direction u.
174 \param v Vector v for which the propjection is to be found
175 \param u Vector specifying the direction
176 \return Vector projection (same type of v)
177 \f[ \vec{proj} = \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
178 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
179 */
180 template <class Vector1, class Vector2>
181 Vector1 ProjVector( const Vector1 & v, const Vector2 & u) {
182 double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
183 if (magU2 == 0) return Vector1(0,0,0);
184 double d = v.Dot(u)/magU2;
185 return Vector1( u.X() * d, u.Y() * d, u.Z() * d);
186 }
187
188 /**
189 Find the vector component of v perpendicular to the given direction of u
190 \param v Vector v for which the perpendicular component is to be found
191 \param u Vector specifying the direction
192 \return Vector component of v which is perpendicular to u
193 \f[ \vec{perp} = \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
194 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
195 */
196 template <class Vector1, class Vector2>
197 inline Vector1 PerpVector( const Vector1 & v, const Vector2 & u) {
198 return v - ProjVector(v,u);
199 }
200
201 /**
202 Find the magnitude square of the vector component of v perpendicular to the given direction of u
203 \param v Vector v for which the perpendicular component is to be found
204 \param u Vector specifying the direction
205 \return square value of the component of v which is perpendicular to u
206 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} |^2 \f]
207 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
208 */
209 template <class Vector1, class Vector2>
210 inline double Perp2( const Vector1 & v, const Vector2 & u) {
211 double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
212 double prjvu = v.Dot(u);
213 double magV2 = v.Dot(v);
214 return magU2 > 0.0 ? magV2-prjvu*prjvu/magU2 : magV2;
215 }
216
217 /**
218 Find the magnitude of the vector component of v perpendicular to the given direction of u
219 \param v Vector v for which the perpendicular component is to be found
220 \param u Vector specifying the direction
221 \return value of the component of v which is perpendicular to u
222 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} | \f]
223 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
224 */
225 template <class Vector1, class Vector2>
226 inline double Perp( const Vector1 & v, const Vector2 & u) {
227 using std::sqrt;
228 return sqrt(Perp2(v,u) );
229 }
230
231
232
233 // Lorentz Vector functions
234
235
236 /**
237 return the invariant mass of two LorentzVector
238 The only requirement on the LorentzVector is that they need to implement the
239 X() , Y(), Z() and E() methods.
240 \param v1 LorenzVector 1
241 \param v2 LorenzVector 2
242 \return invariant mass M
243 \f[ M_{12} = \sqrt{ (\vec{v1} + \vec{v2} ) \cdot (\vec{v1} + \vec{v2} ) } \f]
244 */
245 template <class Vector1, class Vector2>
246 inline typename Vector1::Scalar InvariantMass( const Vector1 & v1, const Vector2 & v2) {
247 typedef typename Vector1::Scalar Scalar;
248 Scalar ee = (v1.E() + v2.E() );
249 Scalar xx = (v1.X() + v2.X() );
250 Scalar yy = (v1.Y() + v2.Y() );
251 Scalar zz = (v1.Z() + v2.Z() );
252 Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
253 using std::sqrt;
254 return mm2 < 0.0 ? -sqrt(-mm2) : sqrt(mm2);
255 // PxPyPzE4D<double> q(xx,yy,zz,ee);
256 // return q.M();
257 //return ( v1 + v2).mag();
258 }
259
260 template <class Vector1, class Vector2>
261 inline typename Vector1::Scalar InvariantMass2( const Vector1 & v1, const Vector2 & v2) {
262 typedef typename Vector1::Scalar Scalar;
263 Scalar ee = (v1.E() + v2.E() );
264 Scalar xx = (v1.X() + v2.X() );
265 Scalar yy = (v1.Y() + v2.Y() );
266 Scalar zz = (v1.Z() + v2.Z() );
267 Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
268 return mm2 ; // < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
269 // PxPyPzE4D<double> q(xx,yy,zz,ee);
270 // return q.M();
271 //return ( v1 + v2).mag();
272 }
273
274 // rotation and transformations
275
276
277#ifndef __CINT__
278 /**
279 rotation along X axis for a generic vector by an Angle alpha
280 returning a new vector.
281 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
282 and SetXYZ methods.
283 */
284 template <class Vector>
285 Vector RotateX(const Vector & v, double alpha) {
286 using std::sin;
287 double sina = sin(alpha);
288 using std::cos;
289 double cosa = cos(alpha);
290 double y2 = v.Y() * cosa - v.Z()*sina;
291 double z2 = v.Z() * cosa + v.Y() * sina;
292 Vector vrot;
293 vrot.SetXYZ(v.X(), y2, z2);
294 return vrot;
295 }
296
297 /**
298 rotation along Y axis for a generic vector by an Angle alpha
299 returning a new vector.
300 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
301 and SetXYZ methods.
302 */
303 template <class Vector>
304 Vector RotateY(const Vector & v, double alpha) {
305 using std::sin;
306 double sina = sin(alpha);
307 using std::cos;
308 double cosa = cos(alpha);
309 double x2 = v.X() * cosa + v.Z() * sina;
310 double z2 = v.Z() * cosa - v.X() * sina;
311 Vector vrot;
312 vrot.SetXYZ(x2, v.Y(), z2);
313 return vrot;
314 }
315
316 /**
317 rotation along Z axis for a generic vector by an Angle alpha
318 returning a new vector.
319 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
320 and SetXYZ methods.
321 */
322 template <class Vector>
323 Vector RotateZ(const Vector & v, double alpha) {
324 using std::sin;
325 double sina = sin(alpha);
326 using std::cos;
327 double cosa = cos(alpha);
328 double x2 = v.X() * cosa - v.Y() * sina;
329 double y2 = v.Y() * cosa + v.X() * sina;
330 Vector vrot;
331 vrot.SetXYZ(x2, y2, v.Z());
332 return vrot;
333 }
334
335
336 /**
337 rotation on a generic vector using a generic rotation class.
338 The only requirement on the vector is that implements the
339 X(), Y(), Z() and SetXYZ methods.
340 The requirement on the rotation matrix is that need to implement the
341 (i,j) operator returning the matrix element with R(0,0) = xx element
342 */
343 template<class Vector, class RotationMatrix>
344 Vector Rotate(const Vector &v, const RotationMatrix & rot) {
345 double xX = v.X();
346 double yY = v.Y();
347 double zZ = v.Z();
348 double x2 = rot(0,0)*xX + rot(0,1)*yY + rot(0,2)*zZ;
349 double y2 = rot(1,0)*xX + rot(1,1)*yY + rot(1,2)*zZ;
350 double z2 = rot(2,0)*xX + rot(2,1)*yY + rot(2,2)*zZ;
351 Vector vrot;
352 vrot.SetXYZ(x2,y2,z2);
353 return vrot;
354 }
355
356 /**
357 Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost
358 The only requirement on the vector is that implements the
359 X(), Y(), Z(), T() and SetXYZT methods.
360 The requirement on the boost vector is that needs to implement the
361 X(), Y() , Z() retorning the vector elements describing the boost
362 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
363 */
364 template <class LVector, class BoostVector>
365 LVector boost(const LVector & v, const BoostVector & b) {
366 double bx = b.X();
367 double by = b.Y();
368 double bz = b.Z();
369 double b2 = bx*bx + by*by + bz*bz;
370 if (b2 >= 1) {
371 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
372 return LVector();
373 }
374 using std::sqrt;
375 double gamma = 1.0 / sqrt(1.0 - b2);
376 double bp = bx*v.X() + by*v.Y() + bz*v.Z();
377 double gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
378 double x2 = v.X() + gamma2*bp*bx + gamma*bx*v.T();
379 double y2 = v.Y() + gamma2*bp*by + gamma*by*v.T();
380 double z2 = v.Z() + gamma2*bp*bz + gamma*bz*v.T();
381 double t2 = gamma*(v.T() + bp);
382 LVector lv;
383 lv.SetXYZT(x2,y2,z2,t2);
384 return lv;
385 }
386
387
388 /**
389 Boost a generic Lorentz Vector class along the X direction with a factor beta
390 The only requirement on the vector is that implements the
391 X(), Y(), Z(), T() and SetXYZT methods.
392 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
393 */
394 template <class LVector, class T>
395 LVector boostX(const LVector & v, T beta) {
396 if (beta >= 1) {
397 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
398 return LVector();
399 }
400 using std::sqrt;
401 T gamma = 1.0/ sqrt(1.0 - beta*beta);
402 typename LVector::Scalar x2 = gamma * v.X() + gamma * beta * v.T();
403 typename LVector::Scalar t2 = gamma * beta * v.X() + gamma * v.T();
404
405 LVector lv;
406 lv.SetXYZT(x2,v.Y(),v.Z(),t2);
407 return lv;
408 }
409
410 /**
411 Boost a generic Lorentz Vector class along the Y direction with a factor beta
412 The only requirement on the vector is that implements the
413 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
414 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
415 */
416 template <class LVector>
417 LVector boostY(const LVector & v, double beta) {
418 if (beta >= 1) {
419 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
420 return LVector();
421 }
422 using std::sqrt;
423 double gamma = 1.0/ sqrt(1.0 - beta*beta);
424 double y2 = gamma * v.Y() + gamma * beta * v.T();
425 double t2 = gamma * beta * v.Y() + gamma * v.T();
426 LVector lv;
427 lv.SetXYZT(v.X(),y2,v.Z(),t2);
428 return lv;
429 }
430
431 /**
432 Boost a generic Lorentz Vector class along the Z direction with a factor beta
433 The only requirement on the vector is that implements the
434 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
435 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
436 */
437 template <class LVector>
438 LVector boostZ(const LVector & v, double beta) {
439 if (beta >= 1) {
440 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
441 return LVector();
442 }
443 using std::sqrt;
444 double gamma = 1.0/ sqrt(1.0 - beta*beta);
445 double z2 = gamma * v.Z() + gamma * beta * v.T();
446 double t2 = gamma * beta * v.Z() + gamma * v.T();
447 LVector lv;
448 lv.SetXYZT(v.X(),v.Y(),z2,t2);
449 return lv;
450 }
451
452#endif
453
454
455
456
457 // MATRIX VECTOR MULTIPLICATION
458 // cannot define an operator * otherwise conflicts with rotations
459 // operations like Rotation3D * vector use Mult
460
461 /**
462 Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
463 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 !!
464 */
465 template<class Matrix, class CoordSystem, class U>
466 inline
469 vret.SetXYZ( m(0,0) * v.x() + m(0,1) * v.y() + m(0,2) * v.z() ,
470 m(1,0) * v.x() + m(1,1) * v.y() + m(1,2) * v.z() ,
471 m(2,0) * v.x() + m(2,1) * v.y() + m(2,2) * v.z() );
472 return vret;
473 }
474
475
476 /**
477 Multiplications of a generic matrices with a generic PositionVector
478 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 !!
479 */
480 template<class Matrix, class CoordSystem, class U>
481 inline
484 pret.SetXYZ( m(0,0) * p.x() + m(0,1) * p.y() + m(0,2) * p.z() ,
485 m(1,0) * p.x() + m(1,1) * p.y() + m(1,2) * p.z() ,
486 m(2,0) * p.x() + m(2,1) * p.y() + m(2,2) * p.z() );
487 return pret;
488 }
489
490
491 /**
492 Multiplications of a generic matrices with a LorentzVector described
493 in any coordinate system.
494 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 !!
495 */
496 // this will not be ambigous with operator*(Scalar, LorentzVector) since that one // Scalar is passed by value
497 template<class CoordSystem, class Matrix>
498 inline
501 vret.SetXYZT( m(0,0)*v.x() + m(0,1)*v.y() + m(0,2)*v.z() + m(0,3)* v.t() ,
502 m(1,0)*v.x() + m(1,1)*v.y() + m(1,2)*v.z() + m(1,3)* v.t() ,
503 m(2,0)*v.x() + m(2,1)*v.y() + m(2,2)*v.z() + m(2,3)* v.t() ,
504 m(3,0)*v.x() + m(3,1)*v.y() + m(3,2)*v.z() + m(3,3)* v.t() );
505 return vret;
506 }
507
508
509
510 // non-template utility functions for all objects
511
512
513 /**
514 Return a phi angle in the interval (0,2*PI]
515 */
516 double Phi_0_2pi(double phi);
517 /**
518 Returns phi angle in the interval (-PI,PI]
519 */
520 double Phi_mpi_pi(double phi);
521
522
523
524 } // end namespace Vector Util
525
526
527
528 } // end namespace Math
529
530} // end namespace ROOT
531
532
533#endif /* ROOT_Math_GenVector_VectorUtil */
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
static const double x2[5]
#define M_PI
Definition Rotated.cxx:105
double acos(double)
double cos(double)
double sqrt(double)
double sin(double)
Class describing a generic displacement vector in 3 dimensions.
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...
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
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...
Class describing a generic position vector (point) in 3 dimensions.
double beta(double x, double y)
Calculates the beta function.
Namespace for new Math classes and functions.
void Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
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:467
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
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:226
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:417
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:124
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
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
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:365
double Angle(const Vector1 &v1, const Vector2 &v2)
Find Angle between two vectors.
Definition VectorUtil.h:167
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:323
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
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:304
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:285
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:140
Vector1 PerpVector(const Vector1 &v, const Vector2 &u)
Find the vector component of v perpendicular to the given direction of u.
Definition VectorUtil.h:197
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:210
double Phi_0_2pi(double phi)
Return a phi angle in the interval (0,2*PI].
Vector Rotate(const Vector &v, const RotationMatrix &rot)
rotation on a generic vector using a generic rotation class.
Definition VectorUtil.h:344
double Phi_mpi_pi(double phi)
Returns phi angle in the interval (-PI,PI].
Vector1::Scalar InvariantMass2(const Vector1 &v1, const Vector2 &v2)
Definition VectorUtil.h:261
Vector1 ProjVector(const Vector1 &v, const Vector2 &u)
Find the projection of v along the given direction u.
Definition VectorUtil.h:181
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:438
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:395
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:246
Rotation3D::Scalar Scalar
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
auto * lv
Definition textalign.C:5
auto * m
Definition textangle.C:8