12 #ifndef ROOT_TLorentzVector
13 #define ROOT_TLorentzVector
29 class TLorentzRotation;
32 class TLorentzVector :
public TObject {
41 typedef Double_t Scalar;
43 enum { kX=0, kY=1, kZ=2, kT=3, kNUM_COORDINATES=4, kSIZE=kNUM_COORDINATES };
48 TLorentzVector(Double_t x, Double_t y, Double_t z, Double_t t);
51 TLorentzVector(
const Double_t * carray);
52 TLorentzVector(
const Float_t * carray);
55 TLorentzVector(
const TVector3 & vector3, Double_t t);
58 TLorentzVector(
const TLorentzVector & lorentzvector);
61 virtual ~TLorentzVector(){};
68 inline Double_t X()
const;
69 inline Double_t Y()
const;
70 inline Double_t Z()
const;
71 inline Double_t T()
const;
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);
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;
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);
94 inline TVector3 Vect()
const ;
97 inline void SetVect(
const TVector3 & vect3);
100 inline Double_t Theta()
const;
101 inline Double_t CosTheta()
const;
102 inline Double_t Phi()
const;
103 inline Double_t Rho()
const;
106 inline void SetTheta(Double_t theta);
107 inline void SetPhi(Double_t phi);
108 inline void SetRho(Double_t rho);
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);
119 inline void GetXYZT(Double_t *carray)
const;
120 inline void GetXYZT(Float_t *carray)
const;
124 Double_t operator () (
int i)
const;
125 inline Double_t operator [] (
int i)
const;
128 Double_t & operator () (
int i);
129 inline Double_t & operator [] (
int i);
132 inline TLorentzVector & operator = (
const TLorentzVector &);
135 inline TLorentzVector operator + (
const TLorentzVector &)
const;
136 inline TLorentzVector & operator += (
const TLorentzVector &);
139 inline TLorentzVector operator - (
const TLorentzVector &)
const;
140 inline TLorentzVector & operator -= (
const TLorentzVector &);
143 inline TLorentzVector operator - ()
const;
146 inline TLorentzVector operator * (Double_t a)
const;
147 inline TLorentzVector & operator *= (Double_t a);
150 inline Bool_t operator == (
const TLorentzVector &)
const;
151 inline Bool_t operator != (
const TLorentzVector &)
const;
154 inline Double_t Perp2()
const;
157 inline Double_t Pt()
const;
158 inline Double_t Perp()
const;
161 inline void SetPerp(Double_t);
164 inline Double_t Perp2(
const TVector3 & v)
const;
167 inline Double_t Pt(
const TVector3 & v)
const;
168 inline Double_t Perp(
const TVector3 & v)
const;
171 inline Double_t Et2()
const;
174 inline Double_t Et()
const;
177 inline Double_t Et2(
const TVector3 &)
const;
180 inline Double_t Et(
const TVector3 &)
const;
183 inline Double_t DeltaPhi(
const TLorentzVector &)
const;
184 inline Double_t DeltaR(
const TLorentzVector &)
const;
185 inline Double_t DrEtaPhi(
const TLorentzVector &)
const;
186 inline TVector2 EtaPhiVector();
188 inline Double_t Angle(
const TVector3 & v)
const;
191 inline Double_t Mag2()
const;
192 inline Double_t M2()
const;
195 inline Double_t Mag()
const;
196 inline Double_t M()
const;
199 inline Double_t Mt2()
const;
202 inline Double_t Mt()
const;
205 inline Double_t Beta()
const;
206 inline Double_t Gamma()
const;
208 inline Double_t Dot(
const TLorentzVector &)
const;
209 inline Double_t operator * (
const TLorentzVector &)
const;
212 inline void SetVectMag(
const TVector3 & spatial, Double_t magnitude);
213 inline void SetVectM(
const TVector3 & spatial, Double_t mass);
216 inline Double_t Plus()
const;
217 inline Double_t Minus()
const;
222 inline TVector3 BoostVector()
const ;
225 void Boost(Double_t, Double_t, Double_t);
226 inline void Boost(
const TVector3 &);
229 Double_t Rapidity()
const;
232 inline Double_t Eta()
const;
233 inline Double_t PseudoRapidity()
const;
236 inline void RotateX(Double_t angle);
239 inline void RotateY(Double_t angle);
242 inline void RotateZ(Double_t angle);
245 inline void RotateUz(
const TVector3 & newUzVector);
248 inline void Rotate(Double_t,
const TVector3 &);
251 inline TLorentzVector & operator *= (
const TRotation &);
252 inline TLorentzVector & Transform(
const TRotation &);
255 TLorentzVector & operator *= (
const TLorentzRotation &);
256 TLorentzVector & Transform(
const TLorentzRotation &);
259 virtual void Print(Option_t *option=
"")
const;
261 ClassDef(TLorentzVector,4)
267 inline TLorentzVector operator * (Double_t a,
const TLorentzVector &);
271 inline Double_t TLorentzVector::X()
const {
return fP.X(); }
272 inline Double_t TLorentzVector::Y()
const {
return fP.Y(); }
273 inline Double_t TLorentzVector::Z()
const {
return fP.Z(); }
274 inline Double_t TLorentzVector::T()
const {
return fE; }
276 inline void TLorentzVector::SetX(Double_t a) { fP.SetX(a); }
277 inline void TLorentzVector::SetY(Double_t a) { fP.SetY(a); }
278 inline void TLorentzVector::SetZ(Double_t a) { fP.SetZ(a); }
279 inline void TLorentzVector::SetT(Double_t a) { fE = a; }
281 inline Double_t TLorentzVector::Px()
const {
return X(); }
282 inline Double_t TLorentzVector::Py()
const {
return Y(); }
283 inline Double_t TLorentzVector::Pz()
const {
return Z(); }
284 inline Double_t TLorentzVector::P()
const {
return fP.Mag(); }
285 inline Double_t TLorentzVector::E()
const {
return T(); }
286 inline Double_t TLorentzVector::Energy()
const {
return T(); }
288 inline void TLorentzVector::SetPx(Double_t a) { SetX(a); }
289 inline void TLorentzVector::SetPy(Double_t a) { SetY(a); }
290 inline void TLorentzVector::SetPz(Double_t a) { SetZ(a); }
291 inline void TLorentzVector::SetE(Double_t a) { SetT(a); }
293 inline TVector3 TLorentzVector::Vect()
const {
return fP; }
295 inline void TLorentzVector::SetVect(
const TVector3 &p) { fP = p; }
297 inline Double_t TLorentzVector::Phi()
const {
301 inline Double_t TLorentzVector::Theta()
const {
305 inline Double_t TLorentzVector::CosTheta()
const {
306 return fP.CosTheta();
310 inline Double_t TLorentzVector::Rho()
const {
314 inline void TLorentzVector::SetTheta(Double_t th) {
318 inline void TLorentzVector::SetPhi(Double_t phi) {
322 inline void TLorentzVector::SetRho(Double_t rho) {
326 inline void TLorentzVector::SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t) {
331 inline void TLorentzVector::SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e) {
332 SetXYZT(px, py, pz, e);
335 inline void TLorentzVector::SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m) {
337 SetXYZT( x, y, z, TMath::Sqrt(x*x+y*y+z*z+m*m) );
339 SetXYZT( x, y, z, TMath::Sqrt( TMath::Max((x*x+y*y+z*z-m*m), 0. ) ) );
342 inline void TLorentzVector::SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m) {
344 SetXYZM(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,m);
347 inline void TLorentzVector::SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e) {
349 SetXYZT(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,e);
352 inline void TLorentzVector::GetXYZT(Double_t *carray)
const {
357 inline void TLorentzVector::GetXYZT(Float_t *carray)
const{
362 inline Double_t & TLorentzVector::operator [] (
int i) {
return (*
this)(i); }
363 inline Double_t TLorentzVector::operator [] (
int i)
const {
return (*
this)(i); }
365 inline TLorentzVector &TLorentzVector::operator = (
const TLorentzVector & q) {
371 inline TLorentzVector TLorentzVector::operator + (
const TLorentzVector & q)
const {
372 return TLorentzVector(fP+q.Vect(), fE+q.T());
375 inline TLorentzVector &TLorentzVector::operator += (
const TLorentzVector & q) {
381 inline TLorentzVector TLorentzVector::operator - (
const TLorentzVector & q)
const {
382 return TLorentzVector(fP-q.Vect(), fE-q.T());
385 inline TLorentzVector &TLorentzVector::operator -= (
const TLorentzVector & q) {
391 inline TLorentzVector TLorentzVector::operator - ()
const {
392 return TLorentzVector(-X(), -Y(), -Z(), -T());
395 inline TLorentzVector& TLorentzVector::operator *= (Double_t a) {
401 inline TLorentzVector TLorentzVector::operator * (Double_t a)
const {
402 return TLorentzVector(a*X(), a*Y(), a*Z(), a*T());
405 inline Bool_t TLorentzVector::operator == (
const TLorentzVector & q)
const {
406 return (Vect() == q.Vect() && T() == q.T());
409 inline Bool_t TLorentzVector::operator != (
const TLorentzVector & q)
const {
410 return (Vect() != q.Vect() || T() != q.T());
413 inline Double_t TLorentzVector::Perp2()
const {
return fP.Perp2(); }
415 inline Double_t TLorentzVector::Perp()
const {
return fP.Perp(); }
417 inline Double_t TLorentzVector::Pt()
const {
return Perp(); }
419 inline void TLorentzVector::SetPerp(Double_t r) {
423 inline Double_t TLorentzVector::Perp2(
const TVector3 &v)
const {
427 inline Double_t TLorentzVector::Perp(
const TVector3 &v)
const {
431 inline Double_t TLorentzVector::Pt(
const TVector3 &v)
const {
435 inline Double_t TLorentzVector::Et2()
const {
436 Double_t pt2 = fP.Perp2();
437 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+Z()*Z());
440 inline Double_t TLorentzVector::Et()
const {
441 Double_t etet = Et2();
442 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
445 inline Double_t TLorentzVector::Et2(
const TVector3 & v)
const {
446 Double_t pt2 = fP.Perp2(v);
447 Double_t pv = fP.Dot(v.Unit());
448 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+pv*pv);
451 inline Double_t TLorentzVector::Et(
const TVector3 & v)
const {
452 Double_t etet = Et2(v);
453 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
456 inline Double_t TLorentzVector::DeltaPhi(
const TLorentzVector & v)
const {
457 return TVector2::Phi_mpi_pi(Phi()-v.Phi());
460 inline Double_t TLorentzVector::Eta()
const {
461 return PseudoRapidity();
463 inline Double_t TLorentzVector::DeltaR(
const TLorentzVector & v)
const {
464 Double_t deta = Eta()-v.Eta();
465 Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
466 return TMath::Sqrt( deta*deta+dphi*dphi );
469 inline Double_t TLorentzVector::DrEtaPhi(
const TLorentzVector & v)
const{
473 inline TVector2 TLorentzVector::EtaPhiVector() {
474 return TVector2 (Eta(),Phi());
478 inline Double_t TLorentzVector::Angle(
const TVector3 &v)
const {
482 inline Double_t TLorentzVector::Mag2()
const {
483 return T()*T() - fP.Mag2();
486 inline Double_t TLorentzVector::Mag()
const {
487 Double_t mm = Mag2();
488 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
491 inline Double_t TLorentzVector::M2()
const {
return Mag2(); }
492 inline Double_t TLorentzVector::M()
const {
return Mag(); }
494 inline Double_t TLorentzVector::Mt2()
const {
495 return E()*E() - Z()*Z();
498 inline Double_t TLorentzVector::Mt()
const {
500 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
503 inline Double_t TLorentzVector::Beta()
const {
504 return fP.Mag() / fE;
507 inline Double_t TLorentzVector::Gamma()
const {
509 return 1.0/TMath::Sqrt(1- b*b);
512 inline void TLorentzVector::SetVectMag(
const TVector3 & spatial, Double_t magnitude) {
513 SetXYZM(spatial.X(), spatial.Y(), spatial.Z(), magnitude);
516 inline void TLorentzVector::SetVectM(
const TVector3 & spatial, Double_t mass) {
517 SetVectMag(spatial, mass);
520 inline Double_t TLorentzVector::Dot(
const TLorentzVector & q)
const {
521 return T()*q.T() - Z()*q.Z() - Y()*q.Y() - X()*q.X();
524 inline Double_t TLorentzVector::operator * (
const TLorentzVector & q)
const {
539 inline Double_t TLorentzVector::Plus()
const {
543 inline Double_t TLorentzVector::Minus()
const {
547 inline TVector3 TLorentzVector::BoostVector()
const {
548 return TVector3(X()/T(), Y()/T(), Z()/T());
551 inline void TLorentzVector::Boost(
const TVector3 & b) {
552 Boost(b.X(), b.Y(), b.Z());
555 inline Double_t TLorentzVector::PseudoRapidity()
const {
556 return fP.PseudoRapidity();
559 inline void TLorentzVector::RotateX(Double_t angle) {
563 inline void TLorentzVector::RotateY(Double_t angle) {
567 inline void TLorentzVector::RotateZ(Double_t angle) {
571 inline void TLorentzVector::RotateUz(
const TVector3 &newUzVector) {
572 fP.RotateUz(newUzVector);
575 inline void TLorentzVector::Rotate(Double_t a,
const TVector3 &v) {
579 inline TLorentzVector &TLorentzVector::operator *= (
const TRotation & m) {
584 inline TLorentzVector &TLorentzVector::Transform(
const TRotation & m) {
589 inline TLorentzVector operator * (Double_t a,
const TLorentzVector & p) {
590 return TLorentzVector(a*p.X(), a*p.Y(), a*p.Z(), a*p.T());
593 inline TLorentzVector::TLorentzVector()
596 inline TLorentzVector::TLorentzVector(Double_t x, Double_t y, Double_t z, Double_t t)
597 : fP(x,y,z), fE(t) {}
599 inline TLorentzVector::TLorentzVector(
const Double_t * x0)
600 : fP(x0), fE(x0[3]) {}
602 inline TLorentzVector::TLorentzVector(
const Float_t * x0)
603 : fP(x0), fE(x0[3]) {}
605 inline TLorentzVector::TLorentzVector(
const TVector3 & p, Double_t e)
608 inline TLorentzVector::TLorentzVector(
const TLorentzVector & p) : TObject(p)
609 , fP(p.Vect()), fE(p.T()) {}
613 inline Double_t TLorentzVector::operator () (
int i)
const
626 Error(
"operator()()",
"bad index (%d) returning 0",i);
631 inline Double_t & TLorentzVector::operator () (
int i)
644 Error(
"operator()()",
"bad index (%d) returning &fE",i);