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