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