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