Logo ROOT  
Reference Guide
TLorentzVector.h
Go to the documentation of this file.
1 // @(#)root/physics:$Id$
2 // Author: Pasha Murat , Peter Malzacher 12/02/99
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #ifndef ROOT_TLorentzVector
13 #define ROOT_TLorentzVector
14 
15 
16 //////////////////////////////////////////////////////////////////////////
17 // //
18 // TLorentzVector //
19 // //
20 // Place holder for real lorentz vector class. //
21 // //
22 //////////////////////////////////////////////////////////////////////////
23 
24 #include "TMath.h"
25 #include "TVector3.h"
26 #include "TRotation.h"
27 
28 
29 class TLorentzRotation;
30 
31 
32 class TLorentzVector : public TObject {
33 
34 private:
35 
36  TVector3 fP; // 3 vector component
37  Double_t fE; // time or energy of (x,y,z,t) or (px,py,pz,e)
38 
39 public:
40 
41  typedef Double_t Scalar; // to be able to use it with the ROOT::Math::VectorUtil functions
42 
43  enum { kX=0, kY=1, kZ=2, kT=3, kNUM_COORDINATES=4, kSIZE=kNUM_COORDINATES };
44  // Safe indexing of the coordinates when using with matrices, arrays, etc.
45 
47 
49  // Constructor giving the components x, y, z, t.
50 
51  TLorentzVector(const Double_t * carray);
52  TLorentzVector(const Float_t * carray);
53  // Constructor from an array, not checked!
54 
55  TLorentzVector(const TVector3 & vector3, Double_t t);
56  // Constructor giving a 3-Vector and a time component.
57 
58  TLorentzVector(const TLorentzVector & lorentzvector);
59  // Copy constructor.
60 
61  virtual ~TLorentzVector(){};
62  // Destructor
63 
64  // inline operator TVector3 () const;
65  // inline operator TVector3 & ();
66  // Conversion (cast) to TVector3.
67 
68  inline Double_t X() const;
69  inline Double_t Y() const;
70  inline Double_t Z() const;
71  inline Double_t T() const;
72  // Get position and time.
73 
74  inline void SetX(Double_t a);
75  inline void SetY(Double_t a);
76  inline void SetZ(Double_t a);
77  inline void SetT(Double_t a);
78  // Set position and time.
79 
80  inline Double_t Px() const;
81  inline Double_t Py() const;
82  inline Double_t Pz() const;
83  inline Double_t P() const;
84  inline Double_t E() const;
85  inline Double_t Energy() const;
86  // Get momentum and energy.
87 
88  inline void SetPx(Double_t a);
89  inline void SetPy(Double_t a);
90  inline void SetPz(Double_t a);
91  inline void SetE(Double_t a);
92  // Set momentum and energy.
93 
94  inline TVector3 Vect() const ;
95  // Get spatial component.
96 
97  inline void SetVect(const TVector3 & vect3);
98  // Set spatial component.
99 
100  inline Double_t Theta() const;
101  inline Double_t CosTheta() const;
102  inline Double_t Phi() const; //returns phi from -pi to pi
103  inline Double_t Rho() const;
104  // Get spatial vector components in spherical coordinate system.
105 
106  inline void SetTheta(Double_t theta);
107  inline void SetPhi(Double_t phi);
108  inline void SetRho(Double_t rho);
109  // Set spatial vector components in spherical coordinate system.
110 
111  inline void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e);
112  inline void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t);
113  inline void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m);
114  inline void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m);
115  inline void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e);
116  // Setters to provide the functionality (but a more meanigful name) of
117  // the previous version eg SetV4... PsetV4...
118 
119  inline void GetXYZT(Double_t *carray) const;
120  inline void GetXYZT(Float_t *carray) const;
121  // Getters into an arry
122  // no checking!
123 
124  Double_t operator () (int i) const;
125  inline Double_t operator [] (int i) const;
126  // Get components by index.
127 
128  Double_t & operator () (int i);
129  inline Double_t & operator [] (int i);
130  // Set components by index.
131 
132  inline TLorentzVector & operator = (const TLorentzVector &);
133  // Assignment.
134 
135  inline TLorentzVector operator + (const TLorentzVector &) const;
136  inline TLorentzVector & operator += (const TLorentzVector &);
137  // Additions.
138 
139  inline TLorentzVector operator - (const TLorentzVector &) const;
140  inline TLorentzVector & operator -= (const TLorentzVector &);
141  // Subtractions.
142 
143  inline TLorentzVector operator - () const;
144  // Unary minus.
145 
146  inline TLorentzVector operator * (Double_t a) const;
148  // Scaling with real numbers.
149 
150  inline Bool_t operator == (const TLorentzVector &) const;
151  inline Bool_t operator != (const TLorentzVector &) const;
152  // Comparisons.
153 
154  inline Double_t Perp2() const;
155  // Transverse component of the spatial vector squared.
156 
157  inline Double_t Pt() const;
158  inline Double_t Perp() const;
159  // Transverse component of the spatial vector (R in cylindrical system).
160 
161  inline void SetPerp(Double_t);
162  // Set the transverse component of the spatial vector.
163 
164  inline Double_t Perp2(const TVector3 & v) const;
165  // Transverse component of the spatial vector w.r.t. given axis squared.
166 
167  inline Double_t Pt(const TVector3 & v) const;
168  inline Double_t Perp(const TVector3 & v) const;
169  // Transverse component of the spatial vector w.r.t. given axis.
170 
171  inline Double_t Et2() const;
172  // Transverse energy squared.
173 
174  inline Double_t Et() const;
175  // Transverse energy.
176 
177  inline Double_t Et2(const TVector3 &) const;
178  // Transverse energy w.r.t. given axis squared.
179 
180  inline Double_t Et(const TVector3 &) const;
181  // Transverse energy w.r.t. given axis.
182 
183  inline Double_t DeltaPhi(const TLorentzVector &) const;
184  inline Double_t DeltaR(const TLorentzVector &, Bool_t useRapidity=kFALSE) const;
185  inline Double_t DrEtaPhi(const TLorentzVector &) const;
186  inline Double_t DrRapidityPhi(const TLorentzVector &) const;
187  inline TVector2 EtaPhiVector();
188 
189  inline Double_t Angle(const TVector3 & v) const;
190  // Angle wrt. another vector.
191 
192  inline Double_t Mag2() const;
193  inline Double_t M2() const;
194  // Invariant mass squared.
195 
196  inline Double_t Mag() const;
197  inline Double_t M() const;
198  // Invariant mass. If mag2() is negative then -sqrt(-mag2()) is returned.
199 
200  inline Double_t Mt2() const;
201  // Transverse mass squared.
202 
203  inline Double_t Mt() const;
204  // Transverse mass.
205 
206  inline Double_t Beta() const;
207  inline Double_t Gamma() const;
208 
209  inline Double_t Dot(const TLorentzVector &) const;
210  inline Double_t operator * (const TLorentzVector &) const;
211  // Scalar product.
212 
213  inline void SetVectMag(const TVector3 & spatial, Double_t magnitude);
214  inline void SetVectM(const TVector3 & spatial, Double_t mass);
215  // Copy spatial coordinates, and set energy = sqrt(mass^2 + spatial^2)
216 
217  inline Double_t Plus() const;
218  inline Double_t Minus() const;
219  // Returns t +/- z.
220  // Related to the positive/negative light-cone component,
221  // which some define this way and others define as (t +/- z)/sqrt(2)
222 
223  inline TVector3 BoostVector() const ;
224  // Returns the spatial components divided by the time component.
225 
227  inline void Boost(const TVector3 &);
228  // Lorentz boost.
229 
230  Double_t Rapidity() const;
231  // Returns the rapidity, i.e. 0.5*ln((E+pz)/(E-pz))
232 
233  inline Double_t Eta() const;
234  inline Double_t PseudoRapidity() const;
235  // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
236 
237  inline void RotateX(Double_t angle);
238  // Rotate the spatial component around the x-axis.
239 
240  inline void RotateY(Double_t angle);
241  // Rotate the spatial component around the y-axis.
242 
243  inline void RotateZ(Double_t angle);
244  // Rotate the spatial component around the z-axis.
245 
246  inline void RotateUz(const TVector3 & newUzVector);
247  // Rotates the reference frame from Uz to newUz (unit vector).
248 
249  inline void Rotate(Double_t, const TVector3 &);
250  // Rotate the spatial component around specified axis.
251 
252  inline TLorentzVector & operator *= (const TRotation &);
253  inline TLorentzVector & Transform(const TRotation &);
254  // Transformation with HepRotation.
255 
258  // Transformation with HepLorenzRotation.
259 
260  virtual void Print(Option_t *option="") const;
261 
262  ClassDef(TLorentzVector,4) // A four vector with (-,-,-,+) metric
263 };
264 
265 
266 //inline TLorentzVector operator * (const TLorentzVector &, Double_t a);
267 // moved to TLorentzVector::operator * (Double_t a)
269 // Scaling LorentzVector with a real number
270 
271 
272 inline Double_t TLorentzVector::X() const { return fP.X(); }
273 inline Double_t TLorentzVector::Y() const { return fP.Y(); }
274 inline Double_t TLorentzVector::Z() const { return fP.Z(); }
275 inline Double_t TLorentzVector::T() const { return fE; }
276 
280 inline void TLorentzVector::SetT(Double_t a) { fE = a; }
281 
282 inline Double_t TLorentzVector::Px() const { return X(); }
283 inline Double_t TLorentzVector::Py() const { return Y(); }
284 inline Double_t TLorentzVector::Pz() const { return Z(); }
285 inline Double_t TLorentzVector::P() const { return fP.Mag(); }
286 inline Double_t TLorentzVector::E() const { return T(); }
287 inline Double_t TLorentzVector::Energy() const { return T(); }
288 
293 
294 inline TVector3 TLorentzVector::Vect() const { return fP; }
295 
296 inline void TLorentzVector::SetVect(const TVector3 &p) { fP = p; }
297 
299  return fP.Phi();
300 }
301 
303  return fP.Theta();
304 }
305 
307  return fP.CosTheta();
308 }
309 
310 
312  return fP.Mag();
313 }
314 
316  fP.SetTheta(th);
317 }
318 
320  fP.SetPhi(phi);
321 }
322 
324  fP.SetMag(rho);
325 }
326 
328  fP.SetXYZ(x, y, z);
329  SetT(t);
330 }
331 
333  SetXYZT(px, py, pz, e);
334 }
335 
337  if ( m >= 0 )
338  SetXYZT( x, y, z, TMath::Sqrt(x*x+y*y+z*z+m*m) );
339  else
340  SetXYZT( x, y, z, TMath::Sqrt( TMath::Max((x*x+y*y+z*z-m*m), 0. ) ) );
341 }
342 
344  pt = TMath::Abs(pt);
345  SetXYZM(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,m);
346 }
347 
349  pt = TMath::Abs(pt);
350  SetXYZT(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,e);
351 }
352 
353 inline void TLorentzVector::GetXYZT(Double_t *carray) const {
354  fP.GetXYZ(carray);
355  carray[3] = fE;
356 }
357 
358 inline void TLorentzVector::GetXYZT(Float_t *carray) const{
359  fP.GetXYZ(carray);
360  carray[3] = fE;
361 }
362 
363 inline Double_t & TLorentzVector::operator [] (int i) { return (*this)(i); }
364 inline Double_t TLorentzVector::operator [] (int i) const { return (*this)(i); }
365 
367  fP = q.Vect();
368  fE = q.T();
369  return *this;
370 }
371 
373  return TLorentzVector(fP+q.Vect(), fE+q.T());
374 }
375 
377  fP += q.Vect();
378  fE += q.T();
379  return *this;
380 }
381 
383  return TLorentzVector(fP-q.Vect(), fE-q.T());
384 }
385 
387  fP -= q.Vect();
388  fE -= q.T();
389  return *this;
390 }
391 
393  return TLorentzVector(-X(), -Y(), -Z(), -T());
394 }
395 
397  fP *= a;
398  fE *= a;
399  return *this;
400 }
401 
403  return TLorentzVector(a*X(), a*Y(), a*Z(), a*T());
404 }
405 
407  return (Vect() == q.Vect() && T() == q.T());
408 }
409 
411  return (Vect() != q.Vect() || T() != q.T());
412 }
413 
414 inline Double_t TLorentzVector::Perp2() const { return fP.Perp2(); }
415 
416 inline Double_t TLorentzVector::Perp() const { return fP.Perp(); }
417 
418 inline Double_t TLorentzVector::Pt() const { return Perp(); }
419 
421  fP.SetPerp(r);
422 }
423 
424 inline Double_t TLorentzVector::Perp2(const TVector3 &v) const {
425  return fP.Perp2(v);
426 }
427 
428 inline Double_t TLorentzVector::Perp(const TVector3 &v) const {
429  return fP.Perp(v);
430 }
431 
432 inline Double_t TLorentzVector::Pt(const TVector3 &v) const {
433  return Perp(v);
434 }
435 
437  Double_t pt2 = fP.Perp2();
438  return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+Z()*Z());
439 }
440 
441 inline Double_t TLorentzVector::Et() const {
442  Double_t etet = Et2();
443  return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
444 }
445 
446 inline Double_t TLorentzVector::Et2(const TVector3 & v) const {
447  Double_t pt2 = fP.Perp2(v);
448  Double_t pv = fP.Dot(v.Unit());
449  return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+pv*pv);
450 }
451 
452 inline Double_t TLorentzVector::Et(const TVector3 & v) const {
453  Double_t etet = Et2(v);
454  return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
455 }
456 
458  return TVector2::Phi_mpi_pi(Phi()-v.Phi());
459 }
460 
462  return PseudoRapidity();
463 }
464 
465 inline Double_t TLorentzVector::DeltaR(const TLorentzVector & v, const Bool_t useRapidity) const {
466  if(useRapidity){
467  Double_t drap = Rapidity()-v.Rapidity();
468  Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
469  return TMath::Sqrt( drap*drap+dphi*dphi );
470  } else {
471  Double_t deta = Eta()-v.Eta();
472  Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
473  return TMath::Sqrt( deta*deta+dphi*dphi );
474  }
475 }
476 
478  return DeltaR(v);
479 }
480 
482  return DeltaR(v, kTRUE);
483 }
484 
486  return TVector2 (Eta(),Phi());
487 }
488 
489 
490 inline Double_t TLorentzVector::Angle(const TVector3 &v) const {
491  return fP.Angle(v);
492 }
493 
495  return T()*T() - fP.Mag2();
496 }
497 
499  Double_t mm = Mag2();
500  return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
501 }
502 
503 inline Double_t TLorentzVector::M2() const { return Mag2(); }
504 inline Double_t TLorentzVector::M() const { return Mag(); }
505 
507  return E()*E() - Z()*Z();
508 }
509 
510 inline Double_t TLorentzVector::Mt() const {
511  Double_t mm = Mt2();
512  return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
513 }
514 
516  return fP.Mag() / fE;
517 }
518 
520  Double_t b = Beta();
521  return 1.0/TMath::Sqrt(1- b*b);
522 }
523 
524 inline void TLorentzVector::SetVectMag(const TVector3 & spatial, Double_t magnitude) {
525  SetXYZM(spatial.X(), spatial.Y(), spatial.Z(), magnitude);
526 }
527 
528 inline void TLorentzVector::SetVectM(const TVector3 & spatial, Double_t mass) {
529  SetVectMag(spatial, mass);
530 }
531 
533  return T()*q.T() - Z()*q.Z() - Y()*q.Y() - X()*q.X();
534 }
535 
537  return Dot(q);
538 }
539 
540 //Member functions Plus() and Minus() return the positive and negative
541 //light-cone components:
542 //
543 // Double_t pcone = v.Plus();
544 // Double_t mcone = v.Minus();
545 //
546 //CAVEAT: The values returned are T{+,-}Z. It is known that some authors
547 //find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus
548 //check what definition is used in the physics you're working in and adapt
549 //your code accordingly.
550 
552  return T() + Z();
553 }
554 
556  return T() - Z();
557 }
558 
560  return TVector3(X()/T(), Y()/T(), Z()/T());
561 }
562 
563 inline void TLorentzVector::Boost(const TVector3 & b) {
564  Boost(b.X(), b.Y(), b.Z());
565 }
566 
568  return fP.PseudoRapidity();
569 }
570 
571 inline void TLorentzVector::RotateX(Double_t angle) {
572  fP.RotateX(angle);
573 }
574 
575 inline void TLorentzVector::RotateY(Double_t angle) {
576  fP.RotateY(angle);
577 }
578 
579 inline void TLorentzVector::RotateZ(Double_t angle) {
580  fP.RotateZ(angle);
581 }
582 
583 inline void TLorentzVector::RotateUz(const TVector3 &newUzVector) {
584  fP.RotateUz(newUzVector);
585 }
586 
588  fP.Rotate(a,v);
589 }
590 
592  fP *= m;
593  return *this;
594 }
595 
597  fP.Transform(m);
598  return *this;
599 }
600 
602  return TLorentzVector(a*p.X(), a*p.Y(), a*p.Z(), a*p.T());
603 }
604 
606  : fP(), fE(0.0) {}
607 
609  : fP(x,y,z), fE(t) {}
610 
612  : fP(x0), fE(x0[3]) {}
613 
615  : fP(x0), fE(x0[3]) {}
616 
618  : fP(p), fE(e) {}
619 
621  , fP(p.Vect()), fE(p.T()) {}
622 
623 
624 
626 {
627  //dereferencing operator const
628  switch(i) {
629  case kX:
630  return fP.X();
631  case kY:
632  return fP.Y();
633  case kZ:
634  return fP.Z();
635  case kT:
636  return fE;
637  default:
638  Error("operator()()", "bad index (%d) returning 0",i);
639  }
640  return 0.;
641 }
642 
644 {
645  //dereferencing operator
646  switch(i) {
647  case kX:
648  return fP.fX;
649  case kY:
650  return fP.fY;
651  case kZ:
652  return fP.fZ;
653  case kT:
654  return fE;
655  default:
656  Error("operator()()", "bad index (%d) returning &fE",i);
657  }
658  return fE;
659 }
660 
661 #endif
TVector3::SetPhi
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition: TVector3.cxx:368
m
auto * m
Definition: textangle.C:8
TLorentzVector::operator*
TLorentzVector operator*(Double_t a) const
Definition: TLorentzVector.h:402
TLorentzVector::SetVectMag
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Definition: TLorentzVector.h:524
TLorentzVector::fP
TVector3 fP
Definition: TLorentzVector.h:36
TVector3
Definition: TVector3.h:22
TLorentzVector::operator+=
TLorentzVector & operator+=(const TLorentzVector &)
Definition: TLorentzVector.h:376
TVector3::CosTheta
Double_t CosTheta() const
Definition: TVector3.h:371
TVector3::Transform
TVector3 & Transform(const TRotation &)
Transform this vector with a TRotation.
Definition: TVector3.cxx:196
kTRUE
const Bool_t kTRUE
Definition: RtypesCore.h:100
e
#define e(i)
Definition: RSha256.hxx:103
TLorentzVector::Dot
Double_t Dot(const TLorentzVector &) const
Definition: TLorentzVector.h:532
TLorentzVector::Perp
Double_t Perp() const
Definition: TLorentzVector.h:416
TVector3::SetXYZ
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:227
Option_t
const char Option_t
Definition: RtypesCore.h:66
TVector3::Y
Double_t Y() const
Definition: TVector3.h:217
TLorentzVector::SetPx
void SetPx(Double_t a)
Definition: TLorentzVector.h:289
TLorentzVector::kSIZE
@ kSIZE
Definition: TLorentzVector.h:43
TMath::Max
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
TLorentzVector::kZ
@ kZ
Definition: TLorentzVector.h:43
TLorentzVector::DeltaPhi
Double_t DeltaPhi(const TLorentzVector &) const
Definition: TLorentzVector.h:457
TLorentzVector::RotateUz
void RotateUz(const TVector3 &newUzVector)
Definition: TLorentzVector.h:583
TLorentzVector::SetPtEtaPhiM
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
Definition: TLorentzVector.h:343
TMath::Cos
Double_t Cos(Double_t)
Definition: TMath.h:643
TLorentzVector::kT
@ kT
Definition: TLorentzVector.h:43
TLorentzVector::RotateX
void RotateX(Double_t angle)
Definition: TLorentzVector.h:571
TVector3::fY
Double_t fY
Definition: TVector3.h:185
TVector3::SetMag
void SetMag(Double_t)
Definition: TVector3.h:376
TLorentzVector::SetT
void SetT(Double_t a)
Definition: TLorentzVector.h:280
r
ROOT::R::TRInterface & r
Definition: Object.C:4
TVector3::RotateZ
void RotateZ(Double_t)
Rotate vector around Z.
Definition: TVector3.cxx:285
TObject::Error
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:893
TMath::Sqrt
Double_t Sqrt(Double_t x)
Definition: TMath.h:691
TLorentzVector::operator[]
Double_t operator[](int i) const
Definition: TLorentzVector.h:364
TLorentzVector::Minus
Double_t Minus() const
Definition: TLorentzVector.h:555
TLorentzVector::PseudoRapidity
Double_t PseudoRapidity() const
Definition: TLorentzVector.h:567
TLorentzVector::operator-
TLorentzVector operator-() const
Definition: TLorentzVector.h:392
Float_t
float Float_t
Definition: RtypesCore.h:57
TVector3::RotateUz
void RotateUz(const TVector3 &)
NewUzVector must be normalized !
Definition: TVector3.cxx:305
TVector3::RotateX
void RotateX(Double_t)
Rotate vector around X.
Definition: TVector3.cxx:263
TLorentzVector::kY
@ kY
Definition: TLorentzVector.h:43
x
Double_t x[n]
Definition: legend1.C:17
TLorentzVector::Et
Double_t Et() const
Definition: TLorentzVector.h:441
TLorentzVector::Rho
Double_t Rho() const
Definition: TLorentzVector.h:311
TMath::Abs
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
TLorentzVector::Perp2
Double_t Perp2() const
Definition: TLorentzVector.h:414
TVector3::SetZ
void SetZ(Double_t)
Definition: TVector3.h:225
TLorentzVector::SetPhi
void SetPhi(Double_t phi)
Definition: TLorentzVector.h:319
TLorentzVector::RotateZ
void RotateZ(Double_t angle)
Definition: TLorentzVector.h:579
TVector3::PseudoRapidity
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition: TVector3.cxx:326
TLorentzVector::P
Double_t P() const
Definition: TLorentzVector.h:285
TVector3::Mag2
Double_t Mag2() const
Definition: TVector3.h:339
Bool_t
bool Bool_t
Definition: RtypesCore.h:63
TLorentzVector::Rapidity
Double_t Rapidity() const
Definition: TLorentzVector.cxx:281
v
@ v
Definition: rootcling_impl.cxx:3664
TLorentzVector::Transform
TLorentzVector & Transform(const TRotation &)
Definition: TLorentzVector.h:596
b
#define b(i)
Definition: RSha256.hxx:100
TLorentzVector::Z
Double_t Z() const
Definition: TLorentzVector.h:274
TLorentzVector::E
Double_t E() const
Definition: TLorentzVector.h:286
TVector3::Rotate
void Rotate(Double_t, const TVector3 &)
Rotate vector.
Definition: TVector3.cxx:296
bool
TVector3.h
TLorentzVector::operator-=
TLorentzVector & operator-=(const TLorentzVector &)
Definition: TLorentzVector.h:386
TLorentzVector::M2
Double_t M2() const
Definition: TLorentzVector.h:503
TLorentzVector::Pz
Double_t Pz() const
Definition: TLorentzVector.h:284
q
float * q
Definition: THbookFile.cxx:89
TVector3::SetPerp
void SetPerp(Double_t)
Definition: TVector3.h:388
TLorentzVector::SetVect
void SetVect(const TVector3 &vect3)
Definition: TLorentzVector.h:296
TVector3::Mag
Double_t Mag() const
Definition: TVector3.h:86
TLorentzVector::EtaPhiVector
TVector2 EtaPhiVector()
Definition: TLorentzVector.h:485
TLorentzRotation
The TLorentzRotation class describes Lorentz transformations including Lorentz boosts and rotations (...
Definition: TLorentzRotation.h:20
TLorentzVector::DrEtaPhi
Double_t DrEtaPhi(const TLorentzVector &) const
Definition: TLorentzVector.h:477
TVector3::fX
Double_t fX
Definition: TVector3.h:185
TLorentzVector::SetPtEtaPhiE
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
Definition: TLorentzVector.h:348
TLorentzVector::SetX
void SetX(Double_t a)
Definition: TLorentzVector.h:277
TLorentzVector::SetZ
void SetZ(Double_t a)
Definition: TLorentzVector.h:279
TLorentzVector::DeltaR
Double_t DeltaR(const TLorentzVector &, Bool_t useRapidity=kFALSE) const
Definition: TLorentzVector.h:465
TVector3::Angle
Double_t Angle(const TVector3 &) const
Return the angle w.r.t. another 3-vector.
Definition: TVector3.cxx:203
TLorentzVector::Rotate
void Rotate(Double_t, const TVector3 &)
Definition: TLorentzVector.h:587
TVector3::Theta
Double_t Theta() const
Return the polar angle.
Definition: TVector3.cxx:244
TLorentzVector::T
Double_t T() const
Definition: TLorentzVector.h:275
TLorentzVector::CosTheta
Double_t CosTheta() const
Definition: TLorentzVector.h:306
TLorentzVector::Print
virtual void Print(Option_t *option="") const
Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations.
Definition: TLorentzVector.cxx:327
TLorentzVector::Vect
TVector3 Vect() const
Definition: TLorentzVector.h:294
TLorentzVector::X
Double_t X() const
Definition: TLorentzVector.h:272
a
auto * a
Definition: textangle.C:12
TLorentzVector::Px
Double_t Px() const
Definition: TLorentzVector.h:282
TLorentzVector::operator==
Bool_t operator==(const TLorentzVector &) const
Definition: TLorentzVector.h:406
kFALSE
const Bool_t kFALSE
Definition: RtypesCore.h:101
TLorentzVector
Definition: TLorentzVector.h:32
TLorentzVector::kNUM_COORDINATES
@ kNUM_COORDINATES
Definition: TLorentzVector.h:43
TVector3::SetX
void SetX(Double_t)
Definition: TVector3.h:223
TLorentzVector::SetXYZT
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Definition: TLorentzVector.h:327
TLorentzVector::M
Double_t M() const
Definition: TLorentzVector.h:504
TVector3::fZ
Double_t fZ
Definition: TVector3.h:185
TLorentzVector::SetPxPyPzE
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: TLorentzVector.h:332
TLorentzVector::SetPy
void SetPy(Double_t a)
Definition: TLorentzVector.h:290
TLorentzVector::Gamma
Double_t Gamma() const
Definition: TLorentzVector.h:519
TLorentzVector::Eta
Double_t Eta() const
Definition: TLorentzVector.h:461
TMath::Sin
Double_t Sin(Double_t)
Definition: TMath.h:639
y
Double_t y[n]
Definition: legend1.C:17
TLorentzVector::fE
Double_t fE
Definition: TLorentzVector.h:37
TLorentzVector::Theta
Double_t Theta() const
Definition: TLorentzVector.h:302
sqrt
double sqrt(double)
operator*
TLorentzVector operator*(Double_t a, const TLorentzVector &)
Definition: TLorentzVector.h:601
TLorentzVector::Mag
Double_t Mag() const
Definition: TLorentzVector.h:498
TLorentzVector::Mt2
Double_t Mt2() const
Definition: TLorentzVector.h:506
TLorentzVector::SetTheta
void SetTheta(Double_t theta)
Definition: TLorentzVector.h:315
TGeant4Unit::mm
static constexpr double mm
Definition: TGeant4SystemOfUnits.h:108
TRotation
The TRotation class describes a rotation of objects of the TVector3 class.
Definition: TRotation.h:20
TVector3::Perp2
Double_t Perp2() const
Definition: TVector3.h:353
TLorentzVector::Y
Double_t Y() const
Definition: TLorentzVector.h:273
TLorentzVector::Plus
Double_t Plus() const
Definition: TLorentzVector.h:551
TLorentzVector::SetVectM
void SetVectM(const TVector3 &spatial, Double_t mass)
Definition: TLorentzVector.h:528
TVector3::Perp
Double_t Perp() const
Return the transverse component (R in cylindrical coordinate system)
Definition: TVector3.cxx:219
TVector3::RotateY
void RotateY(Double_t)
Rotate vector around Y.
Definition: TVector3.cxx:274
TLorentzVector::Boost
void Boost(Double_t, Double_t, Double_t)
Definition: TLorentzVector.cxx:267
TVector3::Phi
Double_t Phi() const
Return the azimuth angle. Returns phi from -pi to pi.
Definition: TVector3.cxx:236
TVector3::GetXYZ
void GetXYZ(Double_t *carray) const
Definition: TVector3.h:233
TLorentzVector::Et2
Double_t Et2() const
Definition: TLorentzVector.h:436
TLorentzVector::Py
Double_t Py() const
Definition: TLorentzVector.h:283
TVector3::Dot
Double_t Dot(const TVector3 &) const
Definition: TVector3.h:331
Double_t
double Double_t
Definition: RtypesCore.h:59
TLorentzVector::Mag2
Double_t Mag2() const
Definition: TLorentzVector.h:494
TLorentzVector::operator=
TLorentzVector & operator=(const TLorentzVector &)
Definition: TLorentzVector.h:366
TLorentzVector::~TLorentzVector
virtual ~TLorentzVector()
Definition: TLorentzVector.h:61
TLorentzVector::operator!=
Bool_t operator!=(const TLorentzVector &) const
Definition: TLorentzVector.h:410
TLorentzVector::TLorentzVector
TLorentzVector()
Definition: TLorentzVector.h:605
TLorentzVector::operator+
TLorentzVector operator+(const TLorentzVector &) const
Definition: TLorentzVector.h:372
TLorentzVector::SetPerp
void SetPerp(Double_t)
Definition: TLorentzVector.h:420
TObject
Mother of all ROOT objects.
Definition: TObject.h:37
TLorentzVector::operator()
Double_t operator()(int i) const
Definition: TLorentzVector.h:625
ClassDef
#define ClassDef(name, id)
Definition: Rtypes.h:325
TLorentzVector::Pt
Double_t Pt() const
Definition: TLorentzVector.h:418
TVector3::Z
Double_t Z() const
Definition: TVector3.h:218
TLorentzVector::BoostVector
TVector3 BoostVector() const
Definition: TLorentzVector.h:559
TLorentzVector::Mt
Double_t Mt() const
Definition: TLorentzVector.h:510
ROOT::Math::Chebyshev::T
double T(double x)
Definition: ChebyshevPol.h:34
TVector2
TVector2 is a general two vector class, which can be used for the description of different vectors in...
Definition: TVector2.h:18
sinh
double sinh(double)
TRotation.h
TLorentzVector::operator*=
TLorentzVector & operator*=(Double_t a)
Definition: TLorentzVector.h:396
TLorentzVector::SetXYZM
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Definition: TLorentzVector.h:336
TLorentzVector::SetPz
void SetPz(Double_t a)
Definition: TLorentzVector.h:291
TLorentzVector::Angle
Double_t Angle(const TVector3 &v) const
Definition: TLorentzVector.h:490
TLorentzVector::kX
@ kX
Definition: TLorentzVector.h:43
TLorentzVector::DrRapidityPhi
Double_t DrRapidityPhi(const TLorentzVector &) const
Definition: TLorentzVector.h:481
TVector3::X
Double_t X() const
Definition: TVector3.h:216
pt
TPaveText * pt
Definition: entrylist_figure1.C:7
TLorentzVector::GetXYZT
void GetXYZT(Double_t *carray) const
Definition: TLorentzVector.h:353
TLorentzVector::Scalar
Double_t Scalar
Definition: TLorentzVector.h:41
TLorentzVector::Phi
Double_t Phi() const
Definition: TLorentzVector.h:298
TLorentzVector::SetE
void SetE(Double_t a)
Definition: TLorentzVector.h:292
TVector3::SetTheta
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition: TVector3.cxx:356
TLorentzVector::Beta
Double_t Beta() const
Definition: TLorentzVector.h:515
TLorentzVector::RotateY
void RotateY(Double_t angle)
Definition: TLorentzVector.h:575
TLorentzVector::Energy
Double_t Energy() const
Definition: TLorentzVector.h:287
TLorentzVector::SetY
void SetY(Double_t a)
Definition: TLorentzVector.h:278
TVector2::Phi_mpi_pi
static Double_t Phi_mpi_pi(Double_t x)
Returns phi angle in the interval [-PI,PI)
Definition: TVector2.cxx:103
TVector3::SetY
void SetY(Double_t)
Definition: TVector3.h:224
TMath.h
TLorentzVector::SetRho
void SetRho(Double_t rho)
Definition: TLorentzVector.h:323