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