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 >
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 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 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 // individual coordinate accessors in various coordinate systems
267
268 /**
269 spatial X component
270 */
271 Scalar Px() const { return fCoordinates.Px(); }
272 Scalar X() const { return fCoordinates.Px(); }
273 /**
274 spatial Y component
275 */
276 Scalar Py() const { return fCoordinates.Py(); }
277 Scalar Y() const { return fCoordinates.Py(); }
278 /**
279 spatial Z component
280 */
281 Scalar Pz() const { return fCoordinates.Pz(); }
282 Scalar Z() const { return fCoordinates.Pz(); }
283 /**
284 return 4-th component (time, or energy for a 4-momentum vector)
285 */
286 Scalar E() const { return fCoordinates.E(); }
287 Scalar T() const { return fCoordinates.E(); }
288 /**
289 return magnitude (mass) squared M2 = T**2 - X**2 - Y**2 - Z**2
290 (we use -,-,-,+ metric)
291 */
292 Scalar M2() const { return fCoordinates.M2(); }
293 /**
294 return magnitude (mass) using the (-,-,-,+) metric.
295 If M2 is negative (space-like vector) a GenVector_exception
296 is suggested and if continuing, - sqrt( -M2) is returned
297 */
298 Scalar M() const { return fCoordinates.M();}
299 /**
300 return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
301 */
302 Scalar R() const { return fCoordinates.R(); }
303 Scalar P() const { return fCoordinates.R(); }
304 /**
305 return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
306 */
307 Scalar P2() const { return P() * P(); }
308 /**
309 return the square of the transverse spatial component ( X**2 + Y**2 )
310 */
311 Scalar Perp2( ) const { return fCoordinates.Perp2();}
312
313 /**
314 return the transverse spatial component sqrt ( X**2 + Y**2 )
315 */
316 Scalar Pt() const { return fCoordinates.Pt(); }
317 Scalar Rho() const { return fCoordinates.Pt(); }
318
319 /**
320 return the transverse mass squared
321 \f[ m_t^2 = E^2 - p{_z}^2 \f]
322 */
323 Scalar Mt2() const { return fCoordinates.Mt2(); }
324
325 /**
326 return the transverse mass
327 \f[ \sqrt{ m_t^2 = E^2 - p{_z}^2} X sign(E^ - p{_z}^2) \f]
328 */
329 Scalar Mt() const { return fCoordinates.Mt(); }
330
331 /**
332 return the transverse energy squared
333 \f[ e_t = \frac{E^2 p_{\perp}^2 }{ |p|^2 } \f]
334 */
335 Scalar Et2() const { return fCoordinates.Et2(); }
336
337 /**
338 return the transverse energy
339 \f[ e_t = \sqrt{ \frac{E^2 p_{\perp}^2 }{ |p|^2 } } X sign(E) \f]
340 */
341 Scalar Et() const { return fCoordinates.Et(); }
342
343 /**
344 azimuthal Angle
345 */
346 Scalar Phi() const { return fCoordinates.Phi();}
347
348 /**
349 polar Angle
350 */
351 Scalar Theta() const { return fCoordinates.Theta(); }
352
353 /**
354 pseudorapidity
355 \f[ \eta = - \ln { \tan { \frac { \theta} {2} } } \f]
356 */
357 Scalar Eta() const { return fCoordinates.Eta(); }
358
359 /**
360 get the spatial components of the Vector in a
361 DisplacementVector based on Cartesian Coordinates
362 */
364 return ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> >( X(), Y(), Z() );
365 }
366
367 // ------ Operations combining two Lorentz vectors ------
368
369 /**
370 scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
371 Enable the product using any other LorentzVector implementing
372 the x(), y() , y() and t() member functions
373 \param q any LorentzVector implementing the x(), y() , z() and t()
374 member functions
375 \return the result of v.q of type according to the base scalar type of v
376 */
377
378 template< class OtherLorentzVector >
379 Scalar Dot(const OtherLorentzVector & q) const {
380 return t()*q.t() - x()*q.x() - y()*q.y() - z()*q.z();
381 }
382
383 /**
384 Self addition with another Vector ( v+= q )
385 Enable the addition with any other LorentzVector
386 \param q any LorentzVector implementing the x(), y() , z() and t()
387 member functions
388 */
389 template< class OtherLorentzVector >
390 inline LorentzVector & operator += ( const OtherLorentzVector & q)
391 {
392 SetXYZT( x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t() );
393 return *this;
394 }
395
396 /**
397 Self subtraction of another Vector from this ( v-= q )
398 Enable the addition with any other LorentzVector
399 \param q any LorentzVector implementing the x(), y() , z() and t()
400 member functions
401 */
402 template< class OtherLorentzVector >
403 LorentzVector & operator -= ( const OtherLorentzVector & q) {
404 SetXYZT( x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t() );
405 return *this;
406 }
407
408 /**
409 addition of two LorentzVectors (v3 = v1 + v2)
410 Enable the addition with any other LorentzVector
411 \param v2 any LorentzVector implementing the x(), y() , z() and t()
412 member functions
413 \return a new LorentzVector of the same type as v1
414 */
415 template<class OtherLorentzVector>
416 LorentzVector operator + ( const OtherLorentzVector & v2) const
417 {
419 v3 += v2;
420 return v3;
421 }
422
423 /**
424 subtraction of two LorentzVectors (v3 = v1 - v2)
425 Enable the subtraction of any other LorentzVector
426 \param v2 any LorentzVector implementing the x(), y() , z() and t()
427 member functions
428 \return a new LorentzVector of the same type as v1
429 */
430 template<class OtherLorentzVector>
431 LorentzVector operator - ( const OtherLorentzVector & v2) const {
433 v3 -= v2;
434 return v3;
435 }
436
437 //--- scaling operations ------
438
439 /**
440 multiplication by a scalar quantity v *= a
441 */
443 fCoordinates.Scale(a);
444 return *this;
445 }
446
447 /**
448 division by a scalar quantity v /= a
449 */
451 fCoordinates.Scale(1/a);
452 return *this;
453 }
454
455 /**
456 product of a LorentzVector by a scalar quantity
457 \param a scalar quantity of type a
458 \return a new mathcoreLorentzVector q = v * a same type as v
459 */
461 LorentzVector tmp(*this);
462 tmp *= a;
463 return tmp;
464 }
465
466 /**
467 Divide a LorentzVector by a scalar quantity
468 \param a scalar quantity of type a
469 \return a new mathcoreLorentzVector q = v / a same type as v
470 */
473 tmp /= a;
474 return tmp;
475 }
476
477 /**
478 Negative of a LorentzVector (q = - v )
479 \return a new LorentzVector with opposite direction and time
480 */
482 //LorentzVector<CoordinateType> v(*this);
483 //v.Negate();
484 return operator*( Scalar(-1) );
485 }
487 return *this;
488 }
489
490 // ---- Relativistic Properties ----
491
492 /**
493 Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)]
494 */
495 Scalar Rapidity() const {
496 // TODO - It would be good to check that E > Pz and use the Throw()
497 // mechanism or at least load a NAN if not.
498 // We should then move the code to a .cpp file.
499 const Scalar ee = E();
500 const Scalar ppz = Pz();
501 using std::log;
502 return Scalar(0.5) * log((ee + ppz) / (ee - ppz));
503 }
504
505 /**
506 Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
507 */
509 // TODO - It would be good to check that E > P and use the Throw()
510 // mechanism or at least load a NAN if not.
511 const Scalar ee = E();
512 const Scalar pp = P();
513 using std::log;
514 return Scalar(0.5) * log((ee + pp) / (ee - pp));
515 }
516
517 /**
518 Determine if momentum-energy can represent a physical massive particle
519 */
520 bool isTimelike( ) const {
521 Scalar ee = E(); Scalar pp = P(); return ee*ee > pp*pp;
522 }
523
524 /**
525 Determine if momentum-energy can represent a massless particle
526 */
527 bool isLightlike( Scalar tolerance
528 = 100*std::numeric_limits<Scalar>::epsilon() ) const {
529 Scalar ee = E(); Scalar pp = P(); Scalar delta = ee-pp;
530 if ( ee==0 ) return pp==0;
531 return delta*delta < tolerance * ee*ee;
532 }
533
534 /**
535 Determine if momentum-energy is spacelike, and represents a tachyon
536 */
537 bool isSpacelike( ) const {
538 Scalar ee = E(); Scalar pp = P(); return ee*ee < pp*pp;
539 }
540
542
543 /**
544 The beta vector for the boost that would bring this vector into
545 its center of mass frame (zero momentum)
546 */
548 if (E() == 0) {
549 if (P() == 0) {
550 return BetaVector();
551 } else {
552 // TODO - should attempt to Throw with msg about
553 // boostVector computed for LorentzVector with t=0
554 return -Vect()/E();
555 }
556 }
557 if (M2() <= 0) {
558 // TODO - should attempt to Throw with msg about
559 // boostVector computed for a non-timelike LorentzVector
560 }
561 return -Vect()/E();
562 }
563
564 /**
565 The beta vector for the boost that would bring this vector into
566 its center of mass frame (zero momentum)
567 */
568 template <class Other4Vector>
569 BetaVector BoostToCM(const Other4Vector& v ) const {
570 Scalar eSum = E() + v.E();
572 if (eSum == 0) {
573 if (vecSum.Mag2() == 0) {
574 return BetaVector();
575 } else {
576 // TODO - should attempt to Throw with msg about
577 // boostToCM computed for two 4-vectors with combined t=0
578 return BetaVector(vecSum/eSum);
579 }
580 // TODO - should attempt to Throw with msg about
581 // boostToCM computed for two 4-vectors with combined e=0
582 }
583 return BetaVector (vecSum * (-1./eSum));
584 }
585
586 //beta and gamma
587
588 /**
589 Return beta scalar value
590 */
591 Scalar Beta() const {
592 if ( E() == 0 ) {
593 if ( P2() == 0)
594 // to avoid Nan
595 return 0;
596 else {
597 GenVector::Throw ("LorentzVector::Beta() - beta computed for LorentzVector with t = 0. Return an Infinite result");
598 return 1./E();
599 }
600 }
601 if ( M2() <= 0 ) {
602 GenVector::Throw ("LorentzVector::Beta() - beta computed for non-timelike LorentzVector . Result is physically meaningless" );
603 }
604 return P() / E();
605 }
606 /**
607 Return Gamma scalar value
608 */
609 Scalar Gamma() const {
610 const Scalar v2 = P2();
611 const Scalar t2 = E() * E();
612 if (E() == 0) {
613 if ( P2() == 0) {
614 return 1;
615 } else {
616 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for LorentzVector with t = 0. Return a zero result");
617
618 }
619 }
620 if ( t2 < v2 ) {
621 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a spacelike LorentzVector. Imaginary result");
622 return 0;
623 }
624 else if ( t2 == v2 ) {
625 GenVector::Throw ("LorentzVector::Gamma() - gamma computed for a lightlike LorentzVector. Infinite result");
626 }
627 using std::sqrt;
628 return Scalar(1) / sqrt(Scalar(1) - v2 / t2);
629 } /* gamma */
630
631
632 // Method providing limited backward name compatibility with CLHEP ----
633
634 Scalar x() const { return fCoordinates.Px(); }
635 Scalar y() const { return fCoordinates.Py(); }
636 Scalar z() const { return fCoordinates.Pz(); }
637 Scalar t() const { return fCoordinates.E(); }
638 Scalar px() const { return fCoordinates.Px(); }
639 Scalar py() const { return fCoordinates.Py(); }
640 Scalar pz() const { return fCoordinates.Pz(); }
641 Scalar e() const { return fCoordinates.E(); }
642 Scalar r() const { return fCoordinates.R(); }
643 Scalar theta() const { return fCoordinates.Theta(); }
644 Scalar phi() const { return fCoordinates.Phi(); }
645 Scalar rho() const { return fCoordinates.Rho(); }
646 Scalar eta() const { return fCoordinates.Eta(); }
647 Scalar pt() const { return fCoordinates.Pt(); }
648 Scalar perp2() const { return fCoordinates.Perp2(); }
649 Scalar mag2() const { return fCoordinates.M2(); }
650 Scalar mag() const { return fCoordinates.M(); }
651 Scalar mt() const { return fCoordinates.Mt(); }
652 Scalar mt2() const { return fCoordinates.Mt2(); }
653
654
655 // Methods requested by CMS ---
656 Scalar energy() const { return fCoordinates.E(); }
657 Scalar mass() const { return fCoordinates.M(); }
658 Scalar mass2() const { return fCoordinates.M2(); }
659
660
661 /**
662 Methods setting a Single-component
663 Work only if the component is one of which the vector is represented.
664 For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
665 */
666 LorentzVector<CoordSystem>& SetE ( Scalar a ) { fCoordinates.SetE (a); return *this; }
668 LorentzVector<CoordSystem>& SetM ( Scalar a ) { fCoordinates.SetM (a); return *this; }
670 LorentzVector<CoordSystem>& SetPt ( Scalar a ) { fCoordinates.SetPt (a); return *this; }
671 LorentzVector<CoordSystem>& SetPx ( Scalar a ) { fCoordinates.SetPx (a); return *this; }
672 LorentzVector<CoordSystem>& SetPy ( Scalar a ) { fCoordinates.SetPy (a); return *this; }
673 LorentzVector<CoordSystem>& SetPz ( Scalar a ) { fCoordinates.SetPz (a); return *this; }
674
675 private:
676
677 CoordSystem fCoordinates; // internal coordinate system
678
679
680 }; // LorentzVector<>
681
682
683
684 // global methods
685
686 /**
687 Scale of a LorentzVector with a scalar quantity a
688 \param a scalar quantity of type a
689 \param v mathcore::LorentzVector based on any coordinate system
690 \return a new mathcoreLorentzVector q = v * a same type as v
691 */
692 template< class CoordSystem >
694 ( const typename LorentzVector<CoordSystem>::Scalar & a,
697 tmp *= a;
698 return tmp;
699 }
700
701 // ------------- I/O to/from streams -------------
702
703 template< class char_t, class traits_t, class Coords >
704 inline
705 std::basic_ostream<char_t,traits_t> &
706 operator << ( std::basic_ostream<char_t,traits_t> & os
707 , LorentzVector<Coords> const & v
708 )
709 {
710 if( !os ) return os;
711
712 typename Coords::Scalar a, b, c, d;
713 v.GetCoordinates(a, b, c, d);
714
717 // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
718 }
719 else {
720 os << detail::get_manip( os, detail::open ) << a
721 << detail::get_manip( os, detail::sep ) << b
722 << detail::get_manip( os, detail::sep ) << c
723 << detail::get_manip( os, detail::sep ) << d
725 }
726
727 return os;
728
729 } // op<< <>()
730
731
732 template< class char_t, class traits_t, class Coords >
733 inline
734 std::basic_istream<char_t,traits_t> &
735 operator >> ( std::basic_istream<char_t,traits_t> & is
737 )
738 {
739 if( !is ) return is;
740
741 typename Coords::Scalar a, b, c, d;
742
745 // TODO: call MF's bitwise-accurate functions on each of a, b, c
746 }
747 else {
749 detail::require_delim( is, detail::sep ); is >> b;
750 detail::require_delim( is, detail::sep ); is >> c;
751 detail::require_delim( is, detail::sep ); is >> d;
753 }
754
755 if( is )
756 v.SetCoordinates(a, b, c, d);
757 return is;
758
759 } // op>> <>()
760
761
762
763 } // end namespace Math
764
765} // end namespace ROOT
766
767#include <sstream>
768namespace cling
769{
770template<typename CoordSystem>
771std::string printValue(const ROOT::Math::LorentzVector<CoordSystem> *v)
772{
773 std::stringstream s;
774 s << *v;
775 return s.str();
776}
777
778} // end namespace cling
779
780#endif
781
782//#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 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(const ForeignLorentzVector &v)
Construct from a foreign 4D vector type, for example, HepLorentzVector Precondition: v must implement...
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
Scalar Py() const
spatial Y component
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 )
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
LorentzVector(const LorentzVector< Coords > &v)
constructor from a LorentzVector expressed in different coordinates, or using a different Scalar type
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
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
#define dest(otri, vertexptr)
Definition triangle.c:1041