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 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 square of the difference in true rapidity (y) and Phi betwen 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 betwen 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 betwen 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#ifndef __CINT__
280 /**
281 rotation along X axis for a generic vector by an Angle alpha
282 returning a new vector.
283 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
284 and SetXYZ methods.
285 */
286 template <class Vector>
287 Vector RotateX(const Vector & v, double alpha) {
288 using std::sin;
289 double sina = sin(alpha);
290 using std::cos;
291 double cosa = cos(alpha);
292 double y2 = v.Y() * cosa - v.Z()*sina;
293 double z2 = v.Z() * cosa + v.Y() * sina;
294 Vector vrot;
295 vrot.SetXYZ(v.X(), y2, z2);
296 return vrot;
297 }
298
299 /**
300 rotation along Y axis for a generic vector by an Angle alpha
301 returning a new vector.
302 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
303 and SetXYZ methods.
304 */
305 template <class Vector>
306 Vector RotateY(const Vector & v, double alpha) {
307 using std::sin;
308 double sina = sin(alpha);
309 using std::cos;
310 double cosa = cos(alpha);
311 double x2 = v.X() * cosa + v.Z() * sina;
312 double z2 = v.Z() * cosa - v.X() * sina;
313 Vector vrot;
314 vrot.SetXYZ(x2, v.Y(), z2);
315 return vrot;
316 }
317
318 /**
319 rotation along Z axis for a generic vector by an Angle alpha
320 returning a new vector.
321 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
322 and SetXYZ methods.
323 */
324 template <class Vector>
325 Vector RotateZ(const Vector & v, double alpha) {
326 using std::sin;
327 double sina = sin(alpha);
328 using std::cos;
329 double cosa = cos(alpha);
330 double x2 = v.X() * cosa - v.Y() * sina;
331 double y2 = v.Y() * cosa + v.X() * sina;
332 Vector vrot;
333 vrot.SetXYZ(x2, y2, v.Z());
334 return vrot;
335 }
336
337
338 /**
339 rotation on a generic vector using a generic rotation class.
340 The only requirement on the vector is that implements the
341 X(), Y(), Z() and SetXYZ methods.
342 The requirement on the rotation matrix is that need to implement the
343 (i,j) operator returning the matrix element with R(0,0) = xx element
344 */
345 template<class Vector, class RotationMatrix>
346 Vector Rotate(const Vector &v, const RotationMatrix & rot) {
347 double xX = v.X();
348 double yY = v.Y();
349 double zZ = v.Z();
350 double x2 = rot(0,0)*xX + rot(0,1)*yY + rot(0,2)*zZ;
351 double y2 = rot(1,0)*xX + rot(1,1)*yY + rot(1,2)*zZ;
352 double z2 = rot(2,0)*xX + rot(2,1)*yY + rot(2,2)*zZ;
353 Vector vrot;
354 vrot.SetXYZ(x2,y2,z2);
355 return vrot;
356 }
357
358 /**
359 Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost
360 The only requirement on the vector is that implements the
361 X(), Y(), Z(), T() and SetXYZT methods.
362 The requirement on the boost vector is that needs to implement the
363 X(), Y() , Z() retorning the vector elements describing the boost
364 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
365 */
366 template <class LVector, class BoostVector>
367 LVector boost(const LVector & v, const BoostVector & b) {
368 double bx = b.X();
369 double by = b.Y();
370 double bz = b.Z();
371 double b2 = bx*bx + by*by + bz*bz;
372 if (b2 >= 1) {
373 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
374 return LVector();
375 }
376 using std::sqrt;
377 double gamma = 1.0 / sqrt(1.0 - b2);
378 double bp = bx*v.X() + by*v.Y() + bz*v.Z();
379 double gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
380 double x2 = v.X() + gamma2*bp*bx + gamma*bx*v.T();
381 double y2 = v.Y() + gamma2*bp*by + gamma*by*v.T();
382 double z2 = v.Z() + gamma2*bp*bz + gamma*bz*v.T();
383 double t2 = gamma*(v.T() + bp);
384 LVector lv;
385 lv.SetXYZT(x2,y2,z2,t2);
386 return lv;
387 }
388
389
390 /**
391 Boost a generic Lorentz Vector class along the X direction with a factor beta
392 The only requirement on the vector is that implements the
393 X(), Y(), Z(), T() and SetXYZT methods.
394 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
395 */
396 template <class LVector, class T>
397 LVector boostX(const LVector & v, T beta) {
398 if (beta >= 1) {
399 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
400 return LVector();
401 }
402 using std::sqrt;
403 T gamma = 1.0/ sqrt(1.0 - beta*beta);
404 typename LVector::Scalar x2 = gamma * v.X() + gamma * beta * v.T();
405 typename LVector::Scalar t2 = gamma * beta * v.X() + gamma * v.T();
406
407 LVector lv;
408 lv.SetXYZT(x2,v.Y(),v.Z(),t2);
409 return lv;
410 }
411
412 /**
413 Boost a generic Lorentz Vector class along the Y direction with a factor beta
414 The only requirement on the vector is that implements the
415 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
416 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
417 */
418 template <class LVector>
419 LVector boostY(const LVector & v, double beta) {
420 if (beta >= 1) {
421 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
422 return LVector();
423 }
424 using std::sqrt;
425 double gamma = 1.0/ sqrt(1.0 - beta*beta);
426 double y2 = gamma * v.Y() + gamma * beta * v.T();
427 double t2 = gamma * beta * v.Y() + gamma * v.T();
428 LVector lv;
429 lv.SetXYZT(v.X(),y2,v.Z(),t2);
430 return lv;
431 }
432
433 /**
434 Boost a generic Lorentz Vector class along the Z direction with a factor beta
435 The only requirement on the vector is that implements the
436 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
437 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
438 */
439 template <class LVector>
440 LVector boostZ(const LVector & v, double beta) {
441 if (beta >= 1) {
442 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
443 return LVector();
444 }
445 using std::sqrt;
446 double gamma = 1.0/ sqrt(1.0 - beta*beta);
447 double z2 = gamma * v.Z() + gamma * beta * v.T();
448 double t2 = gamma * beta * v.Z() + gamma * v.T();
449 LVector lv;
450 lv.SetXYZT(v.X(),v.Y(),z2,t2);
451 return lv;
452 }
453
454#endif
455
456
457
458
459 // MATRIX VECTOR MULTIPLICATION
460 // cannot define an operator * otherwise conflicts with rotations
461 // operations like Rotation3D * vector use Mult
462
463 /**
464 Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
465 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 !!
466 */
467 template<class Matrix, class CoordSystem, class U>
468 inline
471 vret.SetXYZ( m(0,0) * v.x() + m(0,1) * v.y() + m(0,2) * v.z() ,
472 m(1,0) * v.x() + m(1,1) * v.y() + m(1,2) * v.z() ,
473 m(2,0) * v.x() + m(2,1) * v.y() + m(2,2) * v.z() );
474 return vret;
475 }
476
477
478 /**
479 Multiplications of a generic matrices with a generic PositionVector
480 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 !!
481 */
482 template<class Matrix, class CoordSystem, class U>
483 inline
486 pret.SetXYZ( m(0,0) * p.x() + m(0,1) * p.y() + m(0,2) * p.z() ,
487 m(1,0) * p.x() + m(1,1) * p.y() + m(1,2) * p.z() ,
488 m(2,0) * p.x() + m(2,1) * p.y() + m(2,2) * p.z() );
489 return pret;
490 }
491
492
493 /**
494 Multiplications of a generic matrices with a LorentzVector described
495 in any coordinate system.
496 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 !!
497 */
498 // this will not be ambigous with operator*(Scalar, LorentzVector) since that one // Scalar is passed by value
499 template<class CoordSystem, class Matrix>
500 inline
503 vret.SetXYZT( m(0,0)*v.x() + m(0,1)*v.y() + m(0,2)*v.z() + m(0,3)* v.t() ,
504 m(1,0)*v.x() + m(1,1)*v.y() + m(1,2)*v.z() + m(1,3)* v.t() ,
505 m(2,0)*v.x() + m(2,1)*v.y() + m(2,2)*v.z() + m(2,3)* v.t() ,
506 m(3,0)*v.x() + m(3,1)*v.y() + m(3,2)*v.z() + m(3,3)* v.t() );
507 return vret;
508 }
509
510
511
512 // non-template utility functions for all objects
513
514
515 /**
516 Return a phi angle in the interval (0,2*PI]
517 */
518 double Phi_0_2pi(double phi);
519 /**
520 Returns phi angle in the interval (-PI,PI]
521 */
522 double Phi_mpi_pi(double phi);
523
524
525
526 } // end namespace Vector Util
527
528
529
530 } // end namespace Math
531
532} // end namespace ROOT
533
534
535#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:469
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 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:419
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:126
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: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:367
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:325
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: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:306
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:287
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:346
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:440
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:397
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
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
TLine lv
Definition textalign.C:5
TMarker m
Definition textangle.C:8