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