]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TUHKMgen/UHKM/DatabasePDG.h
Moving the classes that belong to the following libraries: STEERBase, ESD, CDB, AOD...
[u/mrichter/AliRoot.git] / TUHKMgen / UHKM / DatabasePDG.h
1 // DatabasePDG stores and handles PDG information
2 // The PDG particle definitions and decay channels are read
3 // in the begining from ASCII files
4
5 /*
6   Copyright   : The FASTMC and SPHMC Collaboration
7   Author      : Ionut Cristian Arsene 
8   Affiliation : Oslo University, Norway & Institute for Space Sciences, Bucharest, Romania
9   e-mail      : i.c.arsene@fys.uio.no
10   Date        : 2007/05/30
11
12   This class is using the particle and decay lists provided by the 
13   THERMINATOR (Computer Physics Communications 174 669 (2006)) and
14   SHARE (Computer Physics Communications 167 229 (2005)) collaborations.
15 */
16
17 #ifndef DATABASEPDG_H
18 #define DATABASEPDG_H
19
20 #include <Rtypes.h>
21 #ifndef PARTICLE_PDG
22 #include "ParticlePDG.h"
23 #endif
24
25 const Int_t kMaxParticles = 1000;
26
27 class DatabasePDG {
28  public:
29   DatabasePDG();
30   ~DatabasePDG();
31
32   // Load the particle PDG information from the particle and decay files
33   Bool_t LoadData();                        
34   
35   // Set particle and decay filenames
36   void SetParticleFilename(Char_t *filename);
37   void SetDecayFilename(Char_t *filename);
38   // Set criteria for using particles. Those particle which do not match
39   // these criteria will be flagged with FALSE in the fStatus array.
40   void SetUseCharmParticles(Bool_t flag);
41   void SetMinimumWidth(Double_t value);
42   void SetMaximumWidth(Double_t value);
43   void SetMinimumMass(Double_t value);
44   void SetMaximumMass(Double_t value);
45   void SetWidthRange(Double_t min, Double_t max);
46   void SetMassRange(Double_t min, Double_t max);
47   
48   // Read a list of pdg codes from a specified file. The corresponding particles
49   // will be flagged as good particles. If the exclusive flag is TRUE than
50   // only this criteria will be used in selecting particles and, in consequence,
51   // all the other particles will be flagged as NOT good. If the exclusive flag
52   // is FALSE than we will take into account all the previous applied criterias
53   // and we will flag as good only particles in this list which match also the mass, width and
54   // charmness criteria.
55   // Note: In order for the exclusive=FALSE to be effective, this function must be called after
56   // calling all the width, mass and charmness criteria functions.
57   void UseThisListOfParticles(Char_t *filename, Bool_t exclusive = kTRUE);
58   
59   const Char_t* GetParticleFilename() const {return fParticleFilename;}
60   const Char_t* GetDecayFilename() const {return fDecayFilename;}
61   Int_t GetNParticles(Bool_t all = kFALSE) const;      // true - no. of all particles; false - no. of good status particles
62   ParticlePDG* GetPDGParticleByIndex(Int_t index) const;
63   Bool_t GetPDGParticleStatusByIndex(Int_t index) const;
64   ParticlePDG* GetPDGParticle(Int_t pdg) const;
65   Bool_t GetPDGParticleStatus(Int_t pdg) const;
66   ParticlePDG* GetPDGParticle(Char_t *name) const;
67   Bool_t GetPDGParticleStatus(Char_t *name) const; 
68   Bool_t GetUseCharmParticles() const {return fUseCharmParticles;};
69   Double_t GetMinimumWidth() const {return fMinimumWidth;};
70   Double_t GetMaximumWidth() const {return fMaximumWidth;};
71   Double_t GetMinimumMass() const {return fMinimumMass;};
72   Double_t GetMaximumMass() const {return fMaximumMass;};
73   void DumpData(Bool_t dumpAll = kFALSE) const; // print the PDG information in the console
74   Int_t CheckImpossibleDecays(Bool_t dump = kFALSE) const;   // print all impossible decays included in the database
75   Bool_t IsChannelAllowed(DecayChannel *channel, Double_t motherMass) const;
76   Int_t GetNAllowedChannels(ParticlePDG *particle, Double_t motherMass) const;
77   void SetStable(Int_t pdg, Bool_t value) {GetPDGParticle(pdg)->SetStable(value);}
78   Bool_t GetStableStatus(Int_t pdg) const {return GetPDGParticle(pdg)->GetStableStatus();}
79
80  private:
81   DatabasePDG(const DatabasePDG&);
82   DatabasePDG& operator=(const DatabasePDG&);
83
84   Int_t fNParticles;                        // no. of particles in database
85   ParticlePDG *fParticles[kMaxParticles];   // array of particle pointers
86   Bool_t fStatus[kMaxParticles];            // status of each particle
87   Char_t fParticleFilename[256];            // particle list filename
88   Char_t fDecayFilename[256];               // decay channels filename
89   Bool_t fUseCharmParticles;                // flag for using (or not) charm particles
90   Double_t fMinimumWidth;                   // minimum allowed width for resonances
91   Double_t fMaximumWidth;                   // maximum allowed width for resonances
92   Double_t fMinimumMass;                    // minimum allowed mass for resonances
93   Double_t fMaximumMass;                    // maximum allowed mass for resonances
94
95   Bool_t LoadParticles();
96   Bool_t LoadDecays();
97   void SortParticles();                     // put the good status particles at the beggining of the list
98 };
99
100 #endif