Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
LorentzVector.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
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 class LorentzVector
12//
13// Created by: moneta at Tue May 31 17:06:09 2005
14// Major mods by: fischler at Wed Jul 20 2005
15//
16// Last update: $Id$
17//
18#ifndef ROOT_MathX_GenVectorX_LorentzVector
19#define ROOT_MathX_GenVectorX_LorentzVector 1
20
24#include "MathX/Vector2D.h"
25
27
29
30using namespace ROOT::ROOT_MATH_ARCH;
31
32#include <cmath>
33#include <string>
34
35namespace ROOT {
36
37namespace ROOT_MATH_ARCH {
38
39//__________________________________________________________________________________________
40/** @ingroup GenVectorX
41
42Class describing a generic LorentzVector in the 4D space-time,
43using the specified coordinate system for the spatial vector part.
44The metric used for the LorentzVector is (-,-,-,+).
45In the case of LorentzVector we don't distinguish the concepts
46of points and displacement vectors as in the 3D case,
47since the main use case for 4D Vectors is to describe the kinematics of
48relativistic particles. A LorentzVector behaves like a
49DisplacementVector in 4D. The Minkowski components could be viewed as
50v and t, or for kinematic 4-vectors, as p and E.
51
52ROOT provides specialisations (and aliases to them) of the ROOT::ROOT_MATH_ARCH::LorentzVector template:
53- ROOT::ROOT_MATH_ARCH::PtEtaPhiMVector based on pt (rho),eta,phi and M (t) coordinates in double precision
54- ROOT::ROOT_MATH_ARCH::PtEtaPhiEVector based on pt (rho),eta,phi and E (t) coordinates in double precision
55- ROOT::ROOT_MATH_ARCH::PxPyPzMVector based on px,py,pz and M (mass) coordinates in double precision
56- ROOT::ROOT_MATH_ARCH::PxPyPzEVector based on px,py,pz and E (energy) coordinates in double precision
57- ROOT::ROOT_MATH_ARCH::XYZTVector based on x,y,z,t coordinates (cartesian) in double precision (same as PxPyPzEVector)
58- ROOT::ROOT_MATH_ARCH::XYZTVectorF based on x,y,z,t coordinates (cartesian) in float precision (same as PxPyPzEVector
59but float)
60Pt (or rho) refers to transverse momentum, whereas eta refers to pseudorapidity. M is mass, E is energy, p is momentum.
61
62@see GenVectorX
63*/
64
65template <class CoordSystem>
67
68public:
69 // ------ ctors ------
70
71 typedef typename CoordSystem::Scalar Scalar;
72 typedef CoordSystem CoordinateType;
73
74 /**
75 default constructor of an empty vector (Px = Py = Pz = E = 0 )
76 */
78
79 /**
80 generic constructors from four scalar values.
81 The association between values and coordinate depends on the
82 coordinate system. For PxPyPzE4D,
83 \param a scalar value (Px)
84 \param b scalar value (Py)
85 \param c scalar value (Pz)
86 \param d scalar value (E)
87 */
88 LorentzVector(const Scalar &a, const Scalar &b, const Scalar &c, const Scalar &d) : fCoordinates(a, b, c, d) {}
89
90 /**
91 constructor from a LorentzVector expressed in different
92 coordinates, or using a different Scalar type
93 */
94 template <class Coords>
98
99 /**
100 Construct from a foreign 4D vector type, for example, HepLorentzVector
101 Precondition: v must implement methods x(), y(), z(), and t()
102 */
103 template <class ForeignLorentzVector,
104 typename = decltype(std::declval<ForeignLorentzVector>().x() + std::declval<ForeignLorentzVector>().y() +
105 std::declval<ForeignLorentzVector>().z() + std::declval<ForeignLorentzVector>().t())>
106 explicit LorentzVector(const ForeignLorentzVector &v) : fCoordinates(PxPyPzE4D<Scalar>(v.x(), v.y(), v.z(), v.t()))
107 {
108 }
109
110#ifdef LATER
111 /**
112 construct from a generic linear algebra vector implementing operator []
113 and with a size of at least 4. This could be also a C array
114 In this case v[0] is the first data member
115 ( Px for a PxPyPzE4D base)
116 \param v LA vector
117 \param index0 index of first vector element (Px)
118 */
119 template <class LAVector>
120 explicit LorentzVector(const LAVector &v, size_t index0)
121 {
122 fCoordinates = CoordSystem(v[index0], v[index0 + 1], v[index0 + 2], v[index0 + 3]);
123 }
124#endif
125
126 // ------ assignment ------
127
128 /**
129 Assignment operator from a lorentz vector of arbitrary type
130 */
131 template <class OtherCoords>
133 {
134 fCoordinates = v.Coordinates();
135 return *this;
136 }
137
138 /**
139 assignment from any other Lorentz vector implementing
140 x(), y(), z() and t()
141 */
142 template <class ForeignLorentzVector,
143 typename = decltype(std::declval<ForeignLorentzVector>().x() + std::declval<ForeignLorentzVector>().y() +
144 std::declval<ForeignLorentzVector>().z() + std::declval<ForeignLorentzVector>().t())>
146 {
147 SetXYZT(v.x(), v.y(), v.z(), v.t());
148 return *this;
149 }
150
151#ifdef LATER
152 /**
153 assign from a generic linear algebra vector implementing operator []
154 and with a size of at least 4
155 In this case v[0] is the first data member
156 ( Px for a PxPyPzE4D base)
157 \param v LA vector
158 \param index0 index of first vector element (Px)
159 */
160 template <class LAVector>
161 LorentzVector &AssignFrom(const LAVector &v, size_t index0 = 0)
162 {
163 fCoordinates.SetCoordinates(v[index0], v[index0 + 1], v[index0 + 2], v[index0 + 3]);
164 return *this;
165 }
166#endif
167
168 // ------ Set, Get, and access coordinate data ------
169
170 /**
171 Retrieve a const reference to the coordinates object
172 */
173 const CoordSystem &Coordinates() const { return fCoordinates; }
174
175 /**
176 Set internal data based on an array of 4 Scalar numbers
177 */
179 {
180 fCoordinates.SetCoordinates(src);
181 return *this;
182 }
183
184 /**
185 Set internal data based on 4 Scalar numbers
186 */
188 {
189 fCoordinates.SetCoordinates(a, b, c, d);
190 return *this;
191 }
192
193 /**
194 Set internal data based on 4 Scalars at *begin to *end
195 */
196 template <class IT>
198 {
199 IT a = begin;
200 IT b = ++begin;
201 IT c = ++begin;
202 IT d = ++begin;
203 (void)end;
204 assert(++begin == end);
205 SetCoordinates(*a, *b, *c, *d);
206 return *this;
207 }
208
209 /**
210 get internal data into 4 Scalar numbers
211 */
212 void GetCoordinates(Scalar &a, Scalar &b, Scalar &c, Scalar &d) const { fCoordinates.GetCoordinates(a, b, c, d); }
213
214 /**
215 get internal data into an array of 4 Scalar numbers
216 */
217 void GetCoordinates(Scalar dest[]) const { fCoordinates.GetCoordinates(dest); }
218
219 /**
220 get internal data into 4 Scalars at *begin to *end
221 */
222 template <class IT>
223 void GetCoordinates(IT begin, IT end) const
224 {
225 IT a = begin;
226 IT b = ++begin;
227 IT c = ++begin;
228 IT d = ++begin;
229 (void)end;
230 assert(++begin == end);
231 GetCoordinates(*a, *b, *c, *d);
232 }
233
234 /**
235 get internal data into 4 Scalars at *begin
236 */
237 template <class IT>
238 void GetCoordinates(IT begin) const
239 {
240 Scalar a, b, c, d = 0;
241 GetCoordinates(a, b, c, d);
242 *begin++ = a;
243 *begin++ = b;
244 *begin++ = c;
245 *begin = d;
246 }
247
248 /**
249 set the values of the vector from the cartesian components (x,y,z,t)
250 (if the vector is held in another coordinates, like (Pt,eta,phi,m)
251 then (x, y, z, t) are converted to that form)
252 */
254 {
255 fCoordinates.SetPxPyPzE(xx, yy, zz, tt);
256 return *this;
257 }
259 {
260 fCoordinates.SetPxPyPzE(xx, yy, zz, ee);
261 return *this;
262 }
263
264 // ------------------- Equality -----------------
265
266 /**
267 Exact equality
268 */
269 bool operator==(const LorentzVector &rhs) const { return fCoordinates == rhs.fCoordinates; }
270 bool operator!=(const LorentzVector &rhs) const { return !(operator==(rhs)); }
271
272 // ------ Individual element access, in various coordinate systems ------
273
274 // individual coordinate accessors in various coordinate systems
275
276 /**
277 spatial X component
278 */
279 Scalar Px() const { return fCoordinates.Px(); }
280 Scalar X() const { return fCoordinates.Px(); }
281 /**
282 spatial Y component
283 */
284 Scalar Py() const { return fCoordinates.Py(); }
285 Scalar Y() const { return fCoordinates.Py(); }
286 /**
287 spatial Z component
288 */
289 Scalar Pz() const { return fCoordinates.Pz(); }
290 Scalar Z() const { return fCoordinates.Pz(); }
291 /**
292 return 4-th component (time, or energy for a 4-momentum vector)
293 */
294 Scalar E() const { return fCoordinates.E(); }
295 Scalar T() const { return fCoordinates.E(); }
296 /**
297 return magnitude (mass) squared M2 = T**2 - X**2 - Y**2 - Z**2
298 (we use -,-,-,+ metric)
299 */
300 Scalar M2() const { return fCoordinates.M2(); }
301 /**
302 return magnitude (mass) using the (-,-,-,+) metric.
303 If M2 is negative (space-like vector) a GenVector_exception
304 is suggested and if continuing, - sqrt( -M2) is returned
305 */
306 Scalar M() const { return fCoordinates.M(); }
307 /**
308 return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
309 */
310 Scalar R() const { return fCoordinates.R(); }
311 Scalar P() const { return fCoordinates.R(); }
312 /**
313 return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
314 */
315 Scalar P2() const { return P() * P(); }
316 /**
317 return the square of the transverse spatial component ( X**2 + Y**2 )
318 */
319 Scalar Perp2() const { return fCoordinates.Perp2(); }
320
321 /**
322 return the transverse spatial component sqrt ( X**2 + Y**2 )
323 */
324 Scalar Pt() const { return fCoordinates.Pt(); }
325 Scalar Rho() const { return fCoordinates.Pt(); }
326
327 /**
328 return the transverse mass squared
329 \f[ m_t^2 = E^2 - p{_z}^2 \f]
330 */
331 Scalar Mt2() const { return fCoordinates.Mt2(); }
332
333 /**
334 return the transverse mass
335 \f[ \sqrt{ m_t^2 = E^2 - p{_z}^2} X sign(E^ - p{_z}^2) \f]
336 */
337 Scalar Mt() const { return fCoordinates.Mt(); }
338
339 /**
340 return the transverse energy squared
341 \f[ e_t = \frac{E^2 p_{\perp}^2 }{ |p|^2 } \f]
342 */
343 Scalar Et2() const { return fCoordinates.Et2(); }
344
345 /**
346 return the transverse energy
347 \f[ e_t = \sqrt{ \frac{E^2 p_{\perp}^2 }{ |p|^2 } } X sign(E) \f]
348 */
349 Scalar Et() const { return fCoordinates.Et(); }
350
351 /**
352 azimuthal Angle
353 */
354 Scalar Phi() const { return fCoordinates.Phi(); }
355
356 /**
357 polar Angle
358 */
359 Scalar Theta() const { return fCoordinates.Theta(); }
360
361 /**
362 pseudorapidity
363 \f[ \eta = - \ln { \tan { \frac { \theta} {2} } } \f]
364 */
365 Scalar Eta() const { return fCoordinates.Eta(); }
366
367 /**
368 deltaRapidity between this and vector v
369 \f[ \Delta R = \sqrt { \Delta \eta ^2 + \Delta \phi ^2 } \f]
370 \param useRapidity true to use Rapidity(), false to use Eta()
371 */
372 template <class OtherLorentzVector>
373 Scalar DeltaR(const OtherLorentzVector &v, const bool useRapidity = false) const
374 {
375 const double delta = useRapidity ? Rapidity() - v.Rapidity() : Eta() - v.Eta();
376 double dphi = Phi() - v.Phi();
377 // convert dphi angle to the interval (-PI,PI]
378 if (dphi > M_PI || dphi <= -M_PI) {
379 if (dphi > 0) {
380 int n = static_cast<int>(dphi / (2 * M_PI) + 0.5);
381 dphi -= (2 * M_PI) * n;
382 } else {
383 int n = static_cast<int>(0.5 - dphi / (2 * M_PI));
384 dphi += (2 * M_PI) * n;
385 }
386 }
387 return math_sqrt(delta * delta + dphi * dphi);
388 }
389
390 /**
391 get the spatial components of the Vector in a
392 DisplacementVector based on Cartesian Coordinates
393 */
398
399 // ------ Operations combining two Lorentz vectors ------
400
401 /**
402 scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
403 Enable the product using any other LorentzVector implementing
404 the x(), y() , y() and t() member functions
405 \param q any LorentzVector implementing the x(), y() , z() and t()
406 member functions
407 \return the result of v.q of type according to the base scalar type of v
408 */
409
410 template <class OtherLorentzVector>
412 {
413 return t() * q.t() - x() * q.x() - y() * q.y() - z() * q.z();
414 }
415
416 /**
417 Self addition with another Vector ( v+= q )
418 Enable the addition with any other LorentzVector
419 \param q any LorentzVector implementing the x(), y() , z() and t()
420 member functions
421 */
422 template <class OtherLorentzVector>
424 {
425 SetXYZT(x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t());
426 return *this;
427 }
428
429 /**
430 Self subtraction of another Vector from this ( v-= q )
431 Enable the addition with any other LorentzVector
432 \param q any LorentzVector implementing the x(), y() , z() and t()
433 member functions
434 */
435 template <class OtherLorentzVector>
437 {
438 SetXYZT(x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t());
439 return *this;
440 }
441
442 /**
443 addition of two LorentzVectors (v3 = v1 + v2)
444 Enable the addition with any other LorentzVector
445 \param v2 any LorentzVector implementing the x(), y() , z() and t()
446 member functions
447 \return a new LorentzVector of the same type as v1
448 */
449 template <class OtherLorentzVector>
451 {
453 v3 += v2;
454 return v3;
455 }
456
457 /**
458 subtraction of two LorentzVectors (v3 = v1 - v2)
459 Enable the subtraction of any other LorentzVector
460 \param v2 any LorentzVector implementing the x(), y() , z() and t()
461 member functions
462 \return a new LorentzVector of the same type as v1
463 */
464 template <class OtherLorentzVector>
466 {
468 v3 -= v2;
469 return v3;
470 }
471
472 //--- scaling operations ------
473
474 /**
475 multiplication by a scalar quantity v *= a
476 */
478 {
479 fCoordinates.Scale(a);
480 return *this;
481 }
482
483 /**
484 division by a scalar quantity v /= a
485 */
487 {
488 fCoordinates.Scale(1 / a);
489 return *this;
490 }
491
492 /**
493 product of a LorentzVector by a scalar quantity
494 \param a scalar quantity of type a
495 \return a new mathcoreLorentzVector q = v * a same type as v
496 */
498 {
499 LorentzVector tmp(*this);
500 tmp *= a;
501 return tmp;
502 }
503
504 /**
505 Divide a LorentzVector by a scalar quantity
506 \param a scalar quantity of type a
507 \return a new mathcoreLorentzVector q = v / a same type as v
508 */
510 {
512 tmp /= a;
513 return tmp;
514 }
515
516 /**
517 Negative of a LorentzVector (q = - v )
518 \return a new LorentzVector with opposite direction and time
519 */
521 {
522 // LorentzVector<CoordinateType> v(*this);
523 // v.Negate();
524 return operator*(Scalar(-1));
525 }
526 LorentzVector operator+() const { return *this; }
527
528 // ---- Relativistic Properties ----
529
530 /**
531 Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)]
532 */
534 {
535 // TODO - It would be good to check that E > Pz and use the Throw()
536 // mechanism or at least load a NAN if not.
537 // We should then move the code to a .cpp file.
538 const Scalar ee = E();
539 const Scalar ppz = Pz();
540 return Scalar(0.5) * math_log((ee + ppz) / (ee - ppz));
541 }
542
543 /**
544 Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
545 */
547 {
548 // TODO - It would be good to check that E > P and use the Throw()
549 // mechanism or at least load a NAN if not.
550 const Scalar ee = E();
551 const Scalar pp = P();
552 return Scalar(0.5) * math_log((ee + pp) / (ee - pp));
553 }
554
555 /**
556 Determine if momentum-energy can represent a physical massive particle
557 */
558 bool isTimelike() const
559 {
560 Scalar ee = E();
561 Scalar pp = P();
562 return ee * ee > pp * pp;
563 }
564
565 /**
566 Determine if momentum-energy can represent a massless particle
567 */
568 bool isLightlike(Scalar tolerance = 100 * std::numeric_limits<Scalar>::epsilon()) const
569 {
570 Scalar ee = E();
571 Scalar pp = P();
572 Scalar delta = ee - pp;
573 if (ee == 0)
574 return pp == 0;
575 return delta * delta < tolerance * ee * ee;
576 }
577
578 /**
579 Determine if momentum-energy is spacelike, and represents a tachyon
580 */
581 bool isSpacelike() const
582 {
583 Scalar ee = E();
584 Scalar pp = P();
585 return ee * ee < pp * pp;
586 }
587
589
590 /**
591 The beta vector for the boost that would bring this vector into
592 its center of mass frame (zero momentum)
593 */
595 {
596 if (E() == 0) {
597 if (P() == 0) {
598 return BetaVector();
599 } else {
600 // TODO - should attempt to Throw with msg about
601 // boostVector computed for LorentzVector with t=0
602 return -Vect() / E();
603 }
604 }
605 if (M2() <= 0) {
606 // TODO - should attempt to Throw with msg about
607 // boostVector computed for a non-timelike LorentzVector
608 }
609 return -Vect() / E();
610 }
611
612 /**
613 The beta vector for the boost that would bring this vector into
614 its center of mass frame (zero momentum)
615 */
616 template <class Other4Vector>
618 {
619 Scalar eSum = E() + v.E();
621 if (eSum == 0) {
622 if (vecSum.Mag2() == 0) {
623 return BetaVector();
624 } else {
625 // TODO - should attempt to Throw with msg about
626 // boostToCM computed for two 4-vectors with combined t=0
627 return BetaVector(vecSum / eSum);
628 }
629 // TODO - should attempt to Throw with msg about
630 // boostToCM computed for two 4-vectors with combined e=0
631 }
632 return BetaVector(vecSum * (-1. / eSum));
633 }
634
635 // beta and gamma
636
637 /**
638 Return beta scalar value
639 */
640 Scalar Beta() const
641 {
642 if (E() == 0) {
643 if (P2() == 0)
644 // to avoid Nan
645 return 0;
646 else {
647#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
649 "LorentzVector::Beta() - beta computed for LorentzVector with t = 0. Return an Infinite result");
650#endif
651 return 1. / E();
652 }
653 }
654 if (M2() <= 0) {
655#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
657 "LorentzVector::Beta() - beta computed for non-timelike LorentzVector . Result is physically meaningless");
658#endif
659 }
660 return P() / E();
661 }
662 /**
663 Return Gamma scalar value
664 */
665 Scalar Gamma() const
666 {
667 const Scalar v2 = P2();
668 const Scalar t2 = E() * E();
669 if (E() == 0) {
670 if (P2() == 0) {
671 return 1;
672 } else {
673#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
675 "LorentzVector::Gamma() - gamma computed for LorentzVector with t = 0. Return a zero result");
676#endif
677 }
678 }
679 if (t2 < v2) {
680#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
681 GenVector_Throw("LorentzVector::Gamma() - gamma computed for a spacelike LorentzVector. Imaginary result");
682#endif
683 return 0;
684 } else if (t2 == v2) {
685#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
686 GenVector_Throw("LorentzVector::Gamma() - gamma computed for a lightlike LorentzVector. Infinite result");
687#endif
688 }
689 return Scalar(1) / math_sqrt(Scalar(1) - v2 / t2);
690 } /* gamma */
691
692 // Method providing limited backward name compatibility with CLHEP ----
693
694 Scalar x() const { return fCoordinates.Px(); }
695 Scalar y() const { return fCoordinates.Py(); }
696 Scalar z() const { return fCoordinates.Pz(); }
697 Scalar t() const { return fCoordinates.E(); }
698 Scalar px() const { return fCoordinates.Px(); }
699 Scalar py() const { return fCoordinates.Py(); }
700 Scalar pz() const { return fCoordinates.Pz(); }
701 Scalar e() const { return fCoordinates.E(); }
702 Scalar r() const { return fCoordinates.R(); }
703 Scalar theta() const { return fCoordinates.Theta(); }
704 Scalar phi() const { return fCoordinates.Phi(); }
705 Scalar rho() const { return fCoordinates.Rho(); }
706 Scalar eta() const { return fCoordinates.Eta(); }
707 Scalar pt() const { return fCoordinates.Pt(); }
708 Scalar perp2() const { return fCoordinates.Perp2(); }
709 Scalar mag2() const { return fCoordinates.M2(); }
710 Scalar mag() const { return fCoordinates.M(); }
711 Scalar mt() const { return fCoordinates.Mt(); }
712 Scalar mt2() const { return fCoordinates.Mt2(); }
713
714 // Methods requested by CMS ---
715 Scalar energy() const { return fCoordinates.E(); }
716 Scalar mass() const { return fCoordinates.M(); }
717 Scalar mass2() const { return fCoordinates.M2(); }
718
719 /**
720 Methods setting a Single-component
721 Work only if the component is one of which the vector is represented.
722 For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
723 */
725 {
726 fCoordinates.SetE(a);
727 return *this;
728 }
730 {
731 fCoordinates.SetEta(a);
732 return *this;
733 }
735 {
736 fCoordinates.SetM(a);
737 return *this;
738 }
740 {
741 fCoordinates.SetPhi(a);
742 return *this;
743 }
745 {
746 fCoordinates.SetPt(a);
747 return *this;
748 }
750 {
751 fCoordinates.SetPx(a);
752 return *this;
753 }
755 {
756 fCoordinates.SetPy(a);
757 return *this;
758 }
760 {
761 fCoordinates.SetPz(a);
762 return *this;
763 }
764
765private:
766 CoordSystem fCoordinates; // internal coordinate system
767
768}; // LorentzVector<>
769
770// global methods
771
772/**
773 Scale of a LorentzVector with a scalar quantity a
774 \param a scalar quantity of type a
775 \param v LorentzVector based on any coordinate system
776 \return a new mathcoreLorentzVector q = v * a same type as v
777 */
778template <class CoordSystem>
781{
783 tmp *= a;
784 return tmp;
785}
786
787/**
788 pair (p+ p-) acoplanarity `alpha = 1 - |phi+ - phi-|/pi`.
789 \param pp p+, LorentzVector based on any coordinate system
790 \param pm p-, LorentzVector based on any coordinate system
791 \return a scalar
792 \ingroup GenVector
793 \see http://doi.org/10.1103/PhysRevLett.121.212301
794*/
795template <class CoordSystem>
798{
799 auto dphi = pp.Phi() - pm.Phi();
800 // convert dphi angle to the interval (-PI,PI]
801 if (dphi > M_PI || dphi <= -M_PI) {
802 if (dphi > 0) {
803 int n = static_cast<int>(dphi / (2 * M_PI) + 0.5);
804 dphi -= (2 * M_PI) * n;
805 } else {
806 int n = static_cast<int>(0.5 - dphi / (2 * M_PI));
807 dphi += (2 * M_PI) * n;
808 }
809 }
810 return 1 - math_fabs(dphi) / M_PI;
811}
812
813/**
814 pair (p+ p-) vectorial asymmetry `Av = |Pt+ - Pt-|/|Pt+ + Pt-|`.
815 In an experimental setting, it reflects a convolution of the experimental resolutions
816 of particle energy and azimuthal angle measurement.
817 \param pp p+, LorentzVector based on any coordinate system
818 \param pm p-, LorentzVector based on any coordinate system
819 \return a scalar. Returns -1 if both momenta are exactly mirrored.
820 \ingroup GenVector
821 \see http://doi.org/10.1103/PhysRevLett.121.212301, https://doi.org/10.1103/PhysRevD.99.093013
822*/
823template <class CoordSystem>
826{
827 ROOT::ROOT_MATH_ARCH::XYVector vp(pp.Px(), pp.Py());
829 auto denom = (vp + vm).R();
830 if (denom == 0.)
831 return -1;
832 return (vp - vm).R() / denom;
833}
834
835/**
836 pair (p+ p-) scalar asymmetry `As = ||Pt+| - |Pt-|/||Pt+| + |Pt-||`.
837 Measures the relative difference in transverse momentum of the pair, e.g. two photons,
838 and would be ideally zero for two back-to-back photons.
839 \param pp p+, LorentzVector based on any coordinate system
840 \param pm p-, LorentzVector based on any coordinate system
841 \return a scalar. Returns 0 if both transverse momenta are zero
842 \ingroup GenVector
843 \see http://doi.org/10.1103/PhysRevLett.121.212301, https://doi.org/10.1103/PhysRevD.99.093013
844*/
845template <class CoordSystem>
848{
849 auto ptp = pp.Pt();
850 auto ptm = pm.Pt();
851 if (ptp == 0 && ptm == 0)
852 return 0;
853 return math_fabs(ptp - ptm) / (ptp + ptm);
854}
855
856#if !defined(ROOT_MATH_SYCL) && !defined(ROOT_MATH_CUDA)
857// ------------- I/O to/from streams -------------
858
859template <class char_t, class traits_t, class Coords>
860inline std::basic_ostream<char_t, traits_t> &
861operator<<(std::basic_ostream<char_t, traits_t> &os, LorentzVector<Coords> const &v)
862{
863 if (!os)
864 return os;
865
866 typename Coords::Scalar a, b, c, d;
867 v.GetCoordinates(a, b, c, d);
868
871 // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
872 } else {
876 }
877
878 return os;
879
880} // op<< <>()
881
882template <class char_t, class traits_t, class Coords>
883inline std::basic_istream<char_t, traits_t> &
884operator>>(std::basic_istream<char_t, traits_t> &is, LorentzVector<Coords> &v)
885{
886 if (!is)
887 return is;
888
889 typename Coords::Scalar a, b, c, d;
890
893 // TODO: call MF's bitwise-accurate functions on each of a, b, c
894 } else {
896 is >> a;
898 is >> b;
900 is >> c;
902 is >> d;
904 }
905
906 if (is)
907 v.SetCoordinates(a, b, c, d);
908 return is;
909
910} // op>> <>()
911#endif
912
913} // end namespace ROOT_MATH_ARCH
914
915} // end namespace ROOT
916
917#include <sstream>
918namespace cling {
919template <typename CoordSystem>
920std::string printValue(const LorentzVector<CoordSystem> *v)
921{
922 std::stringstream s;
923 s << *v;
924 return s.str();
925}
926
927} // end namespace cling
928
929#endif
930
931// #include "MathX/GenVectorX/LorentzVectorOperations.h"
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#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.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t dest
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t src
float * q
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
DisplacementVector3D< Cartesian3D< Scalar > > BetaVector
Scalar Pz() const
spatial Z component
void GetCoordinates(Scalar dest[]) const
get internal data into an array of 4 Scalar numbers
void GetCoordinates(IT begin) const
get internal data into 4 Scalars at *begin
Scalar Px() const
spatial X component
Scalar M2() const
return magnitude (mass) squared M2 = T**2 - X**2 - Y**2 - Z**2 (we use -,-,-,+ metric)
bool operator!=(const LorentzVector &rhs) const
Scalar R() const
return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
Scalar Pt() const
return the transverse spatial component sqrt ( X**2 + Y**2 )
Scalar DeltaR(const OtherLorentzVector &v, const bool useRapidity=false) const
deltaRapidity between this and vector v
LorentzVector & operator/=(Scalar a)
division by a scalar quantity v /= a
Scalar Et() const
return the transverse energy
LorentzVector & operator+=(const OtherLorentzVector &q)
Self addition with another Vector ( v+= q ) Enable the addition with any other LorentzVector.
LorentzVector< CoordSystem > & SetPt(Scalar a)
Scalar Mt2() const
return the transverse mass squared
LorentzVector< CoordSystem > & SetM(Scalar a)
Scalar Phi() const
azimuthal Angle
LorentzVector< CoordSystem > & SetPy(Scalar a)
Scalar Eta() const
pseudorapidity
bool isLightlike(Scalar tolerance=100 *std::numeric_limits< Scalar >::epsilon()) const
Determine if momentum-energy can represent a massless particle.
Scalar E() const
return 4-th component (time, or energy for a 4-momentum vector)
LorentzVector()
default constructor of an empty vector (Px = Py = Pz = E = 0 )
LorentzVector operator+(const OtherLorentzVector &v2) const
addition of two LorentzVectors (v3 = v1 + v2) Enable the addition with any other LorentzVector
LorentzVector & operator-=(const OtherLorentzVector &q)
Self subtraction of another Vector from this ( v-= q ) Enable the addition with any other LorentzVect...
LorentzVector operator-() const
Negative of a LorentzVector (q = - v )
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...
LorentzVector< CoordSystem > & SetEta(Scalar a)
BetaVector BoostToCM(const Other4Vector &v) const
The beta vector for the boost that would bring this vector into its center of mass frame (zero moment...
LorentzVector(const ForeignLorentzVector &v)
Construct from a foreign 4D vector type, for example, HepLorentzVector Precondition: v must implement...
LorentzVector< CoordSystem > & SetPxPyPzE(Scalar xx, Scalar yy, Scalar zz, Scalar ee)
LorentzVector< CoordSystem > & SetPhi(Scalar a)
bool isSpacelike() const
Determine if momentum-energy is spacelike, and represents a tachyon.
Scalar Beta() const
Return beta scalar value.
LorentzVector< CoordSystem > operator/(const Scalar &a) const
Divide a LorentzVector by a scalar quantity.
LorentzVector< CoordSystem > & SetPx(Scalar a)
const CoordSystem & Coordinates() const
Retrieve a const reference to the coordinates object.
void GetCoordinates(Scalar &a, Scalar &b, Scalar &c, Scalar &d) const
get internal data into 4 Scalar numbers
bool isTimelike() const
Determine if momentum-energy can represent a physical massive particle.
LorentzVector< CoordSystem > & SetPz(Scalar a)
LorentzVector< CoordSystem > & SetE(Scalar a)
Methods setting a Single-component Work only if the component is one of which the vector is represent...
Scalar Gamma() const
Return Gamma scalar value.
Scalar Perp2() const
return the square of the transverse spatial component ( X**2 + Y**2 )
LorentzVector< CoordSystem > & SetCoordinates(IT begin, IT end)
Set internal data based on 4 Scalars at *begin to *end.
LorentzVector(const LorentzVector< Coords > &v)
constructor from a LorentzVector expressed in different coordinates, or using a different Scalar type
LorentzVector operator*(const Scalar &a) const
product of a LorentzVector by a scalar quantity
Scalar M() const
return magnitude (mass) using the (-,-,-,+) metric.
LorentzVector< CoordSystem > & SetCoordinates(Scalar a, Scalar b, Scalar c, Scalar d)
Set internal data based on 4 Scalar numbers.
bool operator==(const LorentzVector &rhs) const
Exact equality.
Scalar Dot(const OtherLorentzVector &q) const
scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+) Enable the product using any ot...
Scalar Py() const
spatial Y component
LorentzVector & operator=(const LorentzVector< OtherCoords > &v)
Assignment operator from a lorentz vector of arbitrary type.
LorentzVector operator-(const OtherLorentzVector &v2) const
subtraction of two LorentzVectors (v3 = v1 - v2) Enable the subtraction of any other LorentzVector
BetaVector BoostToCM() const
The beta vector for the boost that would bring this vector into its center of mass frame (zero moment...
LorentzVector(const Scalar &a, const Scalar &b, const Scalar &c, const Scalar &d)
generic constructors from four scalar values.
Scalar P2() const
return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
LorentzVector & operator*=(Scalar a)
multiplication by a scalar quantity v *= a
Scalar Theta() const
polar Angle
Scalar ColinearRapidity() const
Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)].
Scalar Mt() const
return the transverse mass
Scalar Et2() const
return the transverse energy squared
LorentzVector< CoordSystem > & SetCoordinates(const Scalar src[])
Set internal data based on an array of 4 Scalar numbers.
void GetCoordinates(IT begin, IT end) const
get internal data into 4 Scalars at *begin to *end
Scalar Rapidity() const
Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)].
DisplacementVector3D< Cartesian3D< Scalar > > Vect() const
get the spatial components of the Vector in a DisplacementVector based on Cartesian Coordinates
LorentzVector< CoordSystem >::Scalar Acoplanarity(LorentzVector< CoordSystem > const &pp, LorentzVector< CoordSystem > const &pm)
pair (p+ p-) acoplanarity alpha = 1 - |phi+ - phi-|/pi.
LorentzVector< CoordSystem >::Scalar AsymmetryVectorial(LorentzVector< CoordSystem > const &pp, LorentzVector< CoordSystem > const &pm)
pair (p+ p-) vectorial asymmetry Av = |Pt+ - Pt-|/|Pt+ + Pt-|.
LorentzVector< CoordSystem >::Scalar AsymmetryScalar(LorentzVector< CoordSystem > const &pp, LorentzVector< CoordSystem > const &pm)
pair (p+ p-) scalar asymmetry As = ||Pt+| - |Pt-|/||Pt+| + |Pt-||.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
std::basic_istream< char_t, traits_t > & require_delim(std::basic_istream< char_t, traits_t > &is, manip_t m)
char_t get_manip(std::basic_ios< char_t, traits_t > &ios, manip_t m)
Definition GenVectorIO.h:60
void set_manip(std::basic_ios< char_t, traits_t > &ios, manip_t m, char_t ch)
Definition GenVectorIO.h:77
Scalar math_log(Scalar x)
void GenVector_Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
Scalar math_sqrt(Scalar x)
std::ostream & operator<<(std::ostream &os, const AxisAngle &a)
Stream Output and Input.
Definition AxisAngle.cxx:98
std::basic_istream< char_t, traits_t > & operator>>(std::basic_istream< char_t, traits_t > &is, DisplacementVector2D< T, U > &v)
AxisAngle operator*(RotationX const &r1, AxisAngle const &r2)
Multiplication of an axial rotation by an AxisAngle.
Scalar math_fabs(Scalar x)
auto * tt
Definition textangle.C:16