Logo ROOT   6.10/09
Reference Guide
TParticle.h
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 // TParticle: defines equivalent of HEPEVT particle //
14 //////////////////////////////////////////////////////////////////////////
15 
16 #ifndef ROOT_TParticle
17 #define ROOT_TParticle
18 
19 #include "TNamed.h"
20 #include "TAttLine.h"
21 #include "TAtt3D.h"
22 #include "TLorentzVector.h"
23 
24 class TParticlePDG;
25 
26 class TParticle : public TObject, public TAttLine, public TAtt3D {
27 
28 
29 protected:
30 
31  Int_t fPdgCode; // PDG code of the particle
32  Int_t fStatusCode; // generation status code
33  Int_t fMother[2]; // Indices of the mother particles
34  Int_t fDaughter[2]; // Indices of the daughter particles
35  Float_t fWeight; // particle weight
36 
37  Double_t fCalcMass; // Calculated mass
38 
39  Double_t fPx; // x component of momentum
40  Double_t fPy; // y component of momentum
41  Double_t fPz; // z component of momentum
42  Double_t fE; // Energy
43 
44  Double_t fVx; // x of production vertex
45  Double_t fVy; // y of production vertex
46  Double_t fVz; // z of production vertex
47  Double_t fVt; // t of production vertex
48 
49  Double_t fPolarTheta; // Polar angle of polarisation
50  Double_t fPolarPhi; // azymutal angle of polarisation
51 
52  mutable TParticlePDG* fParticlePDG; //! reference to the particle record in PDG database
53  //----------------------------------------------------------------------------
54  // functions
55  //----------------------------------------------------------------------------
56 public:
57  // ****** constructors and destructor
58  TParticle();
59 
60  TParticle(Int_t pdg, Int_t status,
61  Int_t mother1, Int_t mother2,
62  Int_t daughter1, Int_t daughter2,
63  Double_t px, Double_t py, Double_t pz, Double_t etot,
64  Double_t vx, Double_t vy, Double_t vz, Double_t time);
65 
66  TParticle(Int_t pdg, Int_t status,
67  Int_t mother1, Int_t mother2,
68  Int_t daughter1, Int_t daughter2,
69  const TLorentzVector &p,
70  const TLorentzVector &v);
71 
72  TParticle(const TParticle &part);
73 
74  virtual ~TParticle();
75 
76  TParticle& operator=(const TParticle&);
77 
78 // virtual TString* Name () const { return fName.Data(); }
79 // virtual char* GetName() const { return fName.Data(); }
80 
81  Double_t Ek () const { return fE-fCalcMass; }
82  Int_t GetStatusCode () const { return fStatusCode; }
83  Int_t GetPdgCode () const { return fPdgCode; }
84  Int_t GetFirstMother () const { return fMother[0]; }
85  Int_t GetMother (Int_t i) const { return fMother[i]; }
86  Int_t GetSecondMother () const { return fMother[1]; }
87  Bool_t IsPrimary () const { return fMother[0]>-1 ? kFALSE : kTRUE; } //Is this particle primary one?
88  Int_t GetFirstDaughter() const { return fDaughter[0]; }
89  Int_t GetDaughter (Int_t i) const { return fDaughter[i]; }
90  Int_t GetLastDaughter () const { return fDaughter[1]; }
91  Double_t GetCalcMass () const { return fCalcMass; }
92  Double_t GetMass () const;
93  Int_t GetNDaughters () const { return fDaughter[1]>0 ? fDaughter[1]-fDaughter[0]+1 : 0;}
94  Float_t GetWeight () const { return fWeight; }
95  void GetPolarisation(TVector3 &v) const;
96  TParticlePDG* GetPDG (Int_t mode = 0) const;
97  Int_t Beauty () const;
98  Int_t Charm () const;
99  Int_t Strangeness () const;
100  Double_t PhiX () const { return TMath::Pi()+TMath::ATan2(-fPz,-fPy); } // note that PhiX() returns an angle between 0 and 2pi
101  Double_t PhiY () const { return TMath::Pi()+TMath::ATan2(-fPx,-fPz); } // note that PhiY() returns an angle between 0 and 2pi
102  Double_t PhiZ () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); } // note that PhiZ() returns an angle between 0 and 2pi
103  Double_t ThetaX () const { return (fPx==0)?TMath::PiOver2():TMath::ACos(fPx/P()); }
104  Double_t ThetaY () const { return (fPy==0)?TMath::PiOver2():TMath::ACos(fPy/P()); }
105  Double_t ThetaZ () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
106  Double_t GetPolarTheta () const { return fPolarTheta; }
107  Double_t GetPolarPhi () const { return fPolarPhi; }
108  void GetPolarisation (Double_t &theta, Double_t &phi) const { theta = fPolarTheta; phi = fPolarPhi; }
109  void SetPolarTheta (Double_t theta) { fPolarTheta = theta; }
110  void SetPolarPhi (Double_t phi) { fPolarPhi = phi; }
111  void SetPolarisation (Double_t theta, Double_t phi) { fPolarTheta = theta; fPolarPhi = phi; }
112  void Momentum(TLorentzVector &v) const { v.SetPxPyPzE(fPx,fPy,fPz,fE); }
113  void ProductionVertex(TLorentzVector &v) const { v.SetXYZT(fVx,fVy,fVz,fVt); }
114 
115  Double_t Theta(const TParticle &p) // Returns the angle between momenta of particles
116  {
117  Double_t v = P()*p.P();
118  if (v == 0) v = 1; else v = (fPx*p.Px()+fPy*p.Py()+fPz*p.Pz())/v;
119  if (v > 1) v = 1; else if (v < -1) v = -1; // just a precaution
120  return TMath::ACos(v);
121  }
122 
123 // ****** redefine several most oftenly used methods of LORENTZ_VECTOR
124 
125  Double_t Vx () const { return fVx; }
126  Double_t Vy () const { return fVy; }
127  Double_t Vz () const { return fVz; }
128  Double_t T () const { return fVt; }
129  Double_t R () const { return TMath::Sqrt(fVx*fVx+fVy*fVy); } //Radius of production vertex in cylindrical system
130  Double_t Rho () const { return TMath::Sqrt(fVx*fVx+fVy*fVy+fVz*fVz); } //Radius of production vertex in spherical system
131  Double_t Px () const { return fPx; }
132  Double_t Py () const { return fPy; }
133  Double_t Pz () const { return fPz; }
134  Double_t P () const { return TMath::Sqrt(fPx*fPx+fPy*fPy+fPz*fPz); }
135  Double_t Pt () const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
136  Double_t Energy () const { return fE; }
137  Double_t Eta () const
138  {
139  Double_t pmom = P();
140  if (pmom != TMath::Abs(fPz)) return 0.5*TMath::Log((pmom+fPz)/(pmom-fPz));
141  else return 1.e30;
142  }
143  Double_t Y () const
144  {
145  if (fE != TMath::Abs(fPz)) return 0.5*TMath::Log((fE+fPz)/(fE-fPz));
146  else return 1.e30;
147  }
148 
149  Double_t Phi () const { return TMath::Pi()+TMath::ATan2(-fPy,-fPx); } // note that Phi() returns an angle between 0 and 2pi
150  Double_t Theta () const { return (fPz==0)?TMath::PiOver2():TMath::ACos(fPz/P()); }
151 
152  // setters
153 
154  void SetFirstMother (int code) { fMother[0] = code ; }
155  void SetMother (int i, int code) { fMother[i] = code ; }
156  void SetLastMother (int code) { fMother[1] = code ; }
157  void SetFirstDaughter(int code) { fDaughter[0] = code ; }
158  void SetDaughter(int i, int code) { fDaughter[i] = code ; }
159  void SetLastDaughter(int code) { fDaughter[1] = code ; }
160  void SetCalcMass(Double_t mass) { fCalcMass=mass;}
161  void SetPdgCode(Int_t pdg);
162  void SetPolarisation(Double_t polx, Double_t poly, Double_t polz);
163  void SetPolarisation(const TVector3& v) {SetPolarisation(v.X(), v.Y(), v.Z());}
164  void SetStatusCode(int status) {fStatusCode = status;}
165  void SetWeight(Float_t weight = 1) {fWeight = weight; }
166  void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e) {fPx=px; fPy=py; fPz=pz; fE=e;}
167  void SetMomentum(const TLorentzVector& p) {SetMomentum(p.Px(),p.Py(),p.Pz(),p.Energy());}
168  void SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t) {fVx=vx; fVy=vy; fVz=vz; fVt=t;}
169  void SetProductionVertex(const TLorentzVector& v) {SetProductionVertex(v.X(),v.Y(),v.Z(),v.T());}
170 
171  // ****** overloaded functions of TObject
172 
173  virtual void Paint(Option_t *option = "");
174  virtual void Print(Option_t *option = "") const;
175  virtual void Sizeof3D() const;
176  virtual Int_t DistancetoPrimitive(Int_t px, Int_t py);
177  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py);
178  virtual const char *GetName() const;
179  virtual const char *GetTitle() const;
180 
181 
182  ClassDef(TParticle,2) // TParticle vertex particle information
183 };
184 
185 #endif
Double_t Py() const
Definition: TParticle.h:132
Double_t X() const
Definition: TVector3.h:216
void SetPolarisation(const TVector3 &v)
Definition: TParticle.h:163
Int_t Beauty() const
Return beauty quantum number.
Definition: TParticle.cxx:186
Double_t Y() const
Definition: TParticle.h:143
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
void SetDaughter(int i, int code)
Definition: TParticle.h:158
Int_t GetFirstMother() const
Definition: TParticle.h:84
Double_t Z() const
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a primary track.
Definition: TParticle.cxx:215
void SetProductionVertex(Double_t vx, Double_t vy, Double_t vz, Double_t t)
Definition: TParticle.h:168
Double_t ThetaY() const
Definition: TParticle.h:104
Double_t Log(Double_t x)
Definition: TMath.h:649
Double_t GetPolarPhi() const
Definition: TParticle.h:107
Double_t R() const
Definition: TParticle.h:129
Double_t GetCalcMass() const
Definition: TParticle.h:91
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Bool_t IsPrimary() const
Definition: TParticle.h:87
Double_t GetMass() const
Return nominal particle mass from PDG table.
Definition: TParticle.cxx:178
Double_t Phi() const
Definition: TParticle.h:149
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
void SetStatusCode(int status)
Definition: TParticle.h:164
Double_t fPolarPhi
Definition: TParticle.h:50
Use this attribute class when an object should have 3D capabilities.
Definition: TAtt3D.h:19
virtual ~TParticle()
destructor
Definition: TParticle.cxx:171
void ProductionVertex(TLorentzVector &v) const
Definition: TParticle.h:113
Int_t GetSecondMother() const
Definition: TParticle.h:86
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t ThetaX() const
Definition: TParticle.h:103
Int_t GetStatusCode() const
Definition: TParticle.h:82
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
Double_t fPolarTheta
Definition: TParticle.h:49
void SetMother(int i, int code)
Definition: TParticle.h:155
Double_t Energy() const
Definition: TParticle.h:136
Double_t Rho() const
Definition: TParticle.h:130
void SetWeight(Float_t weight=1)
Definition: TParticle.h:165
Float_t GetWeight() const
Definition: TParticle.h:94
Double_t Pz() const
#define ClassDef(name, id)
Definition: Rtypes.h:297
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 Momentum(TLorentzVector &v) const
Definition: TParticle.h:112
Double_t Y() const
Definition: TVector3.h:217
void SetPolarPhi(Double_t phi)
Definition: TParticle.h:110
Int_t GetDaughter(Int_t i) const
Definition: TParticle.h:89
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:581
Int_t fStatusCode
Definition: TParticle.h:32
constexpr Double_t Pi()
Definition: TMath.h:40
Int_t GetNDaughters() const
Definition: TParticle.h:93
Double_t Eta() const
Definition: TParticle.h:137
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
Int_t GetMother(Int_t i) const
Definition: TParticle.h:85
Float_t fWeight
Definition: TParticle.h:35
Double_t Z() const
Definition: TVector3.h:218
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...
virtual const char * GetTitle() const
Return particle title.
Definition: TParticle.cxx:298
TParticle & operator=(const TParticle &)
Equal operator.
Definition: TParticle.cxx:134
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
SVector< double, 2 > v
Definition: Dict.h:5
Double_t Theta(const TParticle &p)
Definition: TParticle.h:115
Int_t fDaughter[2]
Definition: TParticle.h:34
Int_t fPdgCode
Definition: TParticle.h:31
Int_t GetLastDaughter() const
Definition: TParticle.h:90
Double_t fVt
Definition: TParticle.h:47
Double_t Vx() const
Definition: TParticle.h:125
void SetCalcMass(Double_t mass)
Definition: TParticle.h:160
Description of the static properties of a particle.
Definition: TParticlePDG.h:19
Double_t PhiY() const
Definition: TParticle.h:101
Int_t Charm() const
Return charm quantum number.
Definition: TParticle.cxx:194
Double_t Energy() const
Double_t ACos(Double_t)
Definition: TMath.h:572
void SetMomentum(const TLorentzVector &p)
Definition: TParticle.h:167
Double_t GetPolarTheta() const
Definition: TParticle.h:106
Double_t PhiZ() const
Definition: TParticle.h:102
Int_t fMother[2]
Definition: TParticle.h:33
void SetLastMother(int code)
Definition: TParticle.h:156
const Bool_t kFALSE
Definition: RtypesCore.h:92
Double_t Px() const
Definition: TParticle.h:131
Double_t T() const
Definition: TParticle.h:128
void SetLastDaughter(int code)
Definition: TParticle.h:159
TParticle()
reference to the particle record in PDG database
Definition: TParticle.cxx:63
Double_t Theta() const
Definition: TParticle.h:150
void GetPolarisation(Double_t &theta, Double_t &phi) const
Definition: TParticle.h:108
Double_t fPx
Definition: TParticle.h:39
Double_t fE
Definition: TParticle.h:42
double Double_t
Definition: RtypesCore.h:55
Double_t Py() const
TParticlePDG * GetPDG(Int_t mode=0) const
Returns a pointer to the TParticlePDG object using the pdgcode.
Definition: TParticle.cxx:273
Double_t Y() const
void SetMomentum(Double_t px, Double_t py, Double_t pz, Double_t e)
Definition: TParticle.h:166
virtual void Paint(Option_t *option="")
Paint a primary track.
Definition: TParticle.cxx:311
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
Double_t Vz() const
Definition: TParticle.h:127
Double_t Px() const
Double_t Vy() const
Definition: TParticle.h:126
Double_t T() const
Double_t Ek() const
Definition: TParticle.h:81
Int_t GetFirstDaughter() const
Definition: TParticle.h:88
Mother of all ROOT objects.
Definition: TObject.h:37
Double_t Pt() const
Definition: TParticle.h:135
void SetFirstMother(int code)
Definition: TParticle.h:154
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 SetProductionVertex(const TLorentzVector &v)
Definition: TParticle.h:169
void SetPolarisation(Double_t theta, Double_t phi)
Definition: TParticle.h:111
void SetFirstDaughter(int code)
Definition: TParticle.h:157
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TParticle.cxx:249
Double_t ThetaZ() const
Definition: TParticle.h:105
void GetPolarisation(TVector3 &v) const
Return particle polarisation.
Definition: TParticle.cxx:284
Double_t X() const
Int_t GetPdgCode() const
Definition: TParticle.h:83
void SetPolarTheta(Double_t theta)
Definition: TParticle.h:109
Double_t fVz
Definition: TParticle.h:46
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
Double_t fVx
Definition: TParticle.h:44
const Bool_t kTRUE
Definition: RtypesCore.h:91
Double_t Pz() const
Definition: TParticle.h:133
Line Attributes class.
Definition: TAttLine.h:18
Double_t fCalcMass
Definition: TParticle.h:37
constexpr Double_t PiOver2()
Definition: TMath.h:48
Double_t PhiX() const
Definition: TParticle.h:100
Double_t fPy
Definition: TParticle.h:40