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 * Copyright (c) 2005 , LCG ROOT MathLib Team *
7 * *
8 * *
9 **********************************************************************/
10
11// Header file for Vector Utility functions
12//
13// Created by: moneta at Tue May 31 21:10:29 2005
14//
15// Last update: Tue May 31 21:10:29 2005
16//
17#ifndef ROOT_MathX_GenVectorX_VectorUtil
18#define ROOT_MathX_GenVectorX_VectorUtil 1
19
20#include "Math/Math.h"
21
23
25
27
28namespace ROOT {
29
30namespace ROOT_MATH_ARCH {
31
32// utility functions for vector classes
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 GenVectorX
41
42 @see GenVectorX
43 */
44
45namespace VectorUtil {
46
47/// \addtogroup GenVector
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 */
60template <class Vector1, class Vector2>
61inline typename Vector1::Scalar DeltaPhi(const Vector1 &v1, const Vector2 &v2)
62{
63 typename Vector1::Scalar dphi = v2.Phi() - v1.Phi();
64 if (dphi > M_PI) {
65 dphi -= 2.0 * M_PI;
66 } else if (dphi <= -M_PI) {
67 dphi += 2.0 * M_PI;
68 }
69 return dphi;
70}
71
72/**
73 Find square of the difference in pseudorapidity (Eta) and Phi between two generic vectors
74 The only requirements on the Vector classes is that they implement the Phi() and Eta() method
75 \param v1 Vector 1
76 \param v2 Vector 2
77 \return Angle between the two vectors
78 \f[ \Delta R2 = ( \Delta \phi )^2 + ( \Delta \eta )^2 \f]
79 */
80template <class Vector1, class Vector2>
81inline typename Vector1::Scalar DeltaR2(const Vector1 &v1, const Vector2 &v2)
82{
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 */
96template <class Vector1, class Vector2>
97inline typename Vector1::Scalar DeltaR2RapidityPhi(const Vector1 &v1, const Vector2 &v2)
98{
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 */
112template <class Vector1, class Vector2>
113inline typename Vector1::Scalar DeltaR(const Vector1 &v1, const Vector2 &v2)
114{
115 return math_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 */
126template <class Vector1, class Vector2>
127inline typename Vector1::Scalar DeltaRapidityPhi(const Vector1 &v1, const Vector2 &v2)
128{
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 ??
142template <class Vector1, class Vector2>
143double CosTheta(const Vector1 &v1, const Vector2 &v2)
144{
145 double arg;
146 double v1_r2 = v1.X() * v1.X() + v1.Y() * v1.Y() + v1.Z() * v1.Z();
147 double v2_r2 = v2.X() * v2.X() + v2.Y() * v2.Y() + v2.Z() * v2.Z();
148 double ptot2 = v1_r2 * v2_r2;
149 if (ptot2 <= 0) {
150 arg = 0.0;
151 } else {
152 double pdot = v1.X() * v2.X() + v1.Y() * v2.Y() + v1.Z() * v2.Z();
153 arg = pdot / math_sqrt(ptot2);
154 if (arg > 1.0)
155 arg = 1.0;
156 if (arg < -1.0)
157 arg = -1.0;
158 }
159 return arg;
160}
161
162/**
163 Find Angle between two vectors.
164 Use the CosTheta() function
165 \param v1 Vector v1
166 \param v2 Vector v2
167 \return Angle between the two vectors
168 \f[ \theta = \cos ^{-1} \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
169 */
170template <class Vector1, class Vector2>
171inline double Angle(const Vector1 &v1, const Vector2 &v2)
172{
173 return math_acos(CosTheta(v1, v2));
174}
175
176/**
177 Find the projection of v along the given direction u.
178 \param v Vector v for which the propjection is to be found
179 \param u Vector specifying the direction
180 \return Vector projection (same type of v)
181 \f[ \vec{proj} = \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
182 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
183 */
184template <class Vector1, class Vector2>
186{
187 double magU2 = u.X() * u.X() + u.Y() * u.Y() + u.Z() * u.Z();
188 if (magU2 == 0)
189 return Vector1(0, 0, 0);
190 double d = v.Dot(u) / magU2;
191 return Vector1(u.X() * d, u.Y() * d, u.Z() * d);
192}
193
194/**
195 Find the vector component of v perpendicular to the given direction of u
196 \param v Vector v for which the perpendicular component is to be found
197 \param u Vector specifying the direction
198 \return Vector component of v which is perpendicular to u
199 \f[ \vec{perp} = \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
200 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
201 */
202template <class Vector1, class Vector2>
203inline Vector1 PerpVector(const Vector1 &v, const Vector2 &u)
204{
205 return v - ProjVector(v, u);
206}
207
208/**
209 Find the magnitude square of the vector component of v perpendicular to the given direction of u
210 \param v Vector v for which the perpendicular component is to be found
211 \param u Vector specifying the direction
212 \return square value of the component of v which is perpendicular to u
213 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} |^2 \f]
214 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
215 */
216template <class Vector1, class Vector2>
217inline double Perp2(const Vector1 &v, const Vector2 &u)
218{
219 double magU2 = u.X() * u.X() + u.Y() * u.Y() + u.Z() * u.Z();
220 double prjvu = v.Dot(u);
221 double magV2 = v.Dot(v);
222 return magU2 > 0.0 ? magV2 - prjvu * prjvu / magU2 : magV2;
223}
224
225/**
226 Find the magnitude of the vector component of v perpendicular to the given direction of u
227 \param v Vector v for which the perpendicular component is to be found
228 \param u Vector specifying the direction
229 \return value of the component of v which is perpendicular to u
230 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} | \f]
231 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
232 */
233template <class Vector1, class Vector2>
234inline double Perp(const Vector1 &v, const Vector2 &u)
235{
236 return math_sqrt(Perp2(v, u));
237}
238
239// Lorentz Vector functions
240
241/**
242 return the invariant mass of two LorentzVector
243 The only requirement on the LorentzVector is that they need to implement the
244 X() , Y(), Z() and E() methods.
245 \param v1 LorenzVector 1
246 \param v2 LorenzVector 2
247 \return invariant mass M
248 \f[ M_{12} = \sqrt{ (\vec{v1} + \vec{v2} ) \cdot (\vec{v1} + \vec{v2} ) } \f]
249 */
250template <class Vector1, class Vector2>
251inline typename Vector1::Scalar InvariantMass(const Vector1 &v1, const Vector2 &v2)
252{
253 typedef typename Vector1::Scalar Scalar;
254 Scalar ee = (v1.E() + v2.E());
255 Scalar xx = (v1.X() + v2.X());
256 Scalar yy = (v1.Y() + v2.Y());
257 Scalar zz = (v1.Z() + v2.Z());
258 Scalar mm2 = ee * ee - xx * xx - yy * yy - zz * zz;
259 return mm2 < 0.0 ? -math_sqrt(-mm2) : math_sqrt(mm2);
260 // PxPyPzE4D<double> q(xx,yy,zz,ee);
261 // return q.M();
262 // return ( v1 + v2).mag();
263}
264
265/// @brief Returns the square of what InvariantMass(const Vector1&, const Vector2&) would return.
266/// This is mostly useful for speed optimisations, because the expensive sqrt operation is skipped.
267template <class Vector1, class Vector2>
268inline typename Vector1::Scalar InvariantMass2(const Vector1 &v1, const Vector2 &v2)
269{
270 typedef typename Vector1::Scalar Scalar;
271 Scalar ee = (v1.E() + v2.E());
272 Scalar xx = (v1.X() + v2.X());
273 Scalar yy = (v1.Y() + v2.Y());
274 Scalar zz = (v1.Z() + v2.Z());
275 Scalar mm2 = ee * ee - xx * xx - yy * yy - zz * zz;
276 return mm2; // < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
277 // PxPyPzE4D<double> q(xx,yy,zz,ee);
278 // return q.M();
279 // return ( v1 + v2).mag();
280}
281
282// rotation and transformations
283
284#ifndef __CINT__
285/**
286 rotation along X axis for a generic vector by an Angle alpha
287 returning a new vector.
288 The only pre requisite on the Vector is that it has to implement the X(), Y() and Z()
289 and SetXYZ methods.
290 */
291template <class Vector>
292Vector RotateX(const Vector &v, double alpha)
293{
294 if (math_fmod(alpha, 2 * M_PI) == 0.)
295 return v;
296 double sina = math_sin(alpha);
297 double cosa = math_cos(alpha);
298 double y2 = v.Y() * cosa - v.Z() * sina;
299 double z2 = v.Z() * cosa + v.Y() * sina;
300 Vector vrot;
301 vrot.SetXYZ(v.X(), y2, z2);
302 return vrot;
303}
304
305/**
306 rotation along Y axis for a generic vector by an Angle alpha
307 returning a new vector.
308 The only pre requisite on the Vector is that it has to implement the X(), Y() and Z()
309 and SetXYZ methods.
310 */
311template <class Vector>
312Vector RotateY(const Vector &v, double alpha)
313{
314 if (math_fmod(alpha, 2 * M_PI) == 0.)
315 return v;
316 double sina = math_sin(alpha);
317 double cosa = math_cos(alpha);
318 double x2 = v.X() * cosa + v.Z() * sina;
319 double z2 = v.Z() * cosa - v.X() * sina;
320 Vector vrot;
321 vrot.SetXYZ(x2, v.Y(), z2);
322 return vrot;
323}
324
325/**
326 rotation along Z axis for a generic vector by an Angle alpha
327 returning a new vector.
328 The only pre requisite on the Vector is that it has to implement the X(), Y() and Z()
329 and SetXYZ methods.
330 */
331template <class Vector>
332Vector RotateZ(const Vector &v, double alpha)
333{
334 if (math_fmod(alpha, 2 * M_PI) == 0.)
335 return v;
336 double sina = math_sin(alpha);
337 double cosa = math_cos(alpha);
338 double x2 = v.X() * cosa - v.Y() * sina;
339 double y2 = v.Y() * cosa + v.X() * sina;
340 Vector vrot;
341 vrot.SetXYZ(x2, y2, v.Z());
342 return vrot;
343}
344
345/**
346 rotation along a custom axis for a generic vector by an Angle alpha (in rad)
347 returning a new vector.
348 The only pre requisite on the Vector is that it has to implement the X(), Y() and Z()
349 and SetXYZ methods.
350 */
351template <class Vector>
352Vector Rotate(const Vector &v, double alpha, const Vector &axis)
353{
354 if (math_fmod(alpha, 2 * M_PI) == 0.)
355 return v;
356 const double ll = math_sqrt(axis.X() * axis.X() + axis.Y() * axis.Y() + axis.Z() * axis.Z());
357#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
358 if (ll == 0.)
359 GenVector_Throw("Axis Vector has zero magnitude");
360#endif
361 const double sa = math_sin(alpha);
362 const double ca = math_cos(alpha);
363 const double dx = axis.X() / ll;
364 const double dy = axis.Y() / ll;
365 const double dz = axis.Z() / ll;
366 // clang-format off
367 const double rot00 = (1 - ca) * dx * dx + ca , rot01 = (1 - ca) * dx * dy - sa * dz, rot02 = (1 - ca) * dx * dz + sa * dy,
368 rot10 = (1 - ca) * dy * dx + sa * dz, rot11 = (1 - ca) * dy * dy + ca , rot12 = (1 - ca) * dy * dz - sa * dx,
369 rot20 = (1 - ca) * dz * dx - sa * dy, rot21 = (1 - ca) * dz * dy + sa * dx, rot22 = (1 - ca) * dz * dz + ca ;
370 // clang-format on
371 const double xX = v.X();
372 const double yY = v.Y();
373 const double zZ = v.Z();
374 const double x2 = rot00 * xX + rot01 * yY + rot02 * zZ;
375 const double y2 = rot10 * xX + rot11 * yY + rot12 * zZ;
376 const double z2 = rot20 * xX + rot21 * yY + rot22 * zZ;
377 Vector vrot;
378 vrot.SetXYZ(x2, y2, z2);
379 return vrot;
380}
381
382/**
383 rotation on a generic vector using a generic rotation class.
384 The only requirement on the vector is that implements the
385 X(), Y(), Z() and SetXYZ methods.
386 The requirement on the rotation matrix is that need to implement the
387 (i,j) operator returning the matrix element with R(0,0) = xx element
388 */
389template <class Vector, class RotationMatrix>
390Vector Rotate(const Vector &v, const RotationMatrix &rot)
391{
392 double xX = v.X();
393 double yY = v.Y();
394 double zZ = v.Z();
395 double x2 = rot(0, 0) * xX + rot(0, 1) * yY + rot(0, 2) * zZ;
396 double y2 = rot(1, 0) * xX + rot(1, 1) * yY + rot(1, 2) * zZ;
397 double z2 = rot(2, 0) * xX + rot(2, 1) * yY + rot(2, 2) * zZ;
398 Vector vrot;
399 vrot.SetXYZ(x2, y2, z2);
400 return vrot;
401}
402
403/**
404 Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost
405 The only requirement on the vector is that implements the
406 X(), Y(), Z(), T() and SetXYZT methods.
407 The requirement on the boost vector is that needs to implement the
408 X(), Y() , Z() retorning the vector elements describing the boost
409 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
410 */
411template <class LVector, class BoostVector>
412LVector boost(const LVector &v, const BoostVector &b)
413{
414 double bx = b.X();
415 double by = b.Y();
416 double bz = b.Z();
417 double b2 = bx * bx + by * by + bz * bz;
418 if (b2 >= 1) {
419#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
420 GenVector_Throw("Beta Vector supplied to set Boost represents speed >= c");
421#endif
422 return LVector();
423 }
424 double gamma = 1.0 / math_sqrt(1.0 - b2);
425 double bp = bx * v.X() + by * v.Y() + bz * v.Z();
426 double gamma2 = b2 > 0 ? (gamma - 1.0) / b2 : 0.0;
427 double x2 = v.X() + gamma2 * bp * bx + gamma * bx * v.T();
428 double y2 = v.Y() + gamma2 * bp * by + gamma * by * v.T();
429 double z2 = v.Z() + gamma2 * bp * bz + gamma * bz * v.T();
430 double t2 = gamma * (v.T() + bp);
431 LVector lv;
432 lv.SetXYZT(x2, y2, z2, t2);
433 return lv;
434}
435
436/**
437 Boost a generic Lorentz Vector class along the X direction with a factor beta
438 The only requirement on the vector is that implements the
439 X(), Y(), Z(), T() and SetXYZT methods.
440 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
441 */
442template <class LVector, class T>
443LVector boostX(const LVector &v, T beta)
444{
445 if (beta >= 1) {
446#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
447 GenVector_Throw("Beta Vector supplied to set Boost represents speed >= c");
448#endif
449 return LVector();
450 }
451 T gamma = 1.0 / math_sqrt(1.0 - beta * beta);
452 typename LVector::Scalar x2 = gamma * v.X() + gamma * beta * v.T();
453 typename LVector::Scalar t2 = gamma * beta * v.X() + gamma * v.T();
454
455 LVector lv;
456 lv.SetXYZT(x2, v.Y(), v.Z(), t2);
457 return lv;
458}
459
460/**
461 Boost a generic Lorentz Vector class along the Y direction with a factor beta
462 The only requirement on the vector is that implements the
463 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
464 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
465 */
466template <class LVector>
467LVector boostY(const LVector &v, double beta)
468{
469 if (beta >= 1) {
470#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
471 GenVector_Throw("Beta Vector supplied to set Boost represents speed >= c");
472#endif
473 return LVector();
474 }
475 double gamma = 1.0 / math_sqrt(1.0 - beta * beta);
476 double y2 = gamma * v.Y() + gamma * beta * v.T();
477 double t2 = gamma * beta * v.Y() + gamma * v.T();
478 LVector lv;
479 lv.SetXYZT(v.X(), y2, v.Z(), t2);
480 return lv;
481}
482
483/**
484 Boost a generic Lorentz Vector class along the Z direction with a factor beta
485 The only requirement on the vector is that implements the
486 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
487 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
488 */
489template <class LVector>
490LVector boostZ(const LVector &v, double beta)
491{
492 if (beta >= 1) {
493#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
494 GenVector_Throw("Beta Vector supplied to set Boost represents speed >= c");
495#endif
496 return LVector();
497 }
498 double gamma = 1.0 / math_sqrt(1.0 - beta * beta);
499 double z2 = gamma * v.Z() + gamma * beta * v.T();
500 double t2 = gamma * beta * v.Z() + gamma * v.T();
501 LVector lv;
502 lv.SetXYZT(v.X(), v.Y(), z2, t2);
503 return lv;
504}
505
506#endif
507
508// MATRIX VECTOR MULTIPLICATION
509// cannot define an operator * otherwise conflicts with rotations
510// operations like Rotation3D * vector use Mult
511
512/**
513 Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
514 Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no
515 check on the matrix size !!
516 */
517template <class Matrix, class CoordSystem, class U>
519{
521 vret.SetXYZ(m(0, 0) * v.x() + m(0, 1) * v.y() + m(0, 2) * v.z(), m(1, 0) * v.x() + m(1, 1) * v.y() + m(1, 2) * v.z(),
522 m(2, 0) * v.x() + m(2, 1) * v.y() + m(2, 2) * v.z());
523 return vret;
524}
525
526/**
527 Multiplications of a generic matrices with a generic PositionVector
528 Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no
529 check on the matrix size !!
530 */
531template <class Matrix, class CoordSystem, class U>
533{
535 pret.SetXYZ(m(0, 0) * p.x() + m(0, 1) * p.y() + m(0, 2) * p.z(), m(1, 0) * p.x() + m(1, 1) * p.y() + m(1, 2) * p.z(),
536 m(2, 0) * p.x() + m(2, 1) * p.y() + m(2, 2) * p.z());
537 return pret;
538}
539
540/**
541 Multiplications of a generic matrices with a LorentzVector described
542 in any coordinate system.
543 Assume that the matrix implements the operator( i,j) and that it has at least 4 columns and 4 rows. There is no
544 check on the matrix size !!
545 */
546// this will not be ambiguous with operator*(Scalar, LorentzVector) since that one // Scalar is passed by value
547template <class CoordSystem, class Matrix>
549{
551 vret.SetXYZT(m(0, 0) * v.x() + m(0, 1) * v.y() + m(0, 2) * v.z() + m(0, 3) * v.t(),
552 m(1, 0) * v.x() + m(1, 1) * v.y() + m(1, 2) * v.z() + m(1, 3) * v.t(),
553 m(2, 0) * v.x() + m(2, 1) * v.y() + m(2, 2) * v.z() + m(2, 3) * v.t(),
554 m(3, 0) * v.x() + m(3, 1) * v.y() + m(3, 2) * v.z() + m(3, 3) * v.t());
555 return vret;
556}
557
558// non-template utility functions for all objects
559
560/**
561 Return a phi angle in the interval (0,2*PI]
562 */
563double Phi_0_2pi(double phi);
564/**
565 Returns phi angle in the interval (-PI,PI]
566 */
567double Phi_mpi_pi(double phi);
568
569// Close of doxygen group:
570/// @}
571
572} // namespace VectorUtil
573
574} // namespace ROOT_MATH_ARCH
575
576} // end namespace ROOT
577
578#endif /* ROOT_MathX_GenVectorX_VectorUtil */
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define M_PI
Definition Rotated.cxx:105
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
double Phi_mpi_pi(double phi)
Returns phi angle in the interval (-PI,PI].
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:217
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:81
Vector1 PerpVector(const Vector1 &v, const Vector2 &u)
Find the vector component of v perpendicular to the given direction of u.
Definition VectorUtil.h:203
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
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:490
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
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:292
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
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:234
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:518
double Angle(const Vector1 &v1, const Vector2 &v2)
Find Angle between two vectors.
Definition VectorUtil.h:171
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
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:412
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:352
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:443
Vector1::Scalar InvariantMass2(const Vector1 &v1, const Vector2 &v2)
Returns the square of what InvariantMass(const Vector1&, const Vector2&) would return.
Definition VectorUtil.h:268
double Phi_0_2pi(double phi)
Return a phi angle in the interval (0,2*PI].
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:332
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:251
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:467
Vector1 ProjVector(const Vector1 &v, const Vector2 &u)
Find the projection of v along the given direction u.
Definition VectorUtil.h:185
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:312
Scalar math_cos(Scalar x)
Rotation3D::Scalar Scalar
void GenVector_Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
Scalar math_sqrt(Scalar x)
Scalar math_acos(Scalar x)
Scalar math_fmod(Scalar x, Scalar y)
Scalar math_sin(Scalar x)
TLine lv
Definition textalign.C:5
TMarker m
Definition textangle.C:8