Logo ROOT   6.12/07
Reference Guide
TParticle.cxx
Go to the documentation of this file.
1 // @(#)root/eg:$Id$
2 // Author: Rene Brun , Federico Carminati 26/04/99
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TParticle
13  \ingroup eg
14 
15 Description of the dynamic properties of a particle.
16 
17 A dynamic particle class created by event generators and used during
18 the propagation in detectors. The static attributes of a TParticle
19 are described by TParticlePDG.
20 
21 \verbatim
22  Int_t fPdgCode; // PDG code of the particle
23  Int_t fStatusCode; // generation status code
24  Int_t fMother[2]; // Indices of the mother particles
25  Int_t fDaughter[2]; // Indices of the daughter particles
26  Float_t fWeight; // particle weight
27 
28  Double_t fCalcMass; // Calculated mass
29 
30  Double_t fPx; // x component of momentum
31  Double_t fPy; // y component of momentum
32  Double_t fPz; // z component of momentum
33  Double_t fE; // Energy
34 
35  Double_t fVx; // x of production vertex
36  Double_t fVy; // y of production vertex
37  Double_t fVz; // z of production vertex
38  Double_t fVt; // t of production vertex
39 
40  Double_t fPolarTheta; // Polar angle of polarisation
41  Double_t fPolarPhi; // azymutal angle of polarisation
42 
43  TParticlePDG* fParticlePDG; //! reference to the particle record in PDG database
44 \endverbatim
45 */
46 
47 #include "TParticle.h"
48 
49 #include "TBuffer.h"
50 #include "TView.h"
51 #include "TVirtualPad.h"
52 #include "TPolyLine3D.h"
53 #include "TParticlePDG.h"
54 #include "TDatabasePDG.h"
55 #include "TClass.h"
56 #include "X3DBuffer.h"
57 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 ///default constructor
62 
64  fPdgCode(0), fStatusCode(0), fWeight(0),fCalcMass(0), fPx(0), fPy(0),
65  fPz(0), fE(0), fVx(0), fVy(0), fVz(0), fVt(0), fPolarTheta(0), fPolarPhi(0)
66 {
67  fMother[0] = 0;
68  fMother[1] = 0;
69  fDaughter[0] = 0;
70  fDaughter[1] = 0;
71  fParticlePDG = 0;
72 }
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 ///constructor
76 
78  Int_t mother1, Int_t mother2,
79  Int_t daughter1, Int_t daughter2,
80  Double_t px, Double_t py, Double_t pz, Double_t etot,
81  Double_t vx, Double_t vy, Double_t vz, Double_t time):
82  fPdgCode(pdg), fStatusCode(status), fWeight(1.),fPx(px), fPy(py),
83  fPz(pz), fE(etot), fVx(vx), fVy(vy), fVz(vz), fVt(time)
84 {
85  fMother[0] = mother1;
86  fMother[1] = mother2;
87  fDaughter[0] = daughter1;
88  fDaughter[1] = daughter2;
89 
90  SetPolarisation(0,0,0);
91 
92  SetPdgCode(pdg);
93 }
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 ///constructor
97 
99  Int_t mother1, Int_t mother2,
100  Int_t daughter1, Int_t daughter2,
101  const TLorentzVector &p,
102  const TLorentzVector &v) :
103  fPdgCode(pdg), fStatusCode(status), fWeight(1.),fPx(p.Px()), fPy(p.Py()),
104  fPz(p.Pz()), fE(p.E()), fVx(v.X()), fVy(v.Y()), fVz(v.Z()), fVt(v.T())
105 {
106  fMother[0] = mother1;
107  fMother[1] = mother2;
108  fDaughter[0] = daughter1;
109  fDaughter[1] = daughter2;
110 
111  SetPolarisation(0,0,0);
112 
113  SetPdgCode(pdg);
114 }
115 
116 ////////////////////////////////////////////////////////////////////////////////
117 /// copy constructor
118 
121  fWeight(p.fWeight), fCalcMass(p.fCalcMass), fPx(p.fPx), fPy(p.fPy), fPz(p.fPz),
122  fE(p.fE), fVx(p.fVx), fVy(p.fVy), fVz(p.fVz), fVt(p.fVt), fPolarTheta(p.fPolarTheta),
124 {
125  fMother[0]=p.fMother[0];
126  fMother[1]=p.fMother[1];
127  fDaughter[0]=p.fDaughter[0];
128  fDaughter[1]=p.fDaughter[1];
129 }
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// Equal operator
133 
135 {
136  if(this!=&p) {
140  fPdgCode=p.fPdgCode;
142  fMother[0]=p.fMother[0];
143  fMother[1]=p.fMother[1];
144  fDaughter[0]=p.fDaughter[0];
145  fDaughter[1]=p.fDaughter[1];
146  fWeight=p.fWeight;
147 
149 
150  fPx=p.fPx;
151  fPy=p.fPy;
152  fPz=p.fPz;
153  fE=p.fE;
154 
155  fVx=p.fVx;
156  fVy=p.fVy;
157  fVz=p.fVz;
158  fVt=p.fVt;
159 
162 
164  }
165  return *this;
166 }
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 ///destructor
170 
172 {
173 }
174 
175 ////////////////////////////////////////////////////////////////////////////////
176 /// Return nominal particle mass from PDG table.
177 
179 {
180  return GetPDG()->Mass();
181 }
182 
183 ////////////////////////////////////////////////////////////////////////////////
184 /// Return beauty quantum number.
185 
187 {
188  return GetPDG()->Beauty();
189 }
190 
191 ////////////////////////////////////////////////////////////////////////////////
192 /// Return charm quantum number.
193 
195 {
196  return GetPDG()->Charm();
197 }
198 
199 ////////////////////////////////////////////////////////////////////////////////
200 /// Return strangeness quantum number.
201 
203 {
204  return GetPDG()->Strangeness();
205 }
206 
207 ////////////////////////////////////////////////////////////////////////////////
208 /// Compute distance from point px,py to a primary track
209 ///
210 /// Compute the closest distance of approach from point px,py to each segment
211 /// of a track.
212 /// The distance is computed in pixels units.
213 ///ß
214 
216 {
217  const Int_t big = 9999;
218  Float_t xv[3], xe[3], xndc[3];
219  Float_t rmin[3], rmax[3];
220  TView *view = gPad->GetView();
221  if(!view) return big;
222 
223  // compute first and last point in pad coordinates
224  Float_t pmom = this->P();
225  if (pmom == 0) return big;
226  view->GetRange(rmin,rmax);
227  Float_t rbox = rmax[2];
228  xv[0] = fVx;
229  xv[1] = fVy;
230  xv[2] = fVz;
231  xe[0] = xv[0]+rbox*fPx/pmom;
232  xe[1] = xv[1]+rbox*fPy/pmom;
233  xe[2] = xv[2]+rbox*fPz/pmom;
234  view->WCtoNDC(xv, xndc);
235  Float_t x1 = xndc[0];
236  Float_t y1 = xndc[1];
237  view->WCtoNDC(xe, xndc);
238  Float_t x2 = xndc[0];
239  Float_t y2 = xndc[1];
240 
241  return DistancetoLine(px,py,x1,y1,x2,y2);
242 }
243 
244 
245 ////////////////////////////////////////////////////////////////////////////////
246 /// Execute action corresponding to one event
247 ///
248 
250 {
251  gPad->SetCursor(kPointer);
252 }
253 
254 ////////////////////////////////////////////////////////////////////////////////
255 /// Return particle name
256 
257 const char* TParticle::GetName() const {
258  static char def[4] = "XXX";
260  if (ap) return ap->GetName();
261  else return def;
262 }
263 
264 
265 ////////////////////////////////////////////////////////////////////////////////
266 /// Returns a pointer to the TParticlePDG object using the pdgcode.
267 /// - if mode == 0 (default) always get a fresh value for the pointer.
268 /// - if mode != 0 this function returns directly the previously
269 /// computed pointer from a previous call
270 /// One can use mode=1 (faster) when the TParticle object is not part of a
271 /// TClonesArray used in split mode in a Root TTree.
272 
274 {
275  if (!mode || !fParticlePDG) {
277  }
278  return fParticlePDG;
279 }
280 
281 ////////////////////////////////////////////////////////////////////////////////
282 /// Return particle polarisation
283 
285 {
286  if(fPolarTheta == -99 && fPolarPhi == -99)
287  //No polarisation to return
288  v.SetXYZ(0.,0.,0.);
289  else
293 }
294 
295 ////////////////////////////////////////////////////////////////////////////////
296 /// Return particle title
297 
298 const char *TParticle::GetTitle() const
299 {
300  static char def[4] = "XXX";
302  if (ap) return ap->GetTitle();
303  else return def;
304 }
305 
306 ////////////////////////////////////////////////////////////////////////////////
307 ///
308 /// Paint a primary track
309 ///
310 
312 {
313  Float_t rmin[3], rmax[3];
314  static TPolyLine3D *pline = 0;
315  if (!pline) {
316  pline = new TPolyLine3D(2);
317  }
318  Float_t pmom = this->P();
319  if (pmom == 0) return;
320  TView *view = gPad->GetView();
321  if (!view) return;
322  view->GetRange(rmin,rmax);
323  Float_t rbox = rmax[2];
324  pline->SetPoint(0,Vx(), Vy(), Vz());
325  Float_t xend = Vx()+rbox*Px()/pmom;
326  Float_t yend = Vy()+rbox*Py()/pmom;
327  Float_t zend = Vz()+rbox*Pz()/pmom;
328  pline->SetPoint(1, xend, yend, zend);
329  pline->SetLineColor(GetLineColor());
330  pline->SetLineStyle(GetLineStyle());
331  pline->SetLineWidth(GetLineWidth());
332  pline->Paint(option);
333 }
334 
335 ////////////////////////////////////////////////////////////////////////////////
336 ///
337 /// Print the internals of the primary vertex particle
338 ///
339 
341 {
342  Printf("TParticle: %-13s p: %8f %8f %8f Vertex: %8e %8e %8e %5d %5d",
343  GetName(),Px(),Py(),Pz(),Vx(),Vy(),Vz(),
344  fMother[0],fMother[1]);
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Change the PDG code for this particle
349 ///
350 /// Get a new pointer to a TParticlePDG from TDatabasePDG.
351 /// Recompute the mass.
352 
354 {
355  static Int_t nWarnings = 0;
356  fPdgCode = pdg;
358  if (fParticlePDG) {
360  } else {
361  if (nWarnings < 10) {
362  Warning("SetPdgCode","PDG code %d unknown from TDatabasePDG",pdg);
363  nWarnings++;
364  }
365  Double_t a2 = fE*fE -fPx*fPx -fPy*fPy -fPz*fPz;
366  if (a2 >= 0) fCalcMass = TMath::Sqrt(a2);
367  else fCalcMass = -TMath::Sqrt(-a2);
368  }
369 }
370 
371 ////////////////////////////////////////////////////////////////////////////////
372 /// Set particle polarisation
373 
375 {
376  if(polx || poly || polz) {
377  fPolarTheta = TMath::ACos(polz/TMath::Sqrt(polx*polx+poly*poly+polz*polz));
378  fPolarPhi = TMath::Pi()+TMath::ATan2(-poly,-polx);
379  } else {
380  fPolarTheta = -99;
381  fPolarPhi = -99;
382  }
383 }
384 
385 ////////////////////////////////////////////////////////////////////////////////
386 /// Return total X3D size of this primary
387 ///
388 
390 {
391  Float_t pmom = this->P();
392  if (pmom == 0) return;
393  Int_t npoints = 2;
394  gSize3D.numPoints += npoints;
395  gSize3D.numSegs += (npoints-1);
396  gSize3D.numPolys += 0;
397 
398 }
399 
400 ////////////////////////////////////////////////////////////////////////////////
401 /// Stream an object of class TParticle.
402 
403 void TParticle::Streamer(TBuffer &R__b)
404 {
405  if (R__b.IsReading()) {
406  UInt_t R__s, R__c;
407  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
408  if (R__v > 1) {
409  R__b.ReadClassBuffer(TParticle::Class(), this, R__v, R__s, R__c);
411  return;
412  }
413  //====process old versions before automatic schema evolution
414  TObject::Streamer(R__b);
415  TAttLine::Streamer(R__b);
416  R__b >> fPdgCode;
417  R__b >> fStatusCode;
418  R__b.ReadStaticArray(fMother);
420  R__b >> fWeight;
421  R__b >> fCalcMass;
422  R__b >> fPx;
423  R__b >> fPy;
424  R__b >> fPz;
425  R__b >> fE;
426  R__b >> fVx;
427  R__b >> fVy;
428  R__b >> fVz;
429  R__b >> fVt;
430  R__b >> fPolarTheta;
431  R__b >> fPolarPhi;
433  R__b.CheckByteCount(R__s, R__c, TParticle::IsA());
434  //====end of old versions
435 
436  } else {
437  R__b.WriteClassBuffer(TParticle::Class(),this);
438  }
439 }
Double_t Py() const
Definition: TParticle.h:132
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Int_t Beauty() const
Return beauty quantum number.
Definition: TParticle.cxx:186
Double_t Y() const
Definition: TParticle.h:143
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
Bool_t IsReading() const
Definition: TBuffer.h:83
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a primary track.
Definition: TParticle.cxx:215
short Version_t
Definition: RtypesCore.h:61
Int_t Beauty() const
Definition: TParticlePDG.h:77
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Double_t GetMass() const
Return nominal particle mass from PDG table.
Definition: TParticle.cxx:178
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
Double_t fPolarPhi
Definition: TParticle.h:50
See TView3D.
Definition: TView.h:25
Use this attribute class when an object should have 3D capabilities.
Definition: TAtt3D.h:19
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
virtual ~TParticle()
destructor
Definition: TParticle.cxx:171
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
int Int_t
Definition: RtypesCore.h:41
A 3-dimensional polyline.
Definition: TPolyLine3D.h:31
virtual Int_t ReadStaticArray(Bool_t *b)=0
Double_t fPolarTheta
Definition: TParticle.h:49
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
static const double x2[5]
static TDatabasePDG * Instance()
static function
void Class()
Definition: Class.C:29
virtual const char * GetName() const
Return particle name.
Definition: TParticle.cxx:257
virtual void Print(Option_t *option="") const
Print the internals of the primary vertex particle.
Definition: TParticle.cxx:340
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:227
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.h:271
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:580
Int_t fStatusCode
Definition: TParticle.h:32
constexpr Double_t Pi()
Definition: TMath.h:40
Int_t Strangeness() const
Return strangeness quantum number.
Definition: TParticle.cxx:202
TVector3 is a general three vector class, which can be used for the description of different vectors ...
Definition: TVector3.h:22
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
Float_t fWeight
Definition: TParticle.h:35
Double_t P() const
Definition: TParticle.h:134
Double_t fVy
Definition: TParticle.h:45
Double_t fPz
Definition: TParticle.h:41
TLorentzVector is a general four-vector class, which can be used either for the description of positi...
#define gSize3D
Definition: X3DBuffer.h:40
virtual const char * GetTitle() const
Return particle title.
Definition: TParticle.cxx:298
TParticle & operator=(const TParticle &)
Equal operator.
Definition: TParticle.cxx:134
SVector< double, 2 > v
Definition: Dict.h:5
Int_t fDaughter[2]
Definition: TParticle.h:34
Int_t fPdgCode
Definition: TParticle.h:31
Double_t fVt
Definition: TParticle.h:47
Double_t Vx() const
Definition: TParticle.h:125
Int_t Strangeness() const
Definition: TParticlePDG.h:75
unsigned int UInt_t
Definition: RtypesCore.h:42
Description of the static properties of a particle.
Definition: TParticlePDG.h:19
Int_t Charm() const
Return charm quantum number.
Definition: TParticle.cxx:194
Double_t ACos(Double_t)
Definition: TMath.h:571
Int_t Charm() const
Definition: TParticlePDG.h:76
constexpr Double_t E()
Definition: TMath.h:74
Int_t fMother[2]
Definition: TParticle.h:33
#define Printf
Definition: TGeoToOCC.h:18
Double_t Cos(Double_t)
Definition: TMath.h:550
Double_t Px() const
Definition: TParticle.h:131
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
Double_t T() const
Definition: TParticle.h:128
TParticle()
reference to the particle record in PDG database
Definition: TParticle.cxx:63
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Double_t fPx
Definition: TParticle.h:39
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
Double_t fE
Definition: TParticle.h:42
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
TParticlePDG * GetPDG(Int_t mode=0) const
Returns a pointer to the TParticlePDG object using the pdgcode.
Definition: TParticle.cxx:273
virtual void Paint(Option_t *option="")
Paint a TPolyLine3D.
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition: TAttLine.cxx:196
virtual void Paint(Option_t *option="")
Paint a primary track.
Definition: TParticle.cxx:311
Double_t Vz() const
Definition: TParticle.h:127
Double_t Vy() const
Definition: TParticle.h:126
TParticlePDG * GetParticle(Int_t pdgCode) const
Get a pointer to the particle object according to the MC code number.
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:42
Binding & operator=(OUT(*fun)(void))
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Sizeof3D() const
Return total X3D size of this primary.
Definition: TParticle.cxx:389
void SetPdgCode(Int_t pdg)
Change the PDG code for this particle.
Definition: TParticle.cxx:353
TParticlePDG * fParticlePDG
Definition: TParticle.h:52
void SetPolarisation(Double_t theta, Double_t phi)
Definition: TParticle.h:111
virtual void GetRange(Float_t *min, Float_t *max)=0
Double_t Sin(Double_t)
Definition: TMath.h:547
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TParticle.cxx:249
#define gPad
Definition: TVirtualPad.h:285
TPolyLine * pline
Definition: polyline.C:4
void GetPolarisation(TVector3 &v) const
Return particle polarisation.
Definition: TParticle.cxx:284
Double_t Mass() const
Definition: TParticlePDG.h:67
Double_t fVz
Definition: TParticle.h:46
Double_t Sqrt(Double_t x)
Definition: TMath.h:590
Double_t fVx
Definition: TParticle.h:44
Double_t Pz() const
Definition: TParticle.h:133
Line Attributes class.
Definition: TAttLine.h:18
Double_t fCalcMass
Definition: TParticle.h:37
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
Double_t fPy
Definition: TParticle.h:40