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