library: libEG #include "TDatabasePDG.h" |
TDatabasePDG
class description - source file - inheritance tree (.pdf)
public:
TDatabasePDG()
TDatabasePDG(const TDatabasePDG&)
virtual ~TDatabasePDG()
virtual TParticlePDG* AddAntiParticle(const char* Name, Int_t PdgCode)
virtual TParticlePDG* AddParticle(const char* Name, const char* Title, Double_t Mass, Bool_t Stable, Double_t DecayWidth, Double_t Charge, const char* ParticleClass, Int_t PdgCode, Int_t Anti = -1, Int_t TrackingCode = 0)
virtual void Browse(TBrowser* b)
static TClass* Class()
virtual Int_t ConvertGeant3ToPdg(Int_t Geant3Number)
virtual Int_t ConvertIsajetToPdg(Int_t isaNumber)
virtual Int_t ConvertPdgToGeant3(Int_t pdgNumber)
TParticlePDG* GetParticle(Int_t pdgCode) const
TParticlePDG* GetParticle(const char* name) const
TParticleClassPDG* GetParticleClass(const char* name)
static TDatabasePDG* Instance()
virtual TClass* IsA() const
virtual Bool_t IsFolder() const
TDatabasePDG& operator=(const TDatabasePDG&)
const THashList* ParticleList() const
virtual void Print(Option_t* opt = "") const
virtual void ReadPDGTable(const char* filename = "")
virtual void ShowMembers(TMemberInspector& insp, char* parent)
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
virtual Int_t WritePDGTable(const char* filename)
protected:
static TDatabasePDG* fgInstance protect against multiple instances
THashList* fParticleList list of PDG particles
TObjArray* fListOfClasses list of classes (leptons etc.)
Particle database manager class
This manager creates a list of particles which by default is
initialised from with the constants used by PYTHIA6 (plus some
other particles added). See definition and the format of the default
particle list in $ROOTSYS/etc/pdg_table.txt
there are 2 ways of redefining the name of the file containing the
particle properties
1. one can define the name in .rootrc file:
Root.DatabasePDG: $(HOME)/my_pdg_table.txt
2. one can use TDatabasePDG::ReadPDGTable method explicitly:
- TDatabasePDG *pdg = new TDatabasePDG();
- pdg->ReadPDGtable(filename)
See TParticlePDG for the description of a static particle properties.
See TParticle for the description of a dynamic particle particle.
TDatabasePDG(): TNamed("PDGDB","The PDG particle data base")
Create PDG database. Initialization of the DB has to be done via explicit
call to ReadDataBasePDG (also done by GetParticle methods)
~TDatabasePDG()
Cleanup the PDG database.
TDatabasePDG* Instance()
static function
TParticlePDG* AddParticle(const char *name, const char *title,
Double_t mass, Bool_t stable,
Double_t width, Double_t charge,
const char* ParticleClass,
Int_t PDGcode,
Int_t Anti,
Int_t TrackingCode)
Particle definition normal constructor. If the particle is set to be
stable, the decay width parameter does have no meaning and can be set to
any value. The parameters granularity, LowerCutOff and HighCutOff are
used for the construction of the mean free path look up tables. The
granularity will be the number of logwise energy points for which the
mean free path will be calculated.
TParticlePDG* AddAntiParticle(const char* Name, Int_t PdgCode)
assuming particle has already been defined
TParticlePDG* GetParticle(const char *name) const
Get a pointer to the particle object according to the name given
TParticlePDG* GetParticle(Int_t PDGcode) const
Get a pointer to the particle object according to the MC code number
void Print(Option_t *option) const
Print contents of PDG database.
Int_t ConvertGeant3ToPdg(Int_t Geant3number)
Converts Geant3 particle codes to PDG convention. (Geant4 uses
PDG convention already)
Source: BaBar User Guide, Neil I. Geddes,
/*
see Conversion table
*/
with some fixes by PB, marked with (PB) below. Checked against
PDG listings from 2000.
Paul Balm, Nov 19, 2001
Int_t ConvertPdgToGeant3(Int_t pdgNumber)
Converts pdg code to geant3 id
Int_t ConvertIsajetToPdg(Int_t isaNumber)
Converts the ISAJET Particle number into the PDG MC number
void ReadPDGTable(const char *FileName)
read list of particles from a file
if the particle list does not exist, it is created, otherwise
particles are added to the existing list
See $ROOTSYS/etc/pdg_table.txt to see the file format
void Browse(TBrowser* b)
Int_t WritePDGTable(const char * /*filename*/)
write contents of the particle DB into a file
Inline Functions
TParticleClassPDG* GetParticleClass(const char* name)
const THashList* ParticleList() const
Bool_t IsFolder() const
TClass* Class()
TClass* IsA() const
void ShowMembers(TMemberInspector& insp, char* parent)
void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
TDatabasePDG TDatabasePDG(const TDatabasePDG&)
TDatabasePDG& operator=(const TDatabasePDG&)
Author: Pasha Murat 12/02/99
Last update: root/eg:$Name: $:$Id: TDatabasePDG.cxx,v 1.22 2005/09/04 11:42:05 brun Exp $
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.