Logo ROOT   6.18/05
Reference Guide
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
29namespace ROOT {
30
31 namespace Math {
32
33//__________________________________________________________________________________________
34 /**
35 Class describing a generic LorentzVector in the 4D space-time,
36 using the specified coordinate system for the spatial vector part.
37 The metric used for the LorentzVector is (-,-,-,+).
38 In the case of LorentzVector we don't distinguish the concepts
39 of points and displacement vectors as in the 3D case,
40 since the main use case for 4D Vectors is to describe the kinematics of
41 relativistic particles. A LorentzVector behaves like a
42 DisplacementVector in 4D. The Minkowski components could be viewed as
43 v and t, or for kinematic 4-vectors, as p and E.
44
45 ROOT provides specialisations and aliases to them of the ROOT::Math::LorentzVector template:
46 - ROOT::Math::PtEtaPhiMVector based on pt (rho),eta,phi and M (t) coordinates in double precision
47 - ROOT::Math::PtEtaPhiEVector based on pt (rho),eta,phi and E (t) coordinates in double precision
48 - ROOT::Math::PxPyPzMVector based on px,py,pz and M (mass) coordinates in double precision
49 - ROOT::Math::XYZTVector based on x,y,z,t coordinates (cartesian) in double precision
50 - ROOT::Math::XYZTVectorF based on x,y,z,t coordinates (cartesian) in float precision
51
52More details about the GenVector package can be found [here](Vector.html).
53
54
55 @ingroup GenVector
56 */
57 template< class CoordSystem >
59
60 public:
61
62 // ------ ctors ------
63
64 typedef typename CoordSystem::Scalar Scalar;
65 typedef CoordSystem CoordinateType;
66
67 /**
68 default constructor of an empty vector (Px = Py = Pz = E = 0 )
69 */
71
72 /**
73 generic constructors from four scalar values.
74 The association between values and coordinate depends on the
75 coordinate system. For PxPyPzE4D,
76 \param a scalar value (Px)
77 \param b scalar value (Py)
78 \param c scalar value (Pz)
79 \param d scalar value (E)
80 */
82 const Scalar & b,
83 const Scalar & c,
84 const Scalar & d) :
85 fCoordinates(a , b, c, d) { }
86
87 /**
88 constructor from a LorentzVector expressed in different
89 coordinates, or using a different Scalar type
90 */
91 template< class Coords >
94
95 /**
96 Construct from a foreign 4D vector type, for example, HepLorentzVector
97 Precondition: v must implement methods x(), y(), z(), and t()
98 */
99 template<class ForeignLorentzVector>
100 explicit LorentzVector( const ForeignLorentzVector & v) :
101 fCoordinates(PxPyPzE4D<Scalar>( v.x(), v.y(), v.z(), v.t() ) ) { }
102
103#ifdef LATER
104 /**
105 construct from a generic linear algebra vector implementing operator []
106 and with a size of at least 4. This could be also a C array
107 In this case v[0] is the first data member
108 ( Px for a PxPyPzE4D base)
109 \param v LA vector
110 \param index0 index of first vector element (Px)
111 */
112 template< class LAVector >
113 explicit LorentzVector(const LAVector & v, size_t index0 ) {
114 fCoordinates = CoordSystem ( v[index0], v[index0+1], v[index0+2], v[index0+3] );
115 }
116#endif
117
118
119 // ------ assignment ------
120
121 /**
122 Assignment operator from a lorentz vector of arbitrary type
123 */
124 template< class OtherCoords >
126 fCoordinates = v.Coordinates();
127 return *this;
128 }
129
130 /**
131 assignment from any other Lorentz vector implementing
132 x(), y(), z() and t()
133 */
134 template<class ForeignLorentzVector>
135 LorentzVector & operator = ( const ForeignLorentzVector & v) {
136 SetXYZT( v.x(), v.y(), v.z(), v.t() );
137 return *this;
138 }
139
140#ifdef LATER
141 /**
142 assign from a generic linear algebra vector implementing operator []
143 and with a size of at least 4
144 In this case v[0] is the first data member
145 ( Px for a PxPyPzE4D base)
146 \param v LA vector
147 \param index0 index of first vector element (Px)
148 */
149 template< class LAVector >
150 LorentzVector & AssignFrom(const LAVector & v, size_t index0=0 ) {
151 fCoordinates.SetCoordinates( v[index0], v[index0+1], v[index0+2], v[index0+3] );
152 return *this;
153 }
154#endif
155
156 // ------ Set, Get, and access coordinate data ------
157
158 /**
159 Retrieve a const reference to the coordinates object
160 */
161 const CoordSystem & Coordinates() const {
162 return fCoordinates;
163 }
164
165 /**
166 Set internal data based on an array of 4 Scalar numbers
167 */
169 fCoordinates.SetCoordinates(src);
170 return *this;
171 }
172
173 /**
174 Set internal data based on 4 Scalar numbers
175 */
177 fCoordinates.SetCoordinates(a, b, c, d);
178 return *this;
179 }
180
181 /**
182 Set internal data based on 4 Scalars at *begin to *end
183 */
184//#ifdef NDEBUG
185 //this does not compile in CINT
186// template< class IT >
187// LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */ ) {
188// #endif
189 template< class IT >
190#ifndef NDEBUG
192#else
193 LorentzVector<CoordSystem>& SetCoordinates( IT begin, IT /* end */ ) {
194#endif
195 IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
196 assert (++begin==end);
197 SetCoordinates (*a,*b,*c,*d);
198 return *this;
199 }
200
201 /**
202 get internal data into 4 Scalar numbers
203 */
204 void GetCoordinates( Scalar& a, Scalar& b, Scalar& c, Scalar & d ) const
205 { fCoordinates.GetCoordinates(a, b, c, d); }
206
207 /**
208 get internal data into an array of 4 Scalar numbers
209 */
210 void GetCoordinates( Scalar dest[] ) const
211 { fCoordinates.GetCoordinates(dest); }
212
213 /**
214 get internal data into 4 Scalars at *begin to *end
215 */
216 template <class IT>
217#ifndef NDEBUG
218 void GetCoordinates( IT begin, IT end ) const
219#else
220 void GetCoordinates( IT begin, IT /* end */ ) const
221#endif
222 { IT a = begin; IT b = ++begin; IT c = ++begin; IT d = ++begin;
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 // individual coordinate accessors in various coordinate systems
269
270 /**
271 spatial X component
272 */
273 Scalar Px() const { return fCoordinates.Px(); }
274 Scalar X() const { return fCoordinates.Px(); }
275 /**
276 spatial Y component
277 */
278 Scalar Py() const { return fCoordinates.Py(); }
279 Scalar Y() const { return fCoordinates.Py(); }
280 /**
281 spatial Z component
282 */
283 Scalar Pz() const { return fCoordinates.Pz(); }
284 Scalar Z() const { return fCoordinates.Pz(); }
285 /**
286 return 4-th component (time, or energy for a 4-momentum vector)
287 */
288 Scalar E() const { return fCoordinates.E(); }
289 Scalar T() const { return fCoordinates.E(); }
290 /**
291 return magnitude (mass) squared M2 = T**2 - X**2 - Y**2 - Z**2
292 (we use -,-,-,+ metric)
293 */
294 Scalar M2() const { return fCoordinates.M2(); }
295 /**
296 return magnitude (mass) using the (-,-,-,+) metric.
297 If M2 is negative (space-like vector) a GenVector_exception
298 is suggested and if continuing, - sqrt( -M2) is returned
299 */
300 Scalar M() const { return fCoordinates.M();}
301 /**
302 return the spatial (3D) magnitude ( sqrt(X**2 + Y**2 + Z**2) )
303 */
304 Scalar R() const { return fCoordinates.R(); }
305 Scalar P() const { return fCoordinates.R(); }
306 /**
307 return the square of the spatial (3D) magnitude ( X**2 + Y**2 + Z**2 )
308 */
309 Scalar P2() const { return P() * P(); }
310 /**
311 return the square of the transverse spatial component ( X**2 + Y**2 )
312 */
313 Scalar Perp2( ) const { return fCoordinates.Perp2();}
314
315 /**
316 return the transverse spatial component sqrt ( X**2 + Y**2 )
317 */
318 Scalar Pt() const { return fCoordinates.Pt(); }
319 Scalar Rho() const { return fCoordinates.Pt(); }
320
321 /**
322 return the transverse mass squared
323 \f[ m_t^2 = E^2 - p{_z}^2 \f]
324 */
325 Scalar Mt2() const { return fCoordinates.Mt2(); }
326
327 /**
328 return the transverse mass
329 \f[ \sqrt{ m_t^2 = E^2 - p{_z}^2} X sign(E^ - p{_z}^2) \f]
330 */
331 Scalar Mt() const { return fCoordinates.Mt(); }
332
333 /**
334 return the transverse energy squared
335 \f[ e_t = \frac{E^2 p_{\perp}^2 }{ |p|^2 } \f]
336 */
337 Scalar Et2() const { return fCoordinates.Et2(); }
338
339 /**
340 return the transverse energy
341 \f[ e_t = \sqrt{ \frac{E^2 p_{\perp}^2 }{ |p|^2 } } X sign(E) \f]
342 */
343 Scalar Et() const { return fCoordinates.Et(); }
344
345 /**
346 azimuthal Angle
347 */
348 Scalar Phi() const { return fCoordinates.Phi();}
349
350 /**
351 polar Angle
352 */
353 Scalar Theta() const { return fCoordinates.Theta(); }
354
355 /**
356 pseudorapidity
357 \f[ \eta = - \ln { \tan { \frac { \theta} {2} } } \f]
358 */
359 Scalar Eta() const { return fCoordinates.Eta(); }
360
361 /**
362 get the spatial components of the Vector in a
363 DisplacementVector based on Cartesian Coordinates
364 */
366 return ::ROOT::Math::DisplacementVector3D<Cartesian3D<Scalar> >( X(), Y(), Z() );
367 }
368
369 // ------ Operations combining two Lorentz vectors ------
370
371 /**
372 scalar (Dot) product of two LorentzVector vectors (metric is -,-,-,+)
373 Enable the product using any other LorentzVector implementing
374 the x(), y() , y() and t() member functions
375 \param q any LorentzVector implementing the x(), y() , z() and t()
376 member functions
377 \return the result of v.q of type according to the base scalar type of v
378 */
379
380 template< class OtherLorentzVector >
381 Scalar Dot(const OtherLorentzVector & q) const {
382 return t()*q.t() - x()*q.x() - y()*q.y() - z()*q.z();
383 }
384
385 /**
386 Self addition with another Vector ( v+= q )
387 Enable the addition with any other LorentzVector
388 \param q any LorentzVector implementing the x(), y() , z() and t()
389 member functions
390 */
391 template< class OtherLorentzVector >
392 inline LorentzVector & operator += ( const OtherLorentzVector & q)
393 {
394 SetXYZT( x() + q.x(), y() + q.y(), z() + q.z(), t() + q.t() );
395 return *this;
396 }
397
398 /**
399 Self subtraction of another Vector from this ( v-= q )
400 Enable the addition with any other LorentzVector
401 \param q any LorentzVector implementing the x(), y() , z() and t()
402 member functions
403 */
404 template< class OtherLorentzVector >
405 LorentzVector & operator -= ( const OtherLorentzVector & q) {
406 SetXYZT( x() - q.x(), y() - q.y(), z() - q.z(), t() - q.t() );
407 return *this;
408 }
409
410 /**
411 addition of two LorentzVectors (v3 = v1 + v2)
412 Enable the addition with any other LorentzVector
413 \param v2 any LorentzVector implementing the x(), y() , z() and t()
414 member functions
415 \return a new LorentzVector of the same type as v1
416 */
417 template<class OtherLorentzVector>
418 LorentzVector operator + ( const OtherLorentzVector & v2) const
419 {
421 v3 += v2;
422 return v3;
423 }
424
425 /**
426 subtraction of two LorentzVectors (v3 = v1 - v2)
427 Enable the subtraction of any other LorentzVector
428 \param v2 any LorentzVector implementing the x(), y() , z() and t()
429 member functions
430 \return a new LorentzVector of the same type as v1
431 */
432 template<class OtherLorentzVector>
433 LorentzVector operator - ( const OtherLorentzVector & v2) const {
435 v3 -= v2;
436 return v3;
437 }
438
439 //--- scaling operations ------
440
441 /**
442 multiplication by a scalar quantity v *= a
443 */
445 fCoordinates.Scale(a);
446 return *this;
447 }
448
449 /**
450 division by a scalar quantity v /= a
451 */
453 fCoordinates.Scale(1/a);
454 return *this;
455 }
456
457 /**
458 product of a LorentzVector by a scalar quantity
459 \param a scalar quantity of type a
460 \return a new mathcoreLorentzVector q = v * a same type as v
461 */
463 LorentzVector tmp(*this);
464 tmp *= a;
465 return tmp;
466 }
467
468 /**
469 Divide a LorentzVector by a scalar quantity
470 \param a scalar quantity of type a
471 \return a new mathcoreLorentzVector q = v / a same type as v
472 */
475 tmp /= a;
476 return tmp;
477 }
478
479 /**
480 Negative of a LorentzVector (q = - v )
481 \return a new LorentzVector with opposite direction and time
482 */
484 //LorentzVector<CoordinateType> v(*this);
485 //v.Negate();
486 return operator*( Scalar(-1) );
487 }
489 return *this;
490 }
491
492 // ---- Relativistic Properties ----
493
494 /**
495 Rapidity relative to the Z axis: .5 log [(E+Pz)/(E-Pz)]
496 */
497 Scalar Rapidity() const {
498 // TODO - It would be good to check that E > Pz and use the Throw()
499 // mechanism or at least load a NAN if not.
500 // We should then move the code to a .cpp file.
501 const Scalar ee = E();
502 const Scalar ppz = Pz();
503 return Scalar(0.5) * log((ee + ppz) / (ee - ppz));
504 }
505
506 /**
507 Rapidity in the direction of travel: atanh (|P|/E)=.5 log[(E+P)/(E-P)]
508 */
510 // TODO - It would be good to check that E > P and use the Throw()
511 // mechanism or at least load a NAN if not.
512 const Scalar ee = E();
513 const Scalar pp = P();
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
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 return Scalar(1) / sqrt(Scalar(1) - v2 / t2);
628 } /* gamma */
629
630
631 // Method providing limited backward name compatibility with CLHEP ----
632
633 Scalar x() const { return fCoordinates.Px(); }
634 Scalar y() const { return fCoordinates.Py(); }
635 Scalar z() const { return fCoordinates.Pz(); }
636 Scalar t() const { return fCoordinates.E(); }
637 Scalar px() const { return fCoordinates.Px(); }
638 Scalar py() const { return fCoordinates.Py(); }
639 Scalar pz() const { return fCoordinates.Pz(); }
640 Scalar e() const { return fCoordinates.E(); }
641 Scalar r() const { return fCoordinates.R(); }
642 Scalar theta() const { return fCoordinates.Theta(); }
643 Scalar phi() const { return fCoordinates.Phi(); }
644 Scalar rho() const { return fCoordinates.Rho(); }
645 Scalar eta() const { return fCoordinates.Eta(); }
646 Scalar pt() const { return fCoordinates.Pt(); }
647 Scalar perp2() const { return fCoordinates.Perp2(); }
648 Scalar mag2() const { return fCoordinates.M2(); }
649 Scalar mag() const { return fCoordinates.M(); }
650 Scalar mt() const { return fCoordinates.Mt(); }
651 Scalar mt2() const { return fCoordinates.Mt2(); }
652
653
654 // Methods requested by CMS ---
655 Scalar energy() const { return fCoordinates.E(); }
656 Scalar mass() const { return fCoordinates.M(); }
657 Scalar mass2() const { return fCoordinates.M2(); }
658
659
660 /**
661 Methods setting a Single-component
662 Work only if the component is one of which the vector is represented.
663 For example SetE will work for a PxPyPzE Vector but not for a PxPyPzM Vector.
664 */
665 LorentzVector<CoordSystem>& SetE ( Scalar a ) { fCoordinates.SetE (a); return *this; }
667 LorentzVector<CoordSystem>& SetM ( Scalar a ) { fCoordinates.SetM (a); return *this; }
669 LorentzVector<CoordSystem>& SetPt ( Scalar a ) { fCoordinates.SetPt (a); return *this; }
670 LorentzVector<CoordSystem>& SetPx ( Scalar a ) { fCoordinates.SetPx (a); return *this; }
671 LorentzVector<CoordSystem>& SetPy ( Scalar a ) { fCoordinates.SetPy (a); return *this; }
672 LorentzVector<CoordSystem>& SetPz ( Scalar a ) { fCoordinates.SetPz (a); return *this; }
673
674 private:
675
676 CoordSystem fCoordinates; // internal coordinate system
677
678
679 }; // LorentzVector<>
680
681
682
683 // global nethods
684
685 /**
686 Scale of a LorentzVector with a scalar quantity a
687 \param a scalar quantity of typpe a
688 \param v mathcore::LorentzVector based on any coordinate system
689 \return a new mathcoreLorentzVector q = v * a same type as v
690 */
691 template< class CoordSystem >
693 ( const typename LorentzVector<CoordSystem>::Scalar & a,
696 tmp *= a;
697 return tmp;
698 }
699
700 // ------------- I/O to/from streams -------------
701
702 template< class char_t, class traits_t, class Coords >
703 inline
704 std::basic_ostream<char_t,traits_t> &
705 operator << ( std::basic_ostream<char_t,traits_t> & os
706 , LorentzVector<Coords> const & v
707 )
708 {
709 if( !os ) return os;
710
711 typename Coords::Scalar a, b, c, d;
712 v.GetCoordinates(a, b, c, d);
713
716 // TODO: call MF's bitwise-accurate functions on each of a, b, c, d
717 }
718 else {
719 os << detail::get_manip( os, detail::open ) << a
720 << detail::get_manip( os, detail::sep ) << b
721 << detail::get_manip( os, detail::sep ) << c
722 << detail::get_manip( os, detail::sep ) << d
724 }
725
726 return os;
727
728 } // op<< <>()
729
730
731 template< class char_t, class traits_t, class Coords >
732 inline
733 std::basic_istream<char_t,traits_t> &
734 operator >> ( std::basic_istream<char_t,traits_t> & is
736 )
737 {
738 if( !is ) return is;
739
740 typename Coords::Scalar a, b, c, d;
741
744 // TODO: call MF's bitwise-accurate functions on each of a, b, c
745 }
746 else {
748 detail::require_delim( is, detail::sep ); is >> b;
749 detail::require_delim( is, detail::sep ); is >> c;
750 detail::require_delim( is, detail::sep ); is >> d;
752 }
753
754 if( is )
755 v.SetCoordinates(a, b, c, d);
756 return is;
757
758 } // op>> <>()
759
760
761
762 } // end namespace Math
763
764} // end namespace ROOT
765
766#include <sstream>
767namespace cling
768{
769template<typename CoordSystem>
770std::string printValue(const ROOT::Math::LorentzVector<CoordSystem> *v)
771{
772 std::stringstream s;
773 s << *v;
774 return s.str();
775}
776
777} // end namespace cling
778
779#endif
780
781//#include "Math/GenVector/LorentzVectorOperations.h"
782
783
784
SVector< double, 2 > v
Definition: Dict.h:5
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define c(i)
Definition: RSha256.hxx:101
float * q
Definition: THbookFile.cxx:87
double log(double)
Class describing a generic displacement vector in 3 dimensions.
Scalar Mag2() const
Magnitute squared ( r^2 in spherical coordinate)
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
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.
Definition: LorentzVector.h:81
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(const ForeignLorentzVector &v)
Construct from a foreign 4D vector type, for example, HepLorentzVector Precondition: v must implement...
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
Definition: LorentzVector.h:64
LorentzVector(const LorentzVector< Coords > &v)
constructor from a LorentzVector expressed in different coordinates, or using a different Scalar type
Definition: LorentzVector.h:92
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 )
Definition: LorentzVector.h:70
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:42
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)
Definition: GenVectorIO.h:113
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
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)
Rotation3D::Scalar Scalar
Namespace for new ROOT classes and functions.
Definition: StringConv.hxx:21
static constexpr double s
auto * tt
Definition: textangle.C:16
auto * a
Definition: textangle.C:12
REAL epsilon
Definition: triangle.c:617
#define dest(otri, vertexptr)
Definition: triangle.c:1040