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