#include "RConfigure.h"
#include "TROOT.h"
#include "TEnv.h"
#include "THashList.h"
#include "TExMap.h"
#include "TSystem.h"
#include "TDatabasePDG.h"
#include "TDecayChannel.h"
#include "TParticlePDG.h"
#include <stdlib.h>
ClassImp(TDatabasePDG)
TDatabasePDG* TDatabasePDG::fgInstance = 0;
TDatabasePDG::TDatabasePDG(): TNamed("PDGDB","The PDG particle data base")
{
fParticleList = 0;
fPdgMap = 0;
fListOfClasses = 0;
if (fgInstance) {
Warning("TDatabasePDG", "object already instantiated");
} else {
fgInstance = this;
gROOT->GetListOfSpecials()->Add(this);
}
}
TDatabasePDG::~TDatabasePDG()
{
if (fParticleList) {
fParticleList->Delete();
delete fParticleList;
if (fPdgMap) delete fPdgMap;
}
if (fListOfClasses) {
fListOfClasses->Delete();
delete fListOfClasses;
}
gROOT->GetListOfSpecials()->Remove(this);
fgInstance = 0;
}
TDatabasePDG* TDatabasePDG::Instance()
{
return (fgInstance) ? (TDatabasePDG*) fgInstance : new TDatabasePDG();
}
void TDatabasePDG::BuildPdgMap() const
{
fPdgMap = new TExMap(4*TMath::Max(600, fParticleList->GetEntries())/3 + 3);
TIter next(fParticleList);
TParticlePDG *p;
while ((p = (TParticlePDG*)next())) {
fPdgMap->Add((Long_t)p->PdgCode(), (Long_t)p);
}
}
TParticlePDG* TDatabasePDG::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)
{
TParticlePDG* old = GetParticle(PDGcode);
if (old) {
printf(" *** TDatabasePDG::AddParticle: particle with PDGcode=%d already defined\n",PDGcode);
return 0;
}
TParticlePDG* p = new TParticlePDG(name, title, mass, stable, width,
charge, ParticleClass, PDGcode, Anti,
TrackingCode);
fParticleList->Add(p);
if (fPdgMap)
fPdgMap->Add((Long_t)PDGcode, (Long_t)p);
TParticleClassPDG* pclass = GetParticleClass(ParticleClass);
if (!pclass) {
pclass = new TParticleClassPDG(ParticleClass);
fListOfClasses->Add(pclass);
}
pclass->AddParticle(p);
return p;
}
TParticlePDG* TDatabasePDG::AddAntiParticle(const char* Name, Int_t PdgCode)
{
TParticlePDG* old = GetParticle(PdgCode);
if (old) {
printf(" *** TDatabasePDG::AddAntiParticle: can't redefine parameters\n");
return NULL;
}
Int_t pdg_code = abs(PdgCode);
TParticlePDG* p = GetParticle(pdg_code);
if (!p) {
printf(" *** TDatabasePDG::AddAntiParticle: particle with pdg code %d not known\n", pdg_code);
return NULL;
}
TParticlePDG* ap = AddParticle(Name,
Name,
p->Mass(),
1,
p->Width(),
-p->Charge(),
p->ParticleClass(),
PdgCode,
1,
p->TrackingCode());
return ap;
}
TParticlePDG *TDatabasePDG::GetParticle(const char *name) const
{
if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
TParticlePDG *def = (TParticlePDG *)fParticleList->FindObject(name);
return def;
}
TParticlePDG *TDatabasePDG::GetParticle(Int_t PDGcode) const
{
if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
if (fPdgMap == 0) BuildPdgMap();
return (TParticlePDG*) (Long_t)fPdgMap->GetValue((Long_t)PDGcode);
}
void TDatabasePDG::Print(Option_t *option) const
{
if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
TIter next(fParticleList);
TParticlePDG *p;
while ((p = (TParticlePDG *)next())) {
p->Print(option);
}
}
Int_t TDatabasePDG::ConvertGeant3ToPdg(Int_t Geant3number) const
{
/*
see <A href="http://www.slac.stanford.edu/BFROOT/www/Computing/Environment/NewUser/htmlbug/node51.html"> Conversion table</A>
*/
//End_Html
switch(Geant3number) {
case 1 : return 22;
case 25 : return -2112;
case 2 : return -11;
case 26 : return -3122;
case 3 : return 11;
case 27 : return -3222;
case 4 : return 12;
case 28 : return -3212;
case 5 : return -13;
case 29 : return -3112;
case 6 : return 13;
case 30 : return -3322;
case 7 : return 111;
case 31 : return -3312;
case 8 : return 211;
case 32 : return -3334;
case 9 : return -211;
case 33 : return -15;
case 10 : return 130;
case 34 : return 15;
case 11 : return 321;
case 35 : return 411;
case 12 : return -321;
case 36 : return -411;
case 13 : return 2112;
case 37 : return 421;
case 14 : return 2212;
case 38 : return -421;
case 15 : return -2212;
case 39 : return 431;
case 16 : return 310;
case 40 : return -431;
case 17 : return 221;
case 41 : return 4122;
case 18 : return 3122;
case 42 : return 24;
case 19 : return 3222;
case 43 : return -24;
case 20 : return 3212;
case 44 : return 23;
case 21 : return 3112;
case 45 : return 0;
case 22 : return 3322;
case 46 : return 0;
case 23 : return 3312;
case 47 : return 0;
case 24 : return 3334;
case 48 : return 0;
default : return 0;
}
}
Int_t TDatabasePDG::ConvertPdgToGeant3(Int_t pdgNumber) const
{
switch(pdgNumber) {
case 22 : return 1;
case -2112 : return 25;
case -11 : return 2;
case -3122 : return 26;
case 11 : return 3;
case -3222 : return 27;
case 12 : return 4;
case -3212 : return 28;
case -13 : return 5;
case -3112 : return 29;
case 13 : return 6;
case -3322 : return 30;
case 111 : return 7;
case -3312 : return 31;
case 211 : return 8;
case -3334 : return 32;
case -211 : return 9;
case -15 : return 33;
case 130 : return 10;
case 15 : return 34;
case 321 : return 11;
case 411 : return 35;
case -321 : return 12;
case -411 : return 36;
case 2112 : return 13;
case 421 : return 37;
case 2212 : return 14;
case -421 : return 38;
case -2212 : return 15;
case 431 : return 39;
case 310 : return 16;
case -431 : return 40;
case 221 : return 17;
case 4122 : return 41;
case 3122 : return 18;
case 24 : return 42;
case 3222 : return 19;
case -24 : return 43;
case 3212 : return 20;
case 23 : return 44;
case 3112 : return 21;
case 3322 : return 22;
case 3312 : return 23;
case 3334 : return 24;
default : return 0;
}
}
Int_t TDatabasePDG::ConvertIsajetToPdg(Int_t isaNumber) const
{
switch (isaNumber) {
case 1 : return 2;
case -1 : return -2;
case 2 : return 1;
case -2 : return -1;
case 3 : return 3;
case -3 : return -3;
case 4 : return 4;
case -4 : return -4;
case 5 : return 5;
case -5 : return -5;
case 6 : return 6;
case -6 : return -6;
case 9 : return 21;
case 80 : return 24;
case -80 : return -24;
case 90 : return 23;
case 230 : return 311;
case -230 : return -311;
case 330 : return 331;
case 340 : return 0;
case -340 : return 0;
case 440 : return 441;
case 111 : return 113;
case 121 : return 213;
case -121 : return -213;
case 221 : return 223;
case 131 : return 323;
case -131 : return -323;
case 231 : return 313;
case -231 : return -313;
case 331 : return 333;
case -140 : return 421;
case 140 : return -421;
case 141 : return -423;
case -141 : return 423;
case -240 : return -411;
case 240 : return 411;
case 241 : return -413;
case -241 : return 413;
case 341 : return 0;
case -341 : return 0;
case 441 : return 443;
case 250 : return 511;
case -250 : return -511;
case 150 : return 521;
case -150 : return -521;
case 350 : return 531;
case -350 : return -531;
case 351 : return 533;
case -351 : return -533;
case 450 : return 541;
case -450 : return -541;
case 1140 : return 4222;
case -1140 : return -4222;
case 1240 : return 4212;
case -1240 : return -4212;
case 2140 : return 4122;
case -2140 : return -4122;
case 2240 : return 4112;
case -2240 : return -4112;
case 1340 : return 0;
case -1340 : return 0;
case 3140 : return 0;
case -3140 : return 0;
case 2340 : return 0;
case -2340 : return 0;
case 3240 : return 0;
case -3240 : return 0;
case 3340 : return 0;
case -3340 : return 0;
case 1440 : return 0;
case -1440 : return 0;
case 2440 : return 0;
case -2440 : return 0;
case 3440 : return 0;
case -3440 : return 0;
case 1111 : return 2224;
case -1111 : return -2224;
case 1121 : return 2214;
case -1121 : return -2214;
case 1221 : return 2114;
case -1221 : return -2114;
case 2221 : return 1114;
case -2221 : return -1114;
case 1131 : return 3224;
case -1131 : return -3224;
case 1231 : return 3214;
case -1231 : return -3214;
case 2231 : return 3114;
case -2231 : return -3114;
case 1331 : return 3324;
case -1331 : return -3324;
case 2331 : return 3314;
case -2331 : return -3314;
case 3331 : return 3334;
case -3331 : return -3334;
case 1141 : return 0;
case -1141 : return 0;
case 1241 : return 0;
case -1241 : return 0;
case 2241 : return 0;
case -2241 : return 0;
case 1341 : return 0;
case -1341 : return 0;
case 2341 : return 0;
case -2341 : return 0;
case 3341 : return 0;
case -3341 : return 0;
case 1441 : return 0;
case -1441 : return 0;
case 2441 : return 0;
case -2441 : return 0;
case 3441 : return 0;
case -3441 : return 0;
case 4441 : return 0;
case -4441 : return 0;
case 10 : return 22;
case 12 : return 11;
case -12 : return -11;
case 14 : return 13;
case -14 : return -13;
case 16 : return 15;
case -16 : return -15;
case 11 : return 12;
case -11 : return -12;
case 13 : return 14;
case -13 : return -14;
case 15 : return 16;
case -15 : return -16;
case 110 : return 111;
case 120 : return 211;
case -120 : return -211;
case 220 : return 221;
case 130 : return 321;
case -130 : return -321;
case -20 : return 130;
case 20 : return 310;
case 1120 : return 2212;
case -1120 : return -2212;
case 1220 : return 2112;
case -1220 : return -2112;
case 2130 : return 3122;
case -2130 : return -3122;
case 1130 : return 3222;
case -1130 : return -3222;
case 1230 : return 3212;
case -1230 : return -3212;
case 2230 : return 3112;
case -2230 : return -3112;
case 1330 : return 3322;
case -1330 : return -3322;
case 2330 : return 3312;
case -2330 : return -3312;
default : return 0;
}
}
void TDatabasePDG::ReadPDGTable(const char *FileName)
{
if (fParticleList == 0) {
fParticleList = new THashList;
fListOfClasses = new TObjArray;
}
TString default_name;
const char *fn;
if (!FileName[0]) {
#ifdef ROOTETCDIR
default_name.Form("%s/pdg_table.txt", ROOTETCDIR);
#else
default_name.Form("%s/etc/pdg_table.txt", gSystem->Getenv("ROOTSYS"));
#endif
fn = gEnv->GetValue("Root.DatabasePDG", default_name.Data());
} else {
fn = FileName;
}
FILE* file = fopen(fn,"r");
if (file == 0) {
Error("ReadPDGTable","Could not open PDG particle file %s",fn);
return;
}
char c[512];
Int_t class_number, anti, isospin, i3, spin, tracking_code;
Int_t ich, kf, nch, charge;
char name[30], class_name[30];
Double_t mass, width, branching_ratio;
Int_t dau[20];
Int_t idecay, decay_type, flavor, ndau, stable;
Int_t input;
while ( (input=getc(file)) != EOF) {
c[0] = input;
if (c[0] != '#') {
ungetc(c[0],file);
if (fscanf(file,"%i",&ich)) {;}
if (fscanf(file,"%s",name )) {;}
if (fscanf(file,"%i",&kf )) {;}
if (fscanf(file,"%i",&anti )) {;}
if (kf < 0) {
AddAntiParticle(name,kf);
if (fgets(c,200,file)) {;}
} else {
if (fscanf(file,"%i",&class_number)) {;}
if (fscanf(file,"%s",class_name)) {;}
if (fscanf(file,"%i",&charge)) {;}
if (fscanf(file,"%le",&mass)) {;}
if (fscanf(file,"%le",&width)) {;}
if (fscanf(file,"%i",&isospin)) {;}
if (fscanf(file,"%i",&i3)) {;}
if (fscanf(file,"%i",&spin)) {;}
if (fscanf(file,"%i",&flavor)) {;}
if (fscanf(file,"%i",&tracking_code)) {;}
if (fscanf(file,"%i",&nch)) {;}
if (fgets(c,200,file)) {;}
if (width > 1e-10) stable = 0;
else stable = 1;
TParticlePDG* part = AddParticle(name,
name,
mass,
stable,
width,
charge,
class_name,
kf,
anti,
tracking_code);
if (nch) {
ich = 0;
Int_t c_input = 0;
while ( ((c_input=getc(file)) != EOF) && (ich <nch)) {
c[0] = c_input;
if (c[0] != '#') {
ungetc(c[0],file);
if (fscanf(file,"%i",&idecay)) {;}
if (fscanf(file,"%i",&decay_type)) {;}
if (fscanf(file,"%le",&branching_ratio)) {;}
if (fscanf(file,"%i",&ndau)) {;}
for (int idau=0; idau<ndau; idau++) {
if (fscanf(file,"%i",&dau[idau])) {;}
}
if (part) part->AddDecayChannel(decay_type,branching_ratio,ndau,dau);
ich++;
}
if (fgets(c,200,file)) {;}
}
}
}
} else {
if (fgets(c,200,file)) {;}
}
}
TIter it(fParticleList);
Int_t code[20];
TParticlePDG *ap, *p, *daughter;
TDecayChannel *dc;
while ((p = (TParticlePDG*) it.Next())) {
if (p->PdgCode() < 0) {
ap = GetParticle(-p->PdgCode());
if (!ap) continue;
nch = ap->NDecayChannels();
for (ich=0; ich<nch; ich++) {
dc = ap->DecayChannel(ich);
if (!dc) continue;
ndau = dc->NDaughters();
for (int i=0; i<ndau; i++) {
code[i] = dc->DaughterPdgCode(i);
daughter = GetParticle(code[i]);
if (daughter && daughter->AntiParticle()) {
code[i] = -code[i];
}
}
p->AddDecayChannel(dc->MatrixElementCode(),
dc->BranchingRatio(),
dc->NDaughters(),
code);
}
p->SetAntiParticle(ap);
ap->SetAntiParticle(p);
}
}
fclose(file);
return;
}
void TDatabasePDG::Browse(TBrowser* b)
{
if (fListOfClasses ) fListOfClasses->Browse(b);
}
Int_t TDatabasePDG::WritePDGTable(const char *filename)
{
if (fParticleList == 0) {
Error("WritePDGTable","Do not have a valid PDG particle list;"
" consider loading it with ReadPDGTable first.");
return -1;
}
FILE *file = fopen(filename,"w");
if (file == 0) {
Error("WritePDGTable","Could not open PDG particle file %s",filename);
return -1;
}
fprintf(file,"#--------------------------------------------------------------------\n");
fprintf(file,"# i NAME............. KF AP CLASS Q MASS WIDTH 2*I+1 I3 2*S+1 FLVR TrkCod N(dec)\n");
fprintf(file,"#--------------------------------------------------------------------\n");
Int_t nparts=fParticleList->GetEntries();
for(Int_t i=0;i<nparts;++i) {
TParticlePDG *p = dynamic_cast<TParticlePDG*>(fParticleList->At(i));
if(!p) continue;
Int_t ich=i+1;
Int_t kf=p->PdgCode();
fprintf(file,"%5i %-20s %- 6i ", ich, p->GetName(), kf);
Int_t anti=p->AntiParticle() ? 1:0;
if(kf<0) {
for(Int_t j=0;j<nparts;++j) {
TParticlePDG *dummy = dynamic_cast<TParticlePDG*>(fParticleList->At(j));
if(dummy==p->AntiParticle()) {
anti=j+1;
break;
}
}
fprintf(file,"%i 0\n",anti);
continue;
}
fprintf(file,"%i ",anti);
fprintf(file,"%i ",100);
fprintf(file,"%s ",p->ParticleClass());
fprintf(file,"% i ",(Int_t)p->Charge());
fprintf(file,"%.5le ",p->Mass());
fprintf(file,"%.5le ",p->Width());
fprintf(file,"%i ",(Int_t)p->Isospin());
fprintf(file,"%i ",(Int_t)p->I3());
fprintf(file,"%i ",(Int_t)p->Spin());
fprintf(file,"%i ",-1);
fprintf(file,"%i ",p->TrackingCode());
Int_t nch=p->NDecayChannels();
fprintf(file,"%i\n",nch);
if(nch==0) {
continue;
}
fprintf(file,"#----------------------------------------------------------------------\n");
fprintf(file,"# decay type(PY6) BR Nd daughters(codes, then names)\n");
fprintf(file,"#----------------------------------------------------------------------\n");
for(Int_t j=0;j<nch;++j) {
TDecayChannel *dc=p->DecayChannel(j);
if (!dc) continue;
fprintf(file,"%9i ",dc->Number()+1);
fprintf(file,"%3i ",dc->MatrixElementCode());
fprintf(file,"%.5le ",dc->BranchingRatio());
Int_t ndau=dc->NDaughters();
fprintf(file,"%3i ",ndau);
for (int idau=0; idau<ndau; idau++) {
fprintf(file,"%- 6i ",dc->DaughterPdgCode(idau));
}
for (int idau=0; idau<ndau; idau++) {
TParticlePDG *dummy=GetParticle(dc->DaughterPdgCode(idau));
if(dummy)
fprintf(file,"%-10s ",dummy->GetName());
else
fprintf(file,"%-10s ","???");
}
fprintf(file,"\n");
}
}
fclose(file);
return nparts;
}