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
22
24
26
27#include <cmath>
28#include <string>
29
30namespace ROOT {
31
32 namespace Math {
33
34//__________________________________________________________________________________________
35/** @ingroup GenVector
36
37Class describing a generic LorentzVector in the 4D space-time,
38using the specified coordinate system for the spatial vector part.
39The metric used for the LorentzVector is (-,-,-,+).
40In the case of LorentzVector we don't distinguish the concepts
41of points and displacement vectors as in the 3D case,
42since the main use case for 4D Vectors is to describe the kinematics of
43relativistic particles. A LorentzVector behaves like a
44DisplacementVector in 4D. The Minkowski components could be viewed as
45v and t, or for kinematic 4-vectors, as p and E.
46
47ROOT provides specialisations and aliases to them of the ROOT::Math::LorentzVector template:
48- ROOT::Math::PtEtaPhiMVector based on pt (rho),eta,phi and M (t) coordinates in double precision
49- ROOT::Math::PtEtaPhiEVector based on pt (rho),eta,phi and E (t) coordinates in double precision
50- ROOT::Math::PxPyPzMVector based on px,py,pz and M (mass) coordinates in double precision
51- ROOT::Math::PxPyPzEVector based on px,py,pz and E (energy) coordinates in double precision
52- ROOT::Math::XYZTVector based on x,y,z,t coordinates (cartesian) in double precision (same as PxPyPzEVector)
53- ROOT::Math::XYZTVectorF based on x,y,z,t coordinates (cartesian) in float precision (same as PxPyPzEVector but float)
54
55@sa Overview of the @ref GenVector "physics vector library"
56*/
57
58 template< class CoordSystem >
60
61 public:
62
63 // ------ ctors ------
64
65 typedef typename CoordSystem::Scalar Scalar;
66 typedef CoordSystem CoordinateType;
67
68 /**
69 default constructor of an empty vector (Px = Py = Pz = E = 0 )
70 */
72
73 /**
74 generic constructors from four scalar values.
75 The association between values and coordinate depends on the
76 coordinate system. For PxPyPzE4D,
77 \param a scalar value (Px)
78 \param b scalar value (Py)
79 \param c scalar value (Pz)
80 \param d scalar value (E)
81 */
83 const Scalar & b,
84 const Scalar & c,
85 const Scalar & d) :
86 fCoordinates(a , b, c, d) { }
87
88 /**
89 constructor from a LorentzVector expressed in different
90 coordinates, or using a different Scalar type
91 */
92 template< class Coords >
93 explicit constexpr LorentzVector(const LorentzVector<Coords> & v ) :
95
96 /**
97 Construct from a foreign 4D vector type, for example, HepLorentzVector
98 Precondition: v must implement methods x(), y(), z(), and t()
99 */
100 template<class ForeignLorentzVector,
101 typename = decltype(std::declval<ForeignLorentzVector>().x()
102 + std::declval<ForeignLorentzVector>().y()
103 + std::declval<ForeignLorentzVector>().z()
104 + std::declval<ForeignLorentzVector>().t())>
105 explicit constexpr LorentzVector( const ForeignLorentzVector & v) :
106 fCoordinates(PxPyPzE4D<Scalar>( v.x(), v.y(), v.z(), v.t() ) ) { }
107
108#ifdef LATER
109 /**
110 construct from a generic linear algebra vector implementing operator []
111 and with a size of at least 4. This could be also a C array
112 In this case v[0] is the first data member
113 ( Px for a PxPyPzE4D base)
114 \param v LA vector
115 \param index0 index of first vector element (Px)
116 */
117 template< class LAVector >
118 explicit constexpr LorentzVector(const LAVector & v, size_t index0 ) {
119 fCoordinates = CoordSystem ( v[index0], v[index0+1], v[index0+2], v[index0+3] );
120 }
121#endif
122
123
124 // ------ assignment ------
125
126 /**
127 Assignment operator from a lorentz vector of arbitrary type
128 */
129 template< class OtherCoords >
132 return *this;
133 }
134
135 /**
136 assignment from any other Lorentz vector implementing
137 x(), y(), z() and t()
138 */
139 template<class ForeignLorentzVector,
140 typename = decltype(std::declval<ForeignLorentzVector>().x()
141 + std::declval<ForeignLorentzVector>().y()
142 + std::declval<ForeignLorentzVector>().z()
143 + std::declval<ForeignLorentzVector>().t())>
144 LorentzVector & operator = ( const ForeignLorentzVector & v) {
145 SetXYZT( v.x(), v.y(), v.z(), v.t() );
146 return *this;
147 }
148
149#ifdef LATER
150 /**
151 assign from a generic linear algebra vector implementing operator []
152 and with a size of at least 4
153 In this case v[0] is the first data member
154 ( Px for a PxPyPzE4D base)
155 \param v LA vector
156 \param index0 index of first vector element (Px)
157 */
158 template< class LAVector >
159 LorentzVector & AssignFrom(const LAVector & v, size_t index0=0 ) {
160 fCoordinates.SetCoordinates( v[index0], v[index0+1], v[index0+2], v[index0+3] );
161 return *this;
162 }
163#endif
164
165 // ------ Set, Get, and access coordinate data ------
166
167 /**
168 Retrieve a const reference to the coordinates object
169 */
170 const CoordSystem & Coordinates() const {
171 return fCoordinates;
172 }
173
174 /**
175 Set internal data based on an array of 4 Scalar numbers
176 */
178 fCoordinates.SetCoordinates(src);
179 return *this;
180 }
181
182 /**
183 Set internal data based on 4 Scalar numbers
184 */
186 fCoordinates.SetCoordinates(a, b, c, d);
187 return *this;
188 }
189
190 /**
191 Set internal data based on 4 Scalars at *begin to *end
192 */
193 template< class IT >
195 IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
196 (void)end;
197 assert (++begin==end);
198 SetCoordinates (*a,*b,*c,*d);
199 return *this;
200 }
201
202 /**
203 get internal data into 4 Scalar numbers
204 */
205 void GetCoordinates( Scalar& a, Scalar& b, Scalar& c, Scalar & d ) const
206 { fCoordinates.GetCoordinates(a, b, c, d); }
207
208 /**
209 get internal data into an array of 4 Scalar numbers
210 */
211 void GetCoordinates( Scalar dest[] ) const
212 { fCoordinates.GetCoordinates(dest); }
213
214 /**
215 get internal data into 4 Scalars at *begin to *end
216 */
217 template <class IT>
218 void GetCoordinates( IT begin, IT end ) const
219 { IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
220 (void)end;
221 assert (++begin==end);
222 GetCoordinates (*a,*b,*c,*d);
223 }
224
225 /**
226 get internal data into 4 Scalars at *begin
227 */
228 template <class IT>
229 void GetCoordinates( IT begin ) const {
230 Scalar a,b,c,d = 0;
231 GetCoordinates (a,b,c,d);
232 *begin++ = a;
233 *begin++ = b;
234 *begin++ = c;
235 *begin = d;
236 }
237
238 /**
239 set the values of the vector from the cartesian components (x,y,z,t)
240 (if the vector is held in another coordinates, like (Pt,eta,phi,m)
241 then (x, y, z, t) are converted to that form)
242 */
244 fCoordinates.SetPxPyPzE(xx,yy,zz,tt);
245 return *this;
246 }
248 fCoordinates.SetPxPyPzE(xx,yy,zz,ee);
249 return *this;
250 }
251
252 // ------------------- Equality -----------------
253
254 /**
255 Exact equality
256 */
257 bool operator==(const LorentzVector & rhs) const {
258 return fCoordinates==rhs.fCoordinates;
259 }
260 bool operator!= (const LorentzVector & rhs) const {
261 return !(operator==(rhs));
262 }
263
264 // ------ Individual element access, in various coordinate systems ------
265
266 /**
267 dimension
268 */
269 unsigned int Dimension() const
270 {
271 return fDimension;
272 };
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 get the spatial components of the Vector in a
369 DisplacementVector based on Cartesian Coordinates
370 */
372 return ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> >( X(), Y(), Z() );
373 }
374
375 // ------ Operations combining two Lorentz vectors ------
376
377 /**
378 scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
379 Enable the product using any other LorentzVector implementing
380 the x(), y() , y() and t() member functions
381 \param q any LorentzVector implementing the x(), y() , z() and t()
382 member functions
383 \return the result of v.q of type according to the base scalar type of v
384 */
385
386 template< class OtherLorentzVector >
387 Scalar Dot(const OtherLorentzVector & q) const {
388 return t()*q.t() - x()*q.x() - y()*q.y() - z()*q.z();
389 }
390
391 /**
392 Self addition with another Vector ( v+= q )
393 Enable the addition with any other LorentzVector
394 \param q any LorentzVector implementing the x(), y() , z() and t()
395 member functions
396 */
397 template< class OtherLorentzVector >
398 inline LorentzVector & operator += ( const OtherLorentzVector & q)
399 {
400 SetXYZT( x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t() );
401 return *this;
402 }
403
404 /**
405 Self subtraction of another Vector from this ( v-= q )
406 Enable the addition with any other LorentzVector
407 \param q any LorentzVector implementing the x(), y() , z() and t()
408 member functions
409 */
410 template< class OtherLorentzVector >
411 LorentzVector & operator -= ( const OtherLorentzVector & q) {
412 SetXYZT( x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t() );
413 return *this;
414 }
415
416 /**
417 addition of two LorentzVectors (v3 = v1 + v2)
418 Enable the addition with any other LorentzVector
419 \param v2 any LorentzVector implementing the x(), y() , z() and t()
420 member functions
421 \return a new LorentzVector of the same type as v1
422 */
423 template<class OtherLorentzVector>
424 LorentzVector operator + ( const OtherLorentzVector & v2) const
425 {
427 v3 += v2;
428 return v3;
429 }
430
431 /**
432 subtraction of two LorentzVectors (v3 = v1 - v2)
433 Enable the subtraction of any other LorentzVector
434 \param v2 any LorentzVector implementing the x(), y() , z() and t()
435 member functions
436 \return a new LorentzVector of the same type as v1
437 */
438 template<class OtherLorentzVector>
439 LorentzVector operator - ( const OtherLorentzVector & v2) const {
441 v3 -= v2;
442 return v3;
443 }
444
445 //--- scaling operations ------
446
447 /**
448 multiplication by a scalar quantity v *= a
449 */
451 fCoordinates.Scale(a);
452 return *this;
453 }
454
455 /**
456 division by a scalar quantity v /= a
457 */
459 fCoordinates.Scale(1/a);
460 return *this;
461 }
462
463 /**
464 product of a LorentzVector by a scalar quantity
465 \param a scalar quantity of type a
466 \return a new mathcoreLorentzVector q = v * a same type as v
467 */
469 LorentzVector tmp(*this);
470 tmp *= a;
471 return tmp;
472 }
473
474 /**
475 Divide a LorentzVector by a scalar quantity
476 \param a scalar quantity of type a
477 \return a new mathcoreLorentzVector q = v / a same type as v
478 */
481 tmp /= a;
482 return tmp;
483 }
484
485 /**
486 Negative of a LorentzVector (q = - v )
487 \return a new LorentzVector with opposite direction and time
488 */
490 //LorentzVector<CoordinateType> v(*this);
491 //v.Negate();
492 return operator*( Scalar(-1) );
493 }
495 return *this;
496 }
497
498 // ---- Relativistic Properties ----
499
500 /**
501 Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)]
502 */
503 Scalar Rapidity() const {
504 // TODO - It would be good to check that E > Pz and use the Throw()
505 // mechanism or at least load a NAN if not.
506 // We should then move the code to a .cpp file.
507 const Scalar ee = E();
508 const Scalar ppz = Pz();
509 using std::log;
510 return Scalar(0.5) * log((ee + ppz) / (ee - ppz));
511 }
512
513 /**
514 Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
515 */
517 // TODO - It would be good to check that E > P and use the Throw()
518 // mechanism or at least load a NAN if not.
519 const Scalar ee = E();
520 const Scalar pp = P();
521 using std::log;
522 return Scalar(0.5) * log((ee + pp) / (ee - pp));
523 }
524
525 /**
526 Determine if momentum-energy can represent a physical massive particle
527 */
528 bool isTimelike( ) const {
529 Scalar ee = E(); Scalar pp = P(); return ee*ee > pp*pp;
530 }
531
532 /**
533 Determine if momentum-energy can represent a massless particle
534 */
535 bool isLightlike( Scalar tolerance
536 = 100*std::numeric_limits<Scalar>::epsilon() ) const {
537 Scalar ee = E(); Scalar pp = P(); Scalar delta = ee-pp;
538 if ( ee==0 ) return pp==0;
539 return delta*delta < tolerance * ee*ee;
540 }
541
542 /**
543 Determine if momentum-energy is spacelike, and represents a tachyon
544 */
545 bool isSpacelike( ) const {
546 Scalar ee = E(); Scalar pp = P(); return ee*ee < pp*pp;
547 }
548
550
551 /**
552 The beta vector for the boost that would bring this vector into
553 its center of mass frame (zero momentum)
554 */
556 if (E() == 0) {
557 if (P() == 0) {
558 return BetaVector();
559 } else {
560 // TODO - should attempt to Throw with msg about
561 // boostVector computed for LorentzVector with t=0
562 return -Vect()/E();
563 }
564 }
565 if (M2() <= 0) {
566 // TODO - should attempt to Throw with msg about
567 // boostVector computed for a non-timelike LorentzVector
568 }
569 return -Vect()/E();
570 }
571
572 /**
573 The beta vector for the boost that would bring this vector into
574 its center of mass frame (zero momentum)
575 */
576 template <class Other4Vector>
577 BetaVector BoostToCM(const Other4Vector& v ) const {
578 Scalar eSum = E() + v.E();
580 if (eSum == 0) {
581 if (vecSum.Mag2() == 0) {
582 return BetaVector();
583 } else {
584 // TODO - should attempt to Throw with msg about
585 // boostToCM computed for two 4-vectors with combined t=0
586 return BetaVector(vecSum/eSum);
587 }
588 // TODO - should attempt to Throw with msg about
589 // boostToCM computed for two 4-vectors with combined e=0
590 }
591 return BetaVector (vecSum * (-1./eSum));
592 }
593
594 //beta and gamma
595
596 /**
597 Return beta scalar value
598 */
599 Scalar Beta() const {
600 if ( E() == 0 ) {
601 if ( P2() == 0)
602 // to avoid Nan
603 return 0;
604 else {
605 GenVector::Throw ("LorentzVector::Beta() - beta computed for LorentzVector with t = 0. Return an Infinite result");
606 return 1./E();
607 }
608 }
609 if ( M2() <= 0 ) {
610 GenVector::Throw ("LorentzVector::Beta() - beta computed for non-timelike LorentzVector . Result is physically meaningless" );
611 }
612 return P() / E();
613 }
614 /**
615 Return Gamma scalar value
616 */
617 Scalar Gamma() const {
618 const Scalar v2 = P2();
619 const Scalar t2 = E() * E();
620 if (E() == 0) {
621 if ( P2() == 0) {
622 return 1;
623 } else {
624 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for LorentzVector with t = 0. Return a zero result");
625
626 }
627 }
628 if ( t2 < v2 ) {
629 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a spacelike LorentzVector. Imaginary result");
630 return 0;
631 }
632 else if ( t2 == v2 ) {
633 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a lightlike LorentzVector. Infinite result");
634 }
635 using std::sqrt;
636 return Scalar(1) / sqrt(Scalar(1) - v2 / t2);
637 } /* gamma */
638
639
640 // Method providing limited backward name compatibility with CLHEP ----
641
642 Scalar x() const { return fCoordinates.Px(); }
643 Scalar y() const { return fCoordinates.Py(); }
644 Scalar z() const { return fCoordinates.Pz(); }
645 Scalar t() const { return fCoordinates.E(); }
646 Scalar px() const { return fCoordinates.Px(); }
647 Scalar py() const { return fCoordinates.Py(); }
648 Scalar pz() const { return fCoordinates.Pz(); }
649 Scalar e() const { return fCoordinates.E(); }
650 Scalar r() const { return fCoordinates.R(); }
651 Scalar theta() const { return fCoordinates.Theta(); }
652 Scalar phi() const { return fCoordinates.Phi(); }
653 Scalar rho() const { return fCoordinates.Rho(); }
654 Scalar eta() const { return fCoordinates.Eta(); }
655 Scalar pt() const { return fCoordinates.Pt(); }
656 Scalar perp2() const { return fCoordinates.Perp2(); }
657 Scalar mag2() const { return fCoordinates.M2(); }
658 Scalar mag() const { return fCoordinates.M(); }
659 Scalar mt() const { return fCoordinates.Mt(); }
660 Scalar mt2() const { return fCoordinates.Mt2(); }
661
662
663 // Methods requested by CMS ---
664 Scalar energy() const { return fCoordinates.E(); }
665 Scalar mass() const { return fCoordinates.M(); }
666 Scalar mass2() const { return fCoordinates.M2(); }
667
668
669 /**
670 Methods setting a Single-component
671 Work only if the component is one of which the vector is represented.
672 For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
673 */
674 LorentzVector<CoordSystem>& SetE ( Scalar a ) { fCoordinates.SetE (a); return *this; }
676 LorentzVector<CoordSystem>& SetM ( Scalar a ) { fCoordinates.SetM (a); return *this; }
678 LorentzVector<CoordSystem>& SetPt ( Scalar a ) { fCoordinates.SetPt (a); return *this; }
679 LorentzVector<CoordSystem>& SetPx ( Scalar a ) { fCoordinates.SetPx (a); return *this; }
680 LorentzVector<CoordSystem>& SetPy ( Scalar a ) { fCoordinates.SetPy (a); return *this; }
681 LorentzVector<CoordSystem>& SetPz ( Scalar a ) { fCoordinates.SetPz (a); return *this; }
682
683 private:
684
685 CoordSystem fCoordinates; // internal coordinate system
686 static constexpr unsigned int fDimension = CoordinateType::Dimension;
687
688 }; // LorentzVector<>
689
690
691
692 // global methods
693
694 /**
695 Scale of a LorentzVector with a scalar quantity a
696 \param a scalar quantity of type a
697 \param v mathcore::LorentzVector based on any coordinate system
698 \return a new mathcoreLorentzVector q = v * a same type as v
699 */
700 template< class CoordSystem >
702 ( const typename LorentzVector<CoordSystem>::Scalar & a,
705 tmp *= a;
706 return tmp;
707 }
708
709 // ------------- I/O to/from streams -------------
710
711 template< class char_t, class traits_t, class Coords >
712 inline
713 std::basic_ostream<char_t,traits_t> &
714 operator << ( std::basic_ostream<char_t,traits_t> & os
715 , LorentzVector<Coords> const & v
716 )
717 {
718 if( !os ) return os;
719
720 typename Coords::Scalar a, b, c, d;
721 v.GetCoordinates(a, b, c, d);
722
725 // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
726 }
727 else {
728 os << detail::get_manip( os, detail::open ) << a
729 << detail::get_manip( os, detail::sep ) << b
730 << detail::get_manip( os, detail::sep ) << c
731 << detail::get_manip( os, detail::sep ) << d
733 }
734
735 return os;
736
737 } // op<< <>()
738
739
740 template< class char_t, class traits_t, class Coords >
741 inline
742 std::basic_istream<char_t,traits_t> &
743 operator >> ( std::basic_istream<char_t,traits_t> & is
745 )
746 {
747 if( !is ) return is;
748
749 typename Coords::Scalar a, b, c, d;
750
753 // TODO: call MF's bitwise-accurate functions on each of a, b, c
754 }
755 else {
757 detail::require_delim( is, detail::sep ); is >> b;
758 detail::require_delim( is, detail::sep ); is >> c;
759 detail::require_delim( is, detail::sep ); is >> d;
761 }
762
763 if( is )
764 v.SetCoordinates(a, b, c, d);
765 return is;
766
767 } // op>> <>()
768
769
770
771 } // end namespace Math
772
773} // end namespace ROOT
774
775#include <sstream>
776namespace cling
777{
778template<typename CoordSystem>
779std::string printValue(const ROOT::Math::LorentzVector<CoordSystem> *v)
780{
781 std::stringstream s;
782 s << *v;
783 return s.str();
784}
785
786} // end namespace cling
787
788#endif
789
790//#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
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.
Scalar Mag2() const
Magnitute squared ( r^2 in spherical coordinate)
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 )
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)
Class describing a 4D cartesian coordinate system (x, y, z, t coordinates) or momentum-energy vectors...
Definition PxPyPzE4D.h:44
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Namespace for new Math classes and functions.
void Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
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::basic_istream< char_t, traits_t > & operator>>(std::basic_istream< char_t, traits_t > &is, DisplacementVector2D< T, U > &v)
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
auto * tt
Definition textangle.C:16