Logo ROOT  
Reference Guide
TLorentzVector.cxx
Go to the documentation of this file.
1// @(#)root/physics:$Id$
2// Author: Pasha Murat , Peter Malzacher 12/02/99
3// Oct 8 1999: changed Warning to Error and
4// return fX in Double_t & operator()
5// Oct 20 1999: dito in Double_t operator()
6// Jan 25 2000: implemented as (fP,fE) instead of (fX,fY,fZ,fE)
7
8/** \class TLorentzVector
9 \ingroup Physics
10
11## Disclaimer
12TLorentzVector is a legacy class.
13It is slower and worse for serialization than the recommended superior ROOT::Math::LorentzVector.
14ROOT provides specialisations of the ROOT::Math::LorentzVector template which
15are superior from the runtime performance offered, i.e.:
16 - ROOT::Math::PtEtaPhiMVector based on pt (rho),eta,phi and M (t) coordinates in double precision
17 - ROOT::Math::PtEtaPhiEVector based on pt (rho),eta,phi and E (t) coordinates in double precision
18 - ROOT::Math::PxPyPzMVector based on px,py,pz and M (mass) coordinates in double precision
19 - ROOT::Math::PxPyPzEVector based on px,py,pz and E (energy) coordinates in double precision
20 - ROOT::Math::XYZTVector based on x,y,z,t coordinates (cartesian) in double precision (same as PxPyPzEVector)
21 - ROOT::Math::XYZTVectorF based on x,y,z,t coordinates (cartesian) in float precision (same as PxPyPzEVector but float)
22
23More details about the GenVector package can be found [here](Vector.html).
24
25### Description
26TLorentzVector is a general four-vector class, which can be used
27either for the description of position and time (x,y,z,t) or momentum and
28energy (px,py,pz,E).
29
30### Declaration
31TLorentzVector has been implemented as a set a TVector3 and a Double_t variable.
32By default all components are initialized by zero.
33
34~~~ {.cpp}
35 TLorentzVector v1; // initialized by (0., 0., 0., 0.)
36 TLorentzVector v2(1., 1., 1., 1.);
37 TLorentzVector v3(v1);
38 TLorentzVector v4(TVector3(1., 2., 3.),4.);
39~~~
40
41For backward compatibility there are two constructors from an Double_t
42and Float_t C array.
43
44
45### Access to the components
46There are two sets of access functions to the components of a LorentzVector:
47X(), Y(), Z(), T() and Px(),
48Py(), Pz() and E(). Both sets return the same values
49but the first set is more relevant for use where TLorentzVector
50describes a combination of position and time and the second set is more
51relevant where TLorentzVector describes momentum and energy:
52
53~~~ {.cpp}
54 Double_t xx =v.X();
55 ...
56 Double_t tt = v.T();
57
58 Double_t px = v.Px();
59 ...
60 Double_t ee = v.E();
61~~~
62
63The components of TLorentzVector can also accessed by index:
64
65~~~ {.cpp}
66 xx = v(0); or xx = v[0];
67 yy = v(1); yy = v[1];
68 zz = v(2); zz = v[2];
69 tt = v(3); tt = v[3];
70~~~
71
72You can use the Vect() member function to get the vector component
73of TLorentzVector:
74
75~~~ {.cpp}
76 TVector3 p = v.Vect();
77~~~
78
79For setting components also two sets of member functions can be used:
80
81~~~ {.cpp}
82 v.SetX(1.); or v.SetPx(1.);
83 ... ...
84 v.SetT(1.); v.SetE(1.);
85~~~
86
87To set more the one component by one call you can use the SetVect()
88function for the TVector3 part or SetXYZT(), SetPxPyPzE(). For convenience there is
89
90also a SetXYZM():
91
92~~~ {.cpp}
93 v.SetVect(TVector3(1,2,3));
94 v.SetXYZT(x,y,z,t);
95 v.SetPxPyPzE(px,py,pz,e);
96 v.SetXYZM(x,y,z,m); // -> v=(x,y,z,e=Sqrt(x*x+y*y+z*z+m*m))
97~~~
98
99### Vector components in non-cartesian coordinate systems
100There are a couple of member functions to get and set the TVector3
101part of the parameters in
102spherical coordinate systems:
103
104~~~ {.cpp}
105 Double_t m, theta, cost, phi, pp, pp2, ppv2, pp2v2;
106 m = v.Rho();
107 t = v.Theta();
108 cost = v.CosTheta();
109 phi = v.Phi();
110
111 v.SetRho(10.);
112 v.SetTheta(TMath::Pi()*.3);
113 v.SetPhi(TMath::Pi());
114~~~
115
116or get information about the r-coordinate in cylindrical systems:
117
118~~~ {.cpp}
119 Double_t pp, pp2, ppv2, pp2v2;
120 pp = v.Perp(); // get transvers component
121 pp2 = v.Perp2(); // get transverse component squared
122 ppv2 = v.Perp(v1); // get transvers component with
123 // respect to another vector
124 pp2v2 = v.Perp(v1);
125~~~
126
127for convenience there are two more set functions SetPtEtaPhiE(pt,eta,phi,e);
128and SetPtEtaPhiM(pt,eta,phi,m);
129
130### Arithmetic and comparison operators
131The TLorentzVector class provides operators to add, subtract or
132compare four-vectors:
133
134~~~ {.cpp}
135 v3 = -v1;
136 v1 = v2+v3;
137 v1+= v3;
138 v1 = v2 + v3;
139 v1-= v3;
140
141 if (v1 == v2) {...}
142 if(v1 != v3) {...}
143~~~
144
145### Magnitude/Invariant mass, beta, gamma, scalar product
146The scalar product of two four-vectors is calculated with the (-,-,-,+)
147metric,
148
149 i.e. `s = v1*v2 = t1*t2-x1*x2-y1*y2-z1*z2`
150The magnitude squared mag2 of a four-vector is therefore:
151
152~~~ {.cpp}
153 mag2 = v*v = t*t-x*x-y*y-z*z
154~~~
155It mag2 is negative mag = -Sqrt(-mag*mag). The member
156functions are:
157
158~~~ {.cpp}
159 Double_t s, s2;
160 s = v1.Dot(v2); // scalar product
161 s = v1*v2; // scalar product
162 s2 = v.Mag2(); or s2 = v.M2();
163 s = v.Mag(); s = v.M();
164~~~
165
166Since in case of momentum and energy the magnitude has the meaning of
167invariant mass TLorentzVector provides the more meaningful aliases
168M2() and M();
169The member functions Beta() and Gamma() returns
170beta and gamma = 1/Sqrt(1-beta*beta).
171### Lorentz boost
172A boost in a general direction can be parameterised with three parameters
173which can be taken as the components of a three vector b = (bx,by,bz).
174With x = (x,y,z) and gamma = 1/Sqrt(1-beta*beta) (beta being the module of vector b),
175an arbitrary active Lorentz boost transformation (from the rod frame
176to the original frame) can be written as:
177
178~~~ {.cpp}
179 x = x' + (gamma-1)/(beta*beta) * (b*x') * b + gamma * t' * b
180 t = gamma (t'+ b*x').
181~~~
182
183The member function Boost() performs a boost transformation
184from the rod frame to the original frame. BoostVector() returns
185a TVector3 of the spatial components divided by the time component:
186
187~~~ {.cpp}
188 TVector3 b;
189 v.Boost(bx,by,bz);
190 v.Boost(b);
191 b = v.BoostVector(); // b=(x/t,y/t,z/t)
192~~~
193
194### Rotations
195There are four sets of functions to rotate the TVector3 component
196of a TLorentzVector:
197
198#### rotation around axes
199
200~~~ {.cpp}
201 v.RotateX(TMath::Pi()/2.);
202 v.RotateY(.5);
203 v.RotateZ(.99);
204~~~
205
206#### rotation around an arbitrary axis
207 v.Rotate(TMath::Pi()/4., v1); // rotation around v1
208
209#### transformation from rotated frame
210
211~~~ {.cpp}
212 v.RotateUz(direction); // direction must be a unit TVector3
213~~~
214
215#### by TRotation (see TRotation)
216
217~~~ {.cpp}
218 TRotation r;
219 v.Transform(r); or v *= r; // Attention v=M*v
220~~~
221
222### Misc
223
224#### Angle between two vectors
225
226~~~ {.cpp}
227 Double_t a = v1.Angle(v2.Vect()); // get angle between v1 and v2
228~~~
229
230#### Light-cone components
231Member functions Plus() and Minus() return the positive
232and negative light-cone components:
233
234~~~ {.cpp}
235 Double_t pcone = v.Plus();
236 Double_t mcone = v.Minus();
237~~~
238
239CAVEAT: The values returned are T{+,-}Z. It is known that some authors
240find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus
241check what definition is used in the physics you're working in and adapt
242your code accordingly.
243
244#### Transformation by TLorentzRotation
245A general Lorentz transformation see class TLorentzRotation can
246be used by the Transform() member function, the *= or
247* operator of the TLorentzRotation class:
248
249~~~ {.cpp}
250 TLorentzRotation l;
251 v.Transform(l);
252 v = l*v; or v *= l; // Attention v = l*v
253~~~
254*/
255
256#include "TLorentzVector.h"
257
258#include "TBuffer.h"
259#include "TString.h"
260#include "TLorentzRotation.h"
261
263
264
266{
267 //Boost this Lorentz vector
268 Double_t b2 = bx*bx + by*by + bz*bz;
269 Double_t gamma = 1.0 / TMath::Sqrt(1.0 - b2);
270 Double_t bp = bx*X() + by*Y() + bz*Z();
271 Double_t gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
272
273 SetX(X() + gamma2*bp*bx + gamma*bx*T());
274 SetY(Y() + gamma2*bp*by + gamma*by*T());
275 SetZ(Z() + gamma2*bp*bz + gamma*bz*T());
276 SetT(gamma*(T() + bp));
277}
278
280{
281 //return rapidity
282 return 0.5*log( (E()+Pz()) / (E()-Pz()) );
283}
284
286{
287 //multiply this Lorentzvector by m
288 return *this = m.VectorMultiplication(*this);
289}
290
292{
293 //Transform this Lorentzvector
294 return *this = m.VectorMultiplication(*this);
295}
296
297void TLorentzVector::Streamer(TBuffer &R__b)
298{
299 // Stream an object of class TLorentzVector.
300 Double_t x, y, z;
301 UInt_t R__s, R__c;
302 if (R__b.IsReading()) {
303 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
304 if (R__v > 3) {
305 R__b.ReadClassBuffer(TLorentzVector::Class(), this, R__v, R__s, R__c);
306 return;
307 }
308 //====process old versions before automatic schema evolution
309 if (R__v != 2) TObject::Streamer(R__b);
310 R__b >> x;
311 R__b >> y;
312 R__b >> z;
313 fP.SetXYZ(x,y,z);
314 R__b >> fE;
315 R__b.CheckByteCount(R__s, R__c, TLorentzVector::IsA());
316 } else {
318 }
319}
320
321
322////////////////////////////////////////////////////////////////////////////////
323/// Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations
324
326{
327 Printf("(x,y,z,t)=(%f,%f,%f,%f) (P,eta,phi,E)=(%f,%f,%f,%f)",
328 fP.x(),fP.y(),fP.z(),fE,
329 P(),Eta(),Phi(),fE);
330}
void Class()
Definition: Class.C:29
short Version_t
Definition: RtypesCore.h:63
unsigned int UInt_t
Definition: RtypesCore.h:44
double Double_t
Definition: RtypesCore.h:57
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
double log(double)
void Printf(const char *fmt,...)
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition: TBuffer.h:85
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
The TLorentzRotation class describes Lorentz transformations including Lorentz boosts and rotations (...
Double_t Rapidity() const
TLorentzVector & operator*=(Double_t a)
void SetY(Double_t a)
Double_t Y() const
void SetT(Double_t a)
Double_t Pz() const
Double_t X() const
void Boost(Double_t, Double_t, Double_t)
Double_t Eta() const
Double_t P() const
TLorentzVector & Transform(const TRotation &)
virtual void Print(Option_t *option="") const
Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations.
Double_t Phi() const
void SetZ(Double_t a)
Double_t E() const
Double_t Z() const
Double_t T() const
void SetX(Double_t a)
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:227
Double_t z() const
Definition: TVector3.h:215
Double_t x() const
Definition: TVector3.h:213
Double_t y() const
Definition: TVector3.h:214
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
double gamma(double x)
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
auto * m
Definition: textangle.C:8