TDatabasePDG.cxx

Go to the documentation of this file.
00001 // @(#)root/eg:$Id: TDatabasePDG.cxx 36178 2010-10-08 09:43:51Z brun $
00002 // Author: Pasha Murat   12/02/99
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "RConfigure.h"
00013 #include "TROOT.h"
00014 #include "TEnv.h"
00015 #include "THashList.h"
00016 #include "TExMap.h"
00017 #include "TSystem.h"
00018 #include "TDatabasePDG.h"
00019 #include "TDecayChannel.h"
00020 #include "TParticlePDG.h"
00021 #include <stdlib.h>
00022 
00023 
00024 ////////////////////////////////////////////////////////////////////////
00025 //
00026 //  Particle database manager class
00027 //
00028 //  This manager creates a list of particles which by default is
00029 //  initialised from with the constants used by PYTHIA6 (plus some
00030 //  other particles added). See definition and the format of the default
00031 //  particle list in $ROOTSYS/etc/pdg_table.txt
00032 //
00033 //  there are 2 ways of redefining the name of the file containing the
00034 //  particle properties
00035 //
00036 //  1. one can define the name in .rootrc file:
00037 //
00038 //  Root.DatabasePDG: $(HOME)/my_pdg_table.txt
00039 //
00040 //  2. one can use TDatabasePDG::ReadPDGTable method explicitly:
00041 //
00042 //     - TDatabasePDG *pdg = new TDatabasePDG();
00043 //     - pdg->ReadPDGtable(filename)
00044 //
00045 //  See TParticlePDG for the description of a static particle properties.
00046 //  See TParticle    for the description of a dynamic particle particle.
00047 //
00048 ////////////////////////////////////////////////////////////////////////
00049 
00050 ClassImp(TDatabasePDG)
00051 
00052 TDatabasePDG*  TDatabasePDG::fgInstance = 0;
00053 
00054 //______________________________________________________________________________
00055 TDatabasePDG::TDatabasePDG(): TNamed("PDGDB","The PDG particle data base")
00056 {
00057   // Create PDG database. Initialization of the DB has to be done via explicit
00058   // call to ReadDataBasePDG (also done by GetParticle methods)
00059 
00060    fParticleList  = 0;
00061    fPdgMap        = 0;
00062    fListOfClasses = 0;
00063       if (fgInstance) {
00064       Warning("TDatabasePDG", "object already instantiated");
00065    } else {
00066       fgInstance = this;
00067       gROOT->GetListOfSpecials()->Add(this);
00068    }
00069 }
00070 
00071 //______________________________________________________________________________
00072 TDatabasePDG::~TDatabasePDG()
00073 {
00074    // Cleanup the PDG database.
00075 
00076    if (fParticleList) {
00077       fParticleList->Delete();
00078       delete fParticleList;
00079       delete fPdgMap;
00080    }
00081                                 // classes do not own particles...
00082    if (fListOfClasses) delete fListOfClasses;
00083    gROOT->GetListOfSpecials()->Remove(this);
00084    fgInstance = 0;
00085 }
00086 
00087 //______________________________________________________________________________
00088 TDatabasePDG*  TDatabasePDG::Instance()
00089 {
00090    //static function
00091    return (fgInstance) ? (TDatabasePDG*) fgInstance : new TDatabasePDG();
00092 }
00093 
00094 //______________________________________________________________________________
00095 void TDatabasePDG::BuildPdgMap() const
00096 {
00097    // Build fPdgMap mapping pdg-code to particle.
00098    //
00099    // Initial size is set so as to be able to hold at least 600
00100    // particles: 521 in default table, ALICE adds 54 more.
00101    // To be revisited after LHC discovers SUSY.
00102 
00103    fPdgMap = new TExMap(4*TMath::Max(600, fParticleList->GetEntries())/3 + 3);
00104    TIter next(fParticleList);
00105    TParticlePDG *p;
00106    while ((p = (TParticlePDG*)next())) {
00107       fPdgMap->Add((Long_t)p->PdgCode(), (Long_t)p);
00108    }
00109 }
00110 
00111 //______________________________________________________________________________
00112 TParticlePDG* TDatabasePDG::AddParticle(const char *name, const char *title,
00113                                         Double_t mass, Bool_t stable,
00114                                         Double_t width, Double_t charge,
00115                                         const char* ParticleClass,
00116                                         Int_t PDGcode,
00117                                         Int_t Anti,
00118                                         Int_t TrackingCode)
00119 {
00120   //
00121   //  Particle definition normal constructor. If the particle is set to be
00122   //  stable, the decay width parameter does have no meaning and can be set to
00123   //  any value. The parameters granularity, LowerCutOff and HighCutOff are
00124   //  used for the construction of the mean free path look up tables. The
00125   //  granularity will be the number of logwise energy points for which the
00126   //  mean free path will be calculated.
00127   //
00128 
00129    TParticlePDG* old = GetParticle(PDGcode);
00130 
00131    if (old) {
00132       printf(" *** TDatabasePDG::AddParticle: particle with PDGcode=%d already defined\n",PDGcode);
00133       return 0;
00134    }
00135 
00136    TParticlePDG* p = new TParticlePDG(name, title, mass, stable, width,
00137                                      charge, ParticleClass, PDGcode, Anti,
00138                                      TrackingCode);
00139    fParticleList->Add(p);
00140    if (fPdgMap)
00141       fPdgMap->Add((Long_t)PDGcode, (Long_t)p);
00142 
00143    TParticleClassPDG* pclass = GetParticleClass(ParticleClass);
00144 
00145    if (!pclass) {
00146       pclass = new TParticleClassPDG(ParticleClass);
00147       fListOfClasses->Add(pclass);
00148    }
00149 
00150    pclass->AddParticle(p);
00151 
00152    return p;
00153 }
00154 
00155 //______________________________________________________________________________
00156 TParticlePDG* TDatabasePDG::AddAntiParticle(const char* Name, Int_t PdgCode)
00157 {
00158    // assuming particle has already been defined
00159 
00160    TParticlePDG* old = GetParticle(PdgCode);
00161 
00162    if (old) {
00163       printf(" *** TDatabasePDG::AddAntiParticle: can't redefine parameters\n");
00164       return NULL;
00165    }
00166 
00167    Int_t pdg_code  = abs(PdgCode);
00168    TParticlePDG* p = GetParticle(pdg_code);
00169 
00170    TParticlePDG* ap = AddParticle(Name,
00171                                   Name,
00172                                   p->Mass(),
00173                                   1,
00174                                   p->Width(),
00175                                   -p->Charge(),
00176                                   p->ParticleClass(),
00177                                   PdgCode,
00178                                   1,
00179                                   p->TrackingCode());
00180    return ap;
00181 }
00182 
00183 
00184 //______________________________________________________________________________
00185 TParticlePDG *TDatabasePDG::GetParticle(const char *name) const
00186 {
00187    //
00188    //  Get a pointer to the particle object according to the name given
00189    //
00190 
00191    if (fParticleList == 0)  ((TDatabasePDG*)this)->ReadPDGTable();
00192 
00193    TParticlePDG *def = (TParticlePDG *)fParticleList->FindObject(name);
00194 //     if (!def) {
00195 //        Error("GetParticle","No match for %s exists!",name);
00196 //     }
00197    return def;
00198 }
00199 
00200 //______________________________________________________________________________
00201 TParticlePDG *TDatabasePDG::GetParticle(Int_t PDGcode) const
00202 {
00203    //
00204    //  Get a pointer to the particle object according to the MC code number
00205    //
00206 
00207    if (fParticleList == 0)  ((TDatabasePDG*)this)->ReadPDGTable();
00208    if (fPdgMap       == 0)  BuildPdgMap();
00209 
00210    return (TParticlePDG*) (Long_t)fPdgMap->GetValue((Long_t)PDGcode);
00211 }
00212 
00213 //______________________________________________________________________________
00214 void TDatabasePDG::Print(Option_t *option) const
00215 {
00216    // Print contents of PDG database.
00217 
00218    if (fParticleList == 0)  ((TDatabasePDG*)this)->ReadPDGTable();
00219 
00220    TIter next(fParticleList);
00221    TParticlePDG *p;
00222    while ((p = (TParticlePDG *)next())) {
00223       p->Print(option);
00224    }
00225 }
00226 
00227 //______________________________________________________________________________
00228 Int_t TDatabasePDG::ConvertGeant3ToPdg(Int_t Geant3number) {
00229   // Converts Geant3 particle codes to PDG convention. (Geant4 uses
00230   // PDG convention already)
00231   // Source: BaBar User Guide, Neil I. Geddes,
00232   //
00233   //Begin_Html
00234   /*
00235    see <A href="http://www.slac.stanford.edu/BFROOT/www/Computing/Environment/NewUser/htmlbug/node51.html"> Conversion table</A>
00236   */
00237   //End_Html
00238   // with some fixes by PB, marked with (PB) below. Checked against
00239   // PDG listings from 2000.
00240   //
00241   // Paul Balm, Nov 19, 2001
00242 
00243    switch(Geant3number) {
00244 
00245       case 1   : return 22;       // photon
00246       case 25  : return -2112;    // anti-neutron
00247       case 2   : return -11;      // e+
00248       case 26  : return -3122;    // anti-Lambda
00249       case 3   : return 11;       // e-
00250       case 27  : return -3222;    // Sigma-
00251       case 4   : return 12;       // e-neutrino (NB: flavour undefined by Geant)
00252       case 28  : return -3212;    // Sigma0
00253       case 5   : return -13;      // mu+
00254       case 29  : return -3112;    // Sigma+ (PB)*/
00255       case 6   : return 13;       // mu-
00256       case 30  : return -3322;    // Xi0
00257       case 7   : return 111;      // pi0
00258       case 31  : return -3312;    // Xi+
00259       case 8   : return 211;      // pi+
00260       case 32  : return -3334;    // Omega+ (PB)
00261       case 9   : return -211;     // pi-
00262       case 33  : return -15;      // tau+
00263       case 10  : return 130;      // K long
00264       case 34  : return 15;       // tau-
00265       case 11  : return 321;      // K+
00266       case 35  : return 411;      // D+
00267       case 12  : return -321;     // K-
00268       case 36  : return -411;     // D-
00269       case 13  : return 2112;     // n
00270       case 37  : return 421;      // D0
00271       case 14  : return 2212;     // p
00272       case 38  : return -421;     // D0
00273       case 15  : return -2212;    // anti-proton
00274       case 39  : return 431;      // Ds+
00275       case 16  : return 310;      // K short
00276       case 40  : return -431;     // anti Ds-
00277       case 17  : return 221;      // eta
00278       case 41  : return 4122;     // Lamba_c+
00279       case 18  : return 3122;     // Lambda
00280       case 42  : return 24;       // W+
00281       case 19  : return 3222;     // Sigma+
00282       case 43  : return -24;      // W-
00283       case 20  : return 3212;     // Sigma0
00284       case 44  : return 23;       // Z
00285       case 21  : return 3112;     // Sigma-
00286       case 45  : return 0;        // deuteron
00287       case 22  : return 3322;     // Xi0
00288       case 46  : return 0;        // triton
00289       case 23  : return 3312;     // Xi-
00290       case 47  : return 0;        // alpha
00291       case 24  : return 3334;     // Omega- (PB)
00292       case 48  : return 0;        // G nu ? PDG ID 0 is undefined
00293 
00294       default  : return 0;
00295 
00296    }
00297 }
00298 
00299 //______________________________________________________________________________
00300 Int_t TDatabasePDG::ConvertPdgToGeant3(Int_t pdgNumber) {
00301    // Converts pdg code to geant3 id
00302 
00303    switch(pdgNumber) {
00304 
00305       case   22     : return  1;    // photon
00306       case   -2112  : return  25;   // anti-neutron
00307       case   -11    : return  2;    // e+
00308       case   -3122  : return  26;   // anti-Lambda
00309       case   11     : return  3;    // e-
00310       case   -3222  : return  27;   // Sigma-
00311       case   12     : return  4;    // e-neutrino (NB: flavour undefined by Geant)
00312       case   -3212  : return  28;   // Sigma0
00313       case   -13    : return  5;    // mu+
00314       case   -3112  : return  29;   // Sigma+ (PB)*/
00315       case   13     : return  6;    // mu-
00316       case   -3322  : return  30;   // Xi0
00317       case   111    : return  7;    // pi0
00318       case   -3312  : return  31;   // Xi+
00319       case   211    : return  8;    // pi+
00320       case   -3334  : return  32;   // Omega+ (PB)
00321       case   -211   : return  9;    // pi-
00322       case   -15    : return  33;   // tau+
00323       case   130    : return  10;   // K long
00324       case   15     : return  34;   // tau-
00325       case   321    : return  11;   // K+
00326       case   411    : return  35;   // D+
00327       case   -321   : return  12;   // K-
00328       case   -411   : return  36;   // D-
00329       case   2112   : return  13;   // n
00330       case   421    : return  37;   // D0
00331       case   2212   : return  14;   // p
00332       case   -421   : return  38;   // D0
00333       case   -2212  : return  15;   // anti-proton
00334       case   431    : return  39;   // Ds+
00335       case   310    : return  16;   // K short
00336       case   -431   : return  40;   // anti Ds-
00337       case   221    : return  17;   // eta
00338       case   4122   : return  41;   // Lamba_c+
00339       case   3122   : return  18;   // Lambda
00340       case   24     : return  42;   // W+
00341       case   3222   : return  19;   // Sigma+
00342       case   -24    : return  43;   // W-
00343       case   3212   : return  20;   // Sigma0
00344       case   23     : return  44;   // Z
00345       case   3112   : return  21;   // Sigma-
00346       case   3322   : return  22;   // Xi0
00347       case   3312   : return  23;   // Xi-
00348       case   3334   : return  24;   // Omega- (PB)
00349 
00350       default  : return 0;
00351 
00352    }
00353 }
00354 
00355 //______________________________________________________________________________
00356 Int_t TDatabasePDG::ConvertIsajetToPdg(Int_t isaNumber)
00357 {
00358 //
00359 //  Converts the ISAJET Particle number into the PDG MC number
00360 //
00361    switch (isaNumber) {
00362       case     1 : return     2; //     UP        .30000E+00       .67
00363       case    -1 : return    -2; //     UB        .30000E+00      -.67
00364       case     2 : return     1; //     DN        .30000E+00      -.33
00365       case    -2 : return    -1; //     DB        .30000E+00       .33
00366       case     3 : return     3; //     ST        .50000E+00      -.33
00367       case    -3 : return    -3; //     SB        .50000E+00       .33
00368       case     4 : return     4; //     CH        .16000E+01       .67
00369       case    -4 : return    -4; //     CB        .16000E+01      -.67
00370       case     5 : return     5; //     BT        .49000E+01      -.33
00371       case    -5 : return    -5; //     BB        .49000E+01       .33
00372       case     6 : return     6; //     TP        .17500E+03       .67
00373       case    -6 : return    -6; //     TB        .17500E+03      -.67
00374       case     9 : return    21; //     GL       0.               0.00
00375       case    80 : return    24; //     W+        SIN2W=.23       1.00
00376       case   -80 : return   -24; //     W-        SIN2W=.23      -1.00
00377       case    90 : return    23; //     Z0        SIN2W=.23       0.00
00378       case   230 : return   311; //     K0        .49767E+00      0.00
00379       case  -230 : return  -311; //     AK0       .49767E+00      0.00
00380       case   330 : return   331; //     ETAP      .95760E+00      0.00
00381       case   340 : return     0; //     F-        .20300E+01     -1.00
00382       case  -340 : return     0; //     F+        .20300E+01      1.00
00383       case   440 : return   441; //     ETAC      .29760E+01      0.00
00384       case   111 : return   113; //     RHO0      .77000E+00      0.00
00385       case   121 : return   213; //     RHO+      .77000E+00      1.00
00386       case  -121 : return  -213; //     RHO-      .77000E+00     -1.00
00387       case   221 : return   223; //     OMEG      .78260E+00      0.00
00388       case   131 : return   323; //     K*+       .88810E+00      1.00
00389       case  -131 : return  -323; //     K*-       .88810E+00     -1.00
00390       case   231 : return   313; //     K*0       .89220E+00      0.00
00391       case  -231 : return  -313; //     AK*0      .89220E+00      0.00
00392       case   331 : return   333; //     PHI       .10196E+01      0.00
00393       case  -140 : return   421; //     D0
00394       case   140 : return  -421; //     D0 bar
00395       case   141 : return  -423; //     AD*0      .20060E+01      0.00
00396       case  -141 : return   423; //     D*0       .20060E+01      0.00
00397       case  -240 : return  -411; //     D+
00398       case   240 : return   411; //     D-
00399       case   241 : return  -413; //     D*-       .20086E+01     -1.00
00400       case  -241 : return   413; //     D*+       .20086E+01      1.00
00401       case   341 : return     0; //     F*-       .21400E+01     -1.00
00402       case  -341 : return     0; //     F*+       .21400E+01      1.00
00403       case   441 : return   443; //     JPSI      .30970E+01      0.00
00404 
00405                                         // B-mesons, Bc still missing
00406       case   250 : return   511; // B0
00407       case  -250 : return  -511; // B0 bar
00408       case   150 : return   521; // B+
00409       case  -150 : return  -521; // B-
00410       case   350 : return   531; // Bs  0
00411       case  -350 : return  -531; // Bs  bar
00412       case   351 : return   533; // Bs* 0
00413       case  -351 : return  -533; // Bs* bar
00414       case   450 : return   541; // Bc  +
00415       case  -450 : return  -541; // Bc  bar
00416 
00417       case  1140 : return  4222; //     SC++      .24300E+01      2.00
00418       case -1140 : return -4222; //     ASC--     .24300E+01     -2.00
00419       case  1240 : return  4212; //     SC+       .24300E+01      1.00
00420       case -1240 : return -4212; //     ASC-      .24300E+01     -1.00
00421       case  2140 : return  4122; //     LC+       .22600E+01      1.00
00422       case -2140 : return -4122; //     ALC-      .22600E+01     -1.00
00423       case  2240 : return  4112; //     SC0       .24300E+01      0.00
00424       case -2240 : return -4112; //     ASC0      .24300E+01      0.00
00425       case  1340 : return     0; //     USC.      .25000E+01      1.00
00426       case -1340 : return     0; //     AUSC.     .25000E+01     -1.00
00427       case  3140 : return     0; //     SUC.      .24000E+01      1.00
00428       case -3140 : return     0; //     ASUC.     .24000E+01     -1.00
00429       case  2340 : return     0; //     DSC.      .25000E+01      0.00
00430       case -2340 : return     0; //     ADSC.     .25000E+01      0.00
00431       case  3240 : return     0; //     SDC.      .24000E+01      0.00
00432       case -3240 : return     0; //     ASDC.     .24000E+01      0.00
00433       case  3340 : return     0; //     SSC.      .26000E+01      0.00
00434       case -3340 : return     0; //     ASSC.     .26000E+01      0.00
00435       case  1440 : return     0; //     UCC.      .35500E+01      2.00
00436       case -1440 : return     0; //     AUCC.     .35500E+01     -2.00
00437       case  2440 : return     0; //     DCC.      .35500E+01      1.00
00438       case -2440 : return     0; //     ADCC.     .35500E+01     -1.00
00439       case  3440 : return     0; //     SCC.      .37000E+01      1.00
00440       case -3440 : return     0; //     ASCC.     .37000E+01     -1.00
00441       case  1111 : return  2224; //     DL++      .12320E+01      2.00
00442       case -1111 : return -2224; //     ADL--     .12320E+01     -2.00
00443       case  1121 : return  2214; //     DL+       .12320E+01      1.00
00444       case -1121 : return -2214; //     ADL-      .12320E+01     -1.00
00445       case  1221 : return  2114; //     DL0       .12320E+01      0.00
00446       case -1221 : return -2114; //     ADL0      .12320E+01      0.00
00447       case  2221 : return   1114; //     DL-       .12320E+01     -1.00
00448       case -2221 : return -1114; //     ADL+      .12320E+01      1.00
00449       case  1131 : return  3224; //     S*+       .13823E+01      1.00
00450       case -1131 : return -3224; //     AS*-      .13823E+01     -1.00
00451       case  1231 : return  3214; //     S*0       .13820E+01      0.00
00452       case -1231 : return -3214; //     AS*0      .13820E+01      0.00
00453       case  2231 : return  3114; //     S*-       .13875E+01     -1.00
00454       case -2231 : return -3114; //     AS*+      .13875E+01      1.00
00455       case  1331 : return  3324; //     XI*0      .15318E+01      0.00
00456       case -1331 : return -3324; //     AXI*0     .15318E+01      0.00
00457       case  2331 : return  3314; //     XI*-      .15350E+01     -1.00
00458       case -2331 : return -3314; //     AXI*+     .15350E+01      1.00
00459       case  3331 : return  3334; //     OM-       .16722E+01     -1.00
00460       case -3331 : return -3334; //     AOM+      .16722E+01      1.00
00461       case  1141 : return     0; //     UUC*      .26300E+01      2.00
00462       case -1141 : return     0; //     AUUC*     .26300E+01     -2.00
00463       case  1241 : return     0; //     UDC*      .26300E+01      1.00
00464       case -1241 : return     0; //     AUDC*     .26300E+01     -1.00
00465       case  2241 : return     0; //     DDC*      .26300E+01      0.00
00466       case -2241 : return     0; //     ADDC*     .26300E+01      0.00
00467       case  1341 : return     0; //     USC*      .27000E+01      1.00
00468       case -1341 : return     0; //     AUSC*     .27000E+01     -1.00
00469       case  2341 : return     0; //     DSC*      .27000E+01      0.00
00470       case -2341 : return     0; //     ADSC*     .27000E+01      0.00
00471       case  3341 : return     0; //     SSC*      .28000E+01      0.00
00472       case -3341 : return     0; //     ASSC*     .28000E+01      0.00
00473       case  1441 : return     0; //     UCC*      .37500E+01      2.00
00474       case -1441 : return     0; //     AUCC*     .37500E+01     -2.00
00475       case  2441 : return     0; //     DCC*      .37500E+01      1.00
00476       case -2441 : return     0; //     ADCC*     .37500E+01     -1.00
00477       case  3441 : return     0; //     SCC*      .39000E+01      1.00
00478       case -3441 : return     0; //     ASCC*     .39000E+01     -1.00
00479       case  4441 : return     0; //     CCC*      .48000E+01      2.00
00480       case -4441 : return     0; //     ACCC*     .48000E+01     -2.00
00481       case    10 : return    22; // Photon
00482       case    12 : return    11; // Electron
00483       case   -12 : return   -11; // Positron
00484       case    14 : return    13; // Muon-
00485       case   -14 : return   -13; // Muon+
00486       case    16 : return    15; // Tau-
00487       case   -16 : return   -15; // Tau+
00488       case    11 : return    12; // Neutrino e
00489       case   -11 : return   -12; // Anti Neutrino e
00490       case    13 : return    14; // Neutrino Muon
00491       case   -13 : return   -14; // Anti Neutrino Muon
00492       case    15 : return    16; // Neutrino Tau
00493       case   -15 : return   -16; // Anti Neutrino Tau
00494       case   110 : return   111; // Pion0
00495       case   120 : return   211; // Pion+
00496       case  -120 : return  -211; // Pion-
00497       case   220 : return   221; // Eta
00498       case   130 : return   321; // Kaon+
00499       case  -130 : return  -321; // Kaon-
00500       case   -20 : return   130; // Kaon Long
00501       case    20 : return   310; // Kaon Short
00502 
00503                                         // baryons
00504       case  1120 : return  2212; // Proton
00505       case -1120 : return -2212; // Anti Proton
00506       case  1220 : return  2112; // Neutron
00507       case -1220 : return -2112; // Anti Neutron
00508       case  2130 : return  3122; // Lambda
00509       case -2130 : return -3122; // Lambda bar
00510       case  1130 : return  3222; // Sigma+
00511       case -1130 : return -3222; // Sigma bar -
00512       case  1230 : return  3212; // Sigma0
00513       case -1230 : return -3212; // Sigma bar 0
00514       case  2230 : return  3112; // Sigma-
00515       case -2230 : return -3112; // Sigma bar +
00516       case  1330 : return  3322; // Xi0
00517       case -1330 : return -3322; // Xi bar 0
00518       case  2330 : return  3312; // Xi-
00519       case -2330 : return -3312; // Xi bar +
00520       default :    return 0;      // isajet or pdg number does not exist
00521    }
00522 }
00523 
00524 //______________________________________________________________________________
00525 void TDatabasePDG::ReadPDGTable(const char *FileName)
00526 {
00527    // read list of particles from a file
00528    // if the particle list does not exist, it is created, otherwise
00529    // particles are added to the existing list
00530    // See $ROOTSYS/etc/pdg_table.txt to see the file format
00531 
00532    if (fParticleList == 0) {
00533       fParticleList  = new THashList;
00534       fListOfClasses = new TObjArray;
00535    }
00536 
00537    TString default_name;
00538    const char *fn;
00539 
00540    if (strlen(FileName) == 0) {
00541 #ifdef ROOTETCDIR
00542       default_name.Form("%s/pdg_table.txt", ROOTETCDIR);
00543 #else
00544       default_name.Form("%s/etc/pdg_table.txt", gSystem->Getenv("ROOTSYS"));
00545 #endif
00546       fn = gEnv->GetValue("Root.DatabasePDG", default_name.Data());
00547    } else {
00548       fn = FileName;
00549    }
00550 
00551    FILE* file = fopen(fn,"r");
00552    if (file == 0) {
00553       Error("ReadPDGTable","Could not open PDG particle file %s",fn);
00554       return;
00555    }
00556 
00557    char      c[512];
00558    Int_t     class_number, anti, isospin, i3, spin, tracking_code;
00559    Int_t     ich, kf, nch, charge;
00560    char      name[30], class_name[30];
00561    Double_t  mass, width, branching_ratio;
00562    Int_t     dau[20];
00563 
00564    Int_t     idecay, decay_type, flavor, ndau, stable;
00565 
00566    Int_t input;
00567    while ( (input=getc(file)) != EOF) {
00568       c[0] = input;
00569       if (c[0] != '#') {
00570          ungetc(c[0],file);
00571          // read channel number
00572          // coverity [secure_coding : FALSE]
00573          if (fscanf(file,"%i",&ich)) {;}
00574          // coverity [secure_coding : FALSE]
00575          if (fscanf(file,"%s",name  )) {;}
00576          // coverity [secure_coding : FALSE]
00577          if (fscanf(file,"%i",&kf   )) {;}
00578          // coverity [secure_coding : FALSE]
00579          if (fscanf(file,"%i",&anti )) {;}
00580 
00581          if (kf < 0) {
00582             AddAntiParticle(name,kf);
00583             // nothing more on this line
00584             if (fgets(c,200,file)) {;}
00585          } else {
00586             // coverity [secure_coding : FALSE]
00587             if (fscanf(file,"%i",&class_number)) {;}
00588             // coverity [secure_coding : FALSE]
00589             if (fscanf(file,"%s",class_name)) {;}
00590             // coverity [secure_coding : FALSE]
00591             if (fscanf(file,"%i",&charge)) {;}
00592             // coverity [secure_coding : FALSE]
00593             if (fscanf(file,"%le",&mass)) {;}
00594             // coverity [secure_coding : FALSE]
00595             if (fscanf(file,"%le",&width)) {;}
00596             // coverity [secure_coding : FALSE]
00597             if (fscanf(file,"%i",&isospin)) {;}
00598             // coverity [secure_coding : FALSE]
00599             if (fscanf(file,"%i",&i3)) {;}
00600             // coverity [secure_coding : FALSE]
00601             if (fscanf(file,"%i",&spin)) {;}
00602             // coverity [secure_coding : FALSE]
00603             if (fscanf(file,"%i",&flavor)) {;}
00604             // coverity [secure_coding : FALSE]
00605             if (fscanf(file,"%i",&tracking_code)) {;}
00606             // coverity [secure_coding : FALSE]
00607             if (fscanf(file,"%i",&nch)) {;}
00608             // nothing more on this line
00609             if (fgets(c,200,file)) {;}
00610             if (width > 1e-10) stable = 0;
00611             else               stable = 1;
00612 
00613             // create particle
00614 
00615             TParticlePDG* part = AddParticle(name,
00616                                              name,
00617                                              mass,
00618                                              stable,
00619                                              width,
00620                                              charge,
00621                                              class_name,
00622                                              kf,
00623                                              anti,
00624                                              tracking_code);
00625 
00626             if (nch) {
00627                // read in decay channels
00628                ich = 0;
00629                Int_t c_input = 0;
00630                while ( ((c_input=getc(file)) != EOF) && (ich <nch)) {
00631                   c[0] = c_input;
00632                   if (c[0] != '#') {
00633                      ungetc(c[0],file);
00634 
00635                      // coverity [secure_coding : FALSE]
00636                      if (fscanf(file,"%i",&idecay)) {;}
00637                      // coverity [secure_coding : FALSE]
00638                      if (fscanf(file,"%i",&decay_type)) {;}
00639                      // coverity [secure_coding : FALSE]
00640                      if (fscanf(file,"%le",&branching_ratio)) {;}
00641                      // coverity [secure_coding : FALSE]
00642                      if (fscanf(file,"%i",&ndau)) {;}
00643                      for (int idau=0; idau<ndau; idau++) {
00644                         // coverity [secure_coding : FALSE]
00645                         if (fscanf(file,"%i",&dau[idau])) {;}
00646                      }
00647                      // add decay channel
00648 
00649                      if (part) part->AddDecayChannel(decay_type,branching_ratio,ndau,dau);
00650                      ich++;
00651                   }
00652                   // skip end of line
00653                   if (fgets(c,200,file)) {;}
00654                }
00655             }
00656          }
00657       } else {
00658          // skip end of line
00659          if (fgets(c,200,file)) {;}
00660       }
00661    }
00662    // in the end loop over the antiparticles and
00663    // define their decay lists
00664    TIter it(fParticleList);
00665 
00666    Int_t code[20];
00667    TParticlePDG  *ap, *p, *daughter;
00668    TDecayChannel *dc;
00669 
00670    while ((p = (TParticlePDG*) it.Next())) {
00671 
00672       // define decay channels for antiparticles
00673       if (p->PdgCode() < 0) {
00674          ap = GetParticle(-p->PdgCode());
00675          nch = ap->NDecayChannels();
00676          for (ich=0; ich<nch; ich++) {
00677             dc = ap->DecayChannel(ich);
00678             ndau = dc->NDaughters();
00679             for (int i=0; i<ndau; i++) {
00680                // conserve CPT
00681 
00682                code[i] = dc->DaughterPdgCode(i);
00683                daughter = GetParticle(code[i]);
00684                if (daughter->AntiParticle()) {
00685                   // this particle does have an
00686                   // antiparticle
00687                   code[i] = -code[i];
00688                }
00689             }
00690             p->AddDecayChannel(dc->MatrixElementCode(),
00691                                dc->BranchingRatio(),
00692                                dc->NDaughters(),
00693                                code);
00694          }
00695          p->SetAntiParticle(ap);
00696          ap->SetAntiParticle(p);
00697       }
00698    }
00699 
00700    fclose(file);
00701    return;
00702 }
00703 
00704 
00705 //______________________________________________________________________________
00706 void TDatabasePDG::Browse(TBrowser* b)
00707 {
00708    //browse data base
00709    if (fListOfClasses ) fListOfClasses->Browse(b);
00710 }
00711 
00712 
00713 //______________________________________________________________________________
00714 Int_t TDatabasePDG::WritePDGTable(const char *filename)
00715 {
00716    // write contents of the particle DB into a file
00717 
00718    if (fParticleList == 0) {
00719       Error("WritePDGTable","Do not have a valid PDG particle list;"
00720                             " consider loading it with ReadPDGTable first.");
00721       return -1;
00722    }
00723 
00724    FILE *file = fopen(filename,"w");
00725    if (file == 0) {
00726       Error("WritePDGTable","Could not open PDG particle file %s",filename);
00727       return -1;
00728    }
00729 
00730    fprintf(file,"#--------------------------------------------------------------------\n");
00731    fprintf(file,"#    i   NAME.............  KF AP   CLASS      Q        MASS     WIDTH  2*I+1 I3 2*S+1 FLVR TrkCod N(dec)\n");
00732    fprintf(file,"#--------------------------------------------------------------------\n");
00733 
00734    Int_t nparts=fParticleList->GetEntries();
00735    for(Int_t i=0;i<nparts;++i) {
00736       TParticlePDG *p = dynamic_cast<TParticlePDG*>(fParticleList->At(i));
00737       if(!p) continue;
00738 
00739       Int_t ich=i+1;
00740       Int_t kf=p->PdgCode();
00741       fprintf(file,"%5i %-20s %- 6i ", ich, p->GetName(), kf);
00742 
00743       Int_t anti=p->AntiParticle() ? 1:0;
00744       if(kf<0) {
00745          for(Int_t j=0;j<nparts;++j) {
00746             TParticlePDG *dummy = dynamic_cast<TParticlePDG*>(fParticleList->At(j));
00747             if(dummy==p->AntiParticle()) {
00748                anti=j+1;
00749                break;
00750             }
00751          }
00752          fprintf(file,"%i 0\n",anti);
00753          continue;
00754       }
00755 
00756       fprintf(file,"%i ",anti);
00757       fprintf(file,"%i ",100);
00758       fprintf(file,"%s ",p->ParticleClass());
00759       fprintf(file,"% i ",(Int_t)p->Charge());
00760       fprintf(file,"%.5le ",p->Mass());
00761       fprintf(file,"%.5le ",p->Width());
00762       fprintf(file,"%i ",(Int_t)p->Isospin());
00763       fprintf(file,"%i ",(Int_t)p->I3());
00764       fprintf(file,"%i ",(Int_t)p->Spin());
00765       fprintf(file,"%i ",-1);
00766       fprintf(file,"%i ",p->TrackingCode());
00767       Int_t nch=p->NDecayChannels();
00768       fprintf(file,"%i\n",nch);
00769       if(nch==0) {
00770          continue;
00771       }
00772       fprintf(file,"#----------------------------------------------------------------------\n");
00773       fprintf(file,"#    decay  type(PY6)    BR     Nd         daughters(codes, then names)\n");
00774       fprintf(file,"#----------------------------------------------------------------------\n");
00775       for(Int_t j=0;j<nch;++j) {
00776          TDecayChannel *dc=p->DecayChannel(j);
00777 
00778          fprintf(file,"%9i   ",dc->Number()+1);
00779          fprintf(file,"%3i   ",dc->MatrixElementCode());
00780          fprintf(file,"%.5le  ",dc->BranchingRatio());
00781          Int_t ndau=dc->NDaughters();
00782          fprintf(file,"%3i       ",ndau);
00783          for (int idau=0; idau<ndau; idau++) {
00784             fprintf(file,"%- 6i ",dc->DaughterPdgCode(idau));
00785          }
00786          for (int idau=0; idau<ndau; idau++) {
00787             TParticlePDG *dummy=GetParticle(dc->DaughterPdgCode(idau));
00788             if(dummy)
00789                fprintf(file,"%-10s ",dummy->GetName());
00790             else
00791                fprintf(file,"%-10s ","???");
00792          }
00793          fprintf(file,"\n");
00794       }
00795    }
00796    fclose(file);
00797    return nparts;
00798 }

Generated on Tue Jul 5 14:45:15 2011 for ROOT_528-00b_version by  doxygen 1.5.1