ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 //______________________________________________________________________________
9 //*-*-*-*-*-*-*-*-*-*-*-*The Physics Vector package *-*-*-*-*-*-*-*-*-*-*-*
10 //*-* ========================== *
11 //*-* The Physics Vector package consists of five classes: *
12 //*-* - TVector2 *
13 //*-* - TVector3 *
14 //*-* - TRotation *
15 //*-* - TLorentzVector *
16 //*-* - TLorentzRotation *
17 //*-* It is a combination of CLHEPs Vector package written by *
18 //*-* Leif Lonnblad, Andreas Nilsson and Evgueni Tcherniaev *
19 //*-* and a ROOT package written by Pasha Murat. *
20 //*-* for CLHEP see: http://wwwinfo.cern.ch/asd/lhc++/clhep/ *
21 //*-* Adaption to ROOT by Peter Malzacher *
22 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
23 //BEGIN_HTML <!--
24 /* -->
25 <H2>TLorentzVector</H2>
26 <TT>TLorentzVector</TT> is a general four-vector class, which can be used
27 either for the description of position and time (x,y,z,t) or momentum and
28 energy (px,py,pz,E).
29 <BR>&nbsp;
30 <H3>
31 Declaration</H3>
32 <TT>TLorentzVector</TT> has been implemented as a set a <TT>TVector3</TT> and a <TT>Double_t</TT> variable.
33 By default all components are initialized by zero.
34 
35 <P><TT>&nbsp; TLorentzVector v1;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // initialized
36 by (0., 0., 0., 0.)</TT>
37 <BR><TT>&nbsp; TLorentzVector v2(1., 1., 1., 1.);</TT>
38 <BR><TT>&nbsp; TLorentzVector v3(v1);</TT>
39 <BR><TT>&nbsp; TLorentzVector v4(TVector3(1., 2., 3.),4.);</TT>
40 
41 <P>For backward compatibility there are two constructors from an <TT>Double_t</TT>
42 and <TT>Float_t</TT>&nbsp; C array.
43 <BR>&nbsp;
44 
45 <H3>
46 Access to the components</H3>
47 There are two sets of access functions to the components of a<TT> LorentzVector</TT>:
48 <TT>X(</TT>), <TT>Y()</TT>, <TT>Z()</TT>, <TT>T()</TT> and <TT>Px()</TT>,<TT>
49 Py()</TT>, <TT>Pz()</TT> and <TT>E()</TT>. Both sets return the same values
50 but the first set is more relevant for use where <TT>TLorentzVector</TT>
51 describes a combination of position and time and the second set is more
52 relevant where <TT>TLorentzVector</TT> describes momentum and energy:
53 
54 <P><TT>&nbsp; Double_t xx =v.X();</TT>
55 <BR><TT>&nbsp; ...</TT>
56 <BR><TT>&nbsp; Double_t tt = v.T();</TT>
57 
58 <P><TT>&nbsp; Double_t px = v.Px();</TT>
59 <BR><TT>&nbsp; ...</TT>
60 <BR><TT>&nbsp; Double_t ee = v.E();</TT>
61 
62 <P>The components of TLorentzVector can also accessed by index:
63 
64 <P><TT>&nbsp; xx = v(0);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp;&nbsp;
65 xx = v[0];</TT>
66 <BR><TT>&nbsp; yy = v(1);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
67 yy = v[1];</TT>
68 <BR><TT>&nbsp; zz = v(2);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
69 zz = v[2];</TT>
70 <BR><TT>&nbsp; tt = v(3);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
71 tt = v[3];</TT>
72 
73 <P>You can use the <TT>Vect()</TT> member function to get the vector component
74 of <TT>TLorentzVector</TT>:
75 
76 <P>&nbsp; <TT>TVector3 p = v.Vect();</TT>
77 
78 <P>For setting components also two sets of member functions can be used:
79 SetX(),.., SetPx(),..:
80 <BR>&nbsp;
81 <BR><TT>&nbsp; v.SetX(1.);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp;
82 v.SetPx(1.);</TT>
83 <BR><TT>&nbsp; ...&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
84 ...</TT>
85 <BR><TT>&nbsp; v.SetT(1.);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
86 v.SetE(1.);</TT>
87 
88 <P>To set more the one component by one call you can use the <TT>SetVect()</TT>
89 function for the <TT>TVector3</TT> part or <TT>SetXYZT()</TT>, <TT>SetPxPyPzE()</TT>. For convenience there is
90 
91 also a <TT>SetXYZM()</TT>:
92 
93 <P><TT>&nbsp; v.SetVect(TVector3(1,2,3));</TT>
94 <BR><TT>&nbsp; v.SetXYZT(x,y,z,t);</TT>
95 <BR><TT>&nbsp; v.SetPxPyPzE(px,py,pz,e);</TT>
96 <BR><TT>&nbsp; v.SetXYZM(x,y,z,m);&nbsp;&nbsp; //&nbsp;&nbsp; ->&nbsp;
97 v=(x,y,z,e=Sqrt(x*x+y*y+z*z+m*m))</TT>
98 <H3>
99 Vector components in noncartesian coordinate systems</H3>
100 There are a couple of memberfunctions to get and set the <TT>TVector3</TT>
101 part of the parameters in
102 <BR><B>spherical</B> coordinate systems:
103 
104 <P><TT>&nbsp; Double_t m, theta, cost, phi, pp, pp2, ppv2, pp2v2;</TT>
105 <BR><TT>&nbsp; m = v.Rho();</TT>
106 <BR><TT>&nbsp; t = v.Theta();</TT>
107 <BR><TT>&nbsp; cost = v.CosTheta();</TT>
108 <BR><TT>&nbsp; phi = v.Phi();</TT>
109 
110 <P><TT>&nbsp; v.SetRho(10.);</TT>
111 <BR><TT>&nbsp; v.SetTheta(TMath::Pi()*.3);</TT>
112 <BR><TT>&nbsp; v.SetPhi(TMath::Pi());</TT>
113 
114 <P>or get information about the r-coordinate in <B>cylindrical</B> systems:
115 
116 <P><TT>&nbsp; Double_t pp, pp2, ppv2, pp2v2;</TT>
117 <BR><TT>&nbsp; pp = v.Perp();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// get transvers component</TT>
118 <BR><TT>&nbsp; pp2 = v.Perp2();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// get transverse component squared</TT>
119 <BR><TT>&nbsp; ppv2 = v.Perp(v1);&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // get
120 transvers component with</TT>
121 <BR><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// respect to another vector</TT>
122 <BR><TT>&nbsp; pp2v2 = v.Perp(v1);</TT>
123 
124 <P>for convenience there are two more set functions <TT>SetPtEtaPhiE(pt,eta,phi,e);
125 and SetPtEtaPhiM(pt,eta,phi,m);</TT>
126 <H3>
127 Arithmetic and comparison operators</H3>
128 The <TT>TLorentzVecto</TT>r class provides operators to add, subtract or
129 compare four-vectors:
130 
131 <P><TT>&nbsp; v3 = -v1;</TT>
132 <BR><TT>&nbsp; v1 = v2+v3;</TT>
133 <BR><TT>&nbsp; v1+= v3;</TT>
134 <BR><TT>&nbsp; v1 = v2 + v3;</TT>
135 <BR><TT>&nbsp; v1-= v3;</TT>
136 
137 <P><TT>&nbsp; if (v1 == v2) {...}</TT>
138 <BR><TT>&nbsp; if(v1 != v3) {...}</TT>
139 <H3>
140 Magnitude/Invariant mass, beta, gamma, scalar product</H3>
141 The scalar product of two four-vectors is calculated with the (-,-,-,+)
142 metric,
143 <BR>&nbsp;&nbsp; i.e.&nbsp;&nbsp; <TT><B>s = v1*v2 = </B>t1*t2-x1*x2-y1*y2-z1*z2</TT>
144 <BR>The magnitude squared <B><TT>mag2</TT></B> of a four-vector is therfore:
145 <BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; <TT>mag2<B>
146 = v*v = </B>t*t-x*x-y*y-z*z</TT>
147 <BR>It <TT>mag2</TT> is negative <TT>mag = -Sqrt(-mag*mag)</TT>. The member
148 functions are:
149 
150 <P><TT>&nbsp; Double_t s, s2;</TT>
151 <BR><TT>&nbsp; s&nbsp; = v1.Dot(v2);&nbsp;&nbsp;&nbsp;&nbsp; // scalar
152 product</TT>
153 <BR><TT>&nbsp; s&nbsp; = v1*v2;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// scalar product</TT>
154 <BR><TT>&nbsp; s2 = v.Mag2();&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp; s2 = v.M2();</TT>
155 <BR><TT>&nbsp; s&nbsp; = v.Mag();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
156 s&nbsp; = v.M();</TT>
157 
158 <P>Since in case of momentum and energy the magnitude has the meaning of
159 invariant mass <TT>TLorentzVector</TT> provides the more meaningful aliases
160 <TT>M2()</TT> and <TT>M()</TT>;
161 <P>The member functions <TT>Beta()</TT> and <TT>Gamma()</TT> returns
162 <TT>beta</TT> and <tt>gamma = 1/Sqrt(1-beta*beta)</tt>.
163 <H3>
164 Lorentz boost</H3>
165 A boost in a general direction can be parameterized with three parameters
166 which can be taken as the components of a three vector <TT><B>b </B>= (bx,by,bz)</TT>.
167 With
168 <BR><TT><B>&nbsp; x</B> = (x,y,z) and gamma = 1/Sqrt(1-beta*beta)</TT> (beta being the module of vector b),
169 an arbitary active Lorentz boost transformation (from the rod frame
170 to the original frame) can be written as:
171 <BR>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; <TT><B>x</B>
172 = <B>x'</B> + (gamma-1)/(beta*beta) * (<B>b</B>*<B>x'</B>) * <B>b</B> +
173 gamma * t' *<B> b</B></TT>
174 <BR><B>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; </B><TT>t
175 = gamma (t'+ <B>b</B>*<B>x'</B>).</TT>
176 
177 <P>The member function <TT>Boost()</TT> performs a boost transformation
178 from the rod frame to the original frame. <TT>BoostVector()</TT> returns
179 a <TT>TVector3</TT> of the spatial components divided by the time component:
180 
181 <P><TT>&nbsp; TVector3 b;</TT>
182 <BR><TT>&nbsp; v.Boost(bx,by,bz);</TT>
183 <BR><TT>&nbsp; v.Boost(b);</TT>
184 <BR><TT>&nbsp; b = v.BoostVector();&nbsp;&nbsp; // b=(x/t,y/t,z/t)</TT>
185 <H3>
186 Rotations</H3>
187 There are four sets of functions to rotate the <TT>TVector3</TT> component
188 of a <TT>TLorentzVector</TT>:
189 <H5>
190 rotation around axes</H5>
191 <TT>&nbsp; v.RotateX(TMath::Pi()/2.);</TT>
192 <BR><TT>&nbsp; v.RotateY(.5);</TT>
193 <BR><TT>&nbsp; v.RotateZ(.99);</TT>
194 <H5>
195 rotation around an arbitary axis</H5>
196 <TT>&nbsp; v.Rotate(TMath::Pi()/4., v1); // rotation around v1</TT>
197 <H5>
198 transformation from rotated frame</H5>
199 <TT>&nbsp; v.RotateUz(direction); //&nbsp; direction must be a unit TVector3</TT>
200 <H5>
201 by TRotation (see TRotation)</H5>
202 <TT>&nbsp; TRotation r;</TT>
203 <BR><TT>&nbsp; v.Transform(r);&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp;&nbsp;
204 v *= r; // Attention v=M*v</TT>
205 <H3>
206 Misc</H3>
207 
208 <H5>
209 Angle between two vectors</H5>
210 <TT>&nbsp; Double_t a = v1.Angle(v2.Vect());&nbsp; // get angle between v1 and
211 v2</TT>
212 <H5>
213 Light-cone components</H5>
214 Member functions <TT>Plus()</TT> and <TT>Minus(</TT>) return the positive
215 and negative light-cone components:<TT></TT>
216 
217 <P><TT>&nbsp; Double_t pcone = v.Plus();</TT>
218 <BR><TT>&nbsp; Double_t mcone = v.Minus();</TT>
219 <P>CAVEAT: The values returned are T{+,-}Z. It is known that some authors
220 find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus
221 check what definition is used in the physics you're working in and adapt
222 your code accordingly.
223 
224 <H5>
225 Transformation by TLorentzRotation</H5>
226 A general Lorentz transformation see class <TT>TLorentzRotation</TT> can
227 be used by the <TT>Transform()</TT> member function, the <TT>*=</TT> or
228 <TT>*</TT> operator of the <TT>TLorentzRotation</TT> class:<TT></TT>
229 
230 <P><TT>&nbsp; TLorentzRotation l;</TT>
231 <BR><TT>&nbsp; v.Transform(l);</TT>
232 <BR><TT>&nbsp; v = l*v;&nbsp;&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp;&nbsp;
233 v *= l;&nbsp; // Attention v = l*v</TT>
234 
235 <P>
236 <!--*/
237 // -->END_HTML
238 
239 #include "TLorentzVector.h"
240 
241 #include "TBuffer.h"
242 #include "TClass.h"
243 #include "TError.h"
244 #include "TLorentzRotation.h"
245 
247 
249  : fP(), fE(0.0) {}
250 
252  : fP(x,y,z), fE(t) {}
253 
255  : fP(x0), fE(x0[3]) {}
256 
258  : fP(x0), fE(x0[3]) {}
259 
261  : fP(p), fE(e) {}
262 
264  , fP(p.Vect()), fE(p.T()) {}
265 
267 
269 {
270  //dereferencing operator const
271  switch(i) {
272  case kX:
273  case kY:
274  case kZ:
275  return fP(i);
276  case kT:
277  return fE;
278  default:
279  Error("operator()()", "bad index (%d) returning 0",i);
280  }
281  return 0.;
282 }
283 
285 {
286  //dereferencing operator
287  switch(i) {
288  case kX:
289  case kY:
290  case kZ:
291  return fP(i);
292  case kT:
293  return fE;
294  default:
295  Error("operator()()", "bad index (%d) returning &fE",i);
296  }
297  return fE;
298 }
299 
301 {
302  //Boost this Lorentz vector
303  Double_t b2 = bx*bx + by*by + bz*bz;
304  Double_t gamma = 1.0 / TMath::Sqrt(1.0 - b2);
305  Double_t bp = bx*X() + by*Y() + bz*Z();
306  Double_t gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
307 
308  SetX(X() + gamma2*bp*bx + gamma*bx*T());
309  SetY(Y() + gamma2*bp*by + gamma*by*T());
310  SetZ(Z() + gamma2*bp*bz + gamma*bz*T());
311  SetT(gamma*(T() + bp));
312 }
313 
315 {
316  //return rapidity
317  return 0.5*log( (E()+Pz()) / (E()-Pz()) );
318 }
319 
321 {
322  //multiply this Lorentzvector by m
323  return *this = m.VectorMultiplication(*this);
324 }
325 
327 {
328  //Transform this Lorentzvector
329  return *this = m.VectorMultiplication(*this);
330 }
331 
332 void TLorentzVector::Streamer(TBuffer &R__b)
333 {
334  // Stream an object of class TLorentzVector.
335  Double_t x, y, z;
336  UInt_t R__s, R__c;
337  if (R__b.IsReading()) {
338  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
339  if (R__v > 3) {
340  R__b.ReadClassBuffer(TLorentzVector::Class(), this, R__v, R__s, R__c);
341  return;
342  }
343  //====process old versions before automatic schema evolution
344  if (R__v != 2) TObject::Streamer(R__b);
345  R__b >> x;
346  R__b >> y;
347  R__b >> z;
348  fP.SetXYZ(x,y,z);
349  R__b >> fE;
350  R__b.CheckByteCount(R__s, R__c, TLorentzVector::IsA());
351  } else {
353  }
354 }
355 
356 
357 ////////////////////////////////////////////////////////////////////////////////
358 /// Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations
359 
361 {
362  Printf("(x,y,z,t)=(%f,%f,%f,%f) (P,eta,phi,E)=(%f,%f,%f,%f)",
363  fP.x(),fP.y(),fP.z(),fE,
364  P(),Eta(),Phi(),fE);
365 }
Double_t y() const
Definition: TVector3.h:218
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
void Boost(Double_t, Double_t, Double_t)
Double_t z() const
Definition: TVector3.h:219
Bool_t IsReading() const
Definition: TBuffer.h:83
short Version_t
Definition: RtypesCore.h:61
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Double_t Z() const
Double_t P() const
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Double_t Y() const
Double_t X() const
TLorentzVector & operator*=(Double_t a)
TTree * T
Double_t x[n]
Definition: legend1.C:17
virtual ~TLorentzVector()
void Class()
Definition: Class.C:29
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:231
virtual void Print(Option_t *option="") const
Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Float_t z[5]
Definition: Ifit.C:16
TLorentzVector VectorMultiplication(const TLorentzVector &) const
double gamma(double x)
TThread * t[5]
Definition: threadsh1.C:13
Double_t Eta() const
TClass * IsA() const
unsigned int UInt_t
Definition: RtypesCore.h:42
Double_t operator()(int i) const
TMarker * m
Definition: textangle.C:8
void SetX(Double_t a)
Double_t Rapidity() const
Double_t x() const
Definition: TVector3.h:217
#define Printf
Definition: TGeoToOCC.h:18
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
ClassImp(TLorentzVector) TLorentzVector
double Double_t
Definition: RtypesCore.h:55
Double_t E() const
void SetT(Double_t a)
Double_t y[n]
Definition: legend1.C:17
void SetY(Double_t a)
Double_t Phi() const
Mother of all ROOT objects.
Definition: TObject.h:58
void SetZ(Double_t a)
Double_t Pz() const
TLorentzVector & Transform(const TRotation &)
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Double_t T() const
double log(double)
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0