ROOT  6.06/08
Reference Guide
TDatabasePDG.cxx
Go to the documentation of this file.
1 // @(#)root/eg:$Id$
2 // Author: Pasha Murat 12/02/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 #include "RConfigure.h"
13 #include "TROOT.h"
14 #include "TEnv.h"
15 #include "THashList.h"
16 #include "TExMap.h"
17 #include "TSystem.h"
18 #include "TDatabasePDG.h"
19 #include "TDecayChannel.h"
20 #include "TParticlePDG.h"
21 #include <stdlib.h>
22 
23 
24 ////////////////////////////////////////////////////////////////////////
25 //
26 // Particle database manager class
27 //
28 // This manager creates a list of particles which by default is
29 // initialised from with the constants used by PYTHIA6 (plus some
30 // other particles added). See definition and the format of the default
31 // particle list in $ROOTSYS/etc/pdg_table.txt
32 //
33 // there are 2 ways of redefining the name of the file containing the
34 // particle properties
35 //
36 // 1. one can define the name in .rootrc file:
37 //
38 // Root.DatabasePDG: $(HOME)/my_pdg_table.txt
39 //
40 // 2. one can use TDatabasePDG::ReadPDGTable method explicitly:
41 //
42 // - TDatabasePDG *pdg = new TDatabasePDG();
43 // - pdg->ReadPDGtable(filename)
44 //
45 // See TParticlePDG for the description of a static particle properties.
46 // See TParticle for the description of a dynamic particle particle.
47 //
48 ////////////////////////////////////////////////////////////////////////
49 
51 
52 TDatabasePDG* TDatabasePDG::fgInstance = 0;
53 
54 ////////////////////////////////////////////////////////////////////////////////
55 /// Create PDG database. Initialization of the DB has to be done via explicit
56 /// call to ReadDataBasePDG (also done by GetParticle methods)
57 
58 TDatabasePDG::TDatabasePDG(): TNamed("PDGDB","The PDG particle data base")
59 {
60  fParticleList = 0;
61  fPdgMap = 0;
62  fListOfClasses = 0;
63  if (fgInstance) {
64  Warning("TDatabasePDG", "object already instantiated");
65  } else {
66  fgInstance = this;
67  gROOT->GetListOfSpecials()->Add(this);
68  }
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// Cleanup the PDG database.
73 
75 {
76  if (fParticleList) {
78  delete fParticleList; // this deletes all objects in the list
79  if (fPdgMap) delete fPdgMap;
80  }
81  // classes do not own particles...
82  if (fListOfClasses) {
84  delete fListOfClasses;
85  }
86  gROOT->GetListOfSpecials()->Remove(this);
87  fgInstance = 0;
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 ///static function
92 
94 {
95  return (fgInstance) ? (TDatabasePDG*) fgInstance : new TDatabasePDG();
96 }
97 
98 ////////////////////////////////////////////////////////////////////////////////
99 /// Build fPdgMap mapping pdg-code to particle.
100 ///
101 /// Initial size is set so as to be able to hold at least 600
102 /// particles: 521 in default table, ALICE adds 54 more.
103 /// To be revisited after LHC discovers SUSY.
104 
106 {
107  fPdgMap = new TExMap(4*TMath::Max(600, fParticleList->GetEntries())/3 + 3);
108  TIter next(fParticleList);
109  TParticlePDG *p;
110  while ((p = (TParticlePDG*)next())) {
111  fPdgMap->Add((Long_t)p->PdgCode(), (Long_t)p);
112  }
113 }
114 
115 ////////////////////////////////////////////////////////////////////////////////
116 ///
117 /// Particle definition normal constructor. If the particle is set to be
118 /// stable, the decay width parameter does have no meaning and can be set to
119 /// any value. The parameters granularity, LowerCutOff and HighCutOff are
120 /// used for the construction of the mean free path look up tables. The
121 /// granularity will be the number of logwise energy points for which the
122 /// mean free path will be calculated.
123 ///
124 
125 TParticlePDG* TDatabasePDG::AddParticle(const char *name, const char *title,
126  Double_t mass, Bool_t stable,
127  Double_t width, Double_t charge,
128  const char* ParticleClass,
129  Int_t PDGcode,
130  Int_t Anti,
131  Int_t TrackingCode)
132 {
133  TParticlePDG* old = GetParticle(PDGcode);
134 
135  if (old) {
136  printf(" *** TDatabasePDG::AddParticle: particle with PDGcode=%d already defined\n",PDGcode);
137  return 0;
138  }
139 
140  TParticlePDG* p = new TParticlePDG(name, title, mass, stable, width,
141  charge, ParticleClass, PDGcode, Anti,
142  TrackingCode);
143  fParticleList->Add(p);
144  if (fPdgMap)
145  fPdgMap->Add((Long_t)PDGcode, (Long_t)p);
146 
147  TParticleClassPDG* pclass = GetParticleClass(ParticleClass);
148 
149  if (!pclass) {
150  pclass = new TParticleClassPDG(ParticleClass);
151  fListOfClasses->Add(pclass);
152  }
153 
154  pclass->AddParticle(p);
155 
156  return p;
157 }
158 
159 ////////////////////////////////////////////////////////////////////////////////
160 /// assuming particle has already been defined
161 
163 {
164  TParticlePDG* old = GetParticle(PdgCode);
165 
166  if (old) {
167  printf(" *** TDatabasePDG::AddAntiParticle: can't redefine parameters\n");
168  return NULL;
169  }
170 
171  Int_t pdg_code = abs(PdgCode);
172  TParticlePDG* p = GetParticle(pdg_code);
173 
174  if (!p) {
175  printf(" *** TDatabasePDG::AddAntiParticle: particle with pdg code %d not known\n", pdg_code);
176  return NULL;
177  }
178 
179  TParticlePDG* ap = AddParticle(Name,
180  Name,
181  p->Mass(),
182  1,
183  p->Width(),
184  -p->Charge(),
185  p->ParticleClass(),
186  PdgCode,
187  1,
188  p->TrackingCode());
189  return ap;
190 }
191 
192 
193 ////////////////////////////////////////////////////////////////////////////////
194 ///
195 /// Get a pointer to the particle object according to the name given
196 ///
197 
199 {
200  if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
201 
203 // if (!def) {
204 // Error("GetParticle","No match for %s exists!",name);
205 // }
206  return def;
207 }
208 
209 ////////////////////////////////////////////////////////////////////////////////
210 ///
211 /// Get a pointer to the particle object according to the MC code number
212 ///
213 
215 {
216  if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
217  if (fPdgMap == 0) BuildPdgMap();
218 
219  return (TParticlePDG*) (Long_t)fPdgMap->GetValue((Long_t)PDGcode);
220 }
221 
222 ////////////////////////////////////////////////////////////////////////////////
223 /// Print contents of PDG database.
224 
225 void TDatabasePDG::Print(Option_t *option) const
226 {
227  if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
228 
229  TIter next(fParticleList);
230  TParticlePDG *p;
231  while ((p = (TParticlePDG *)next())) {
232  p->Print(option);
233  }
234 }
235 
236 ////////////////////////////////////////////////////////////////////////////////
237 /// Converts Geant3 particle codes to PDG convention. (Geant4 uses
238 /// PDG convention already)
239 /// Source: BaBar User Guide, Neil I. Geddes,
240 ///
241 /// See <A href="http://www.slac.stanford.edu/BFROOT/www/Computing/Environment/NewUser/htmlbug/node51.html">Conversion table</A>
242 /// with some fixes by PB, marked with (PB) below. Checked against
243 /// PDG listings from 2000.
244 ///
245 /// Paul Balm, Nov 19, 2001
246 
248 {
249  switch(Geant3number) {
250 
251  case 1 : return 22; // photon
252  case 25 : return -2112; // anti-neutron
253  case 2 : return -11; // e+
254  case 26 : return -3122; // anti-Lambda
255  case 3 : return 11; // e-
256  case 27 : return -3222; // Sigma-
257  case 4 : return 12; // e-neutrino (NB: flavour undefined by Geant)
258  case 28 : return -3212; // Sigma0
259  case 5 : return -13; // mu+
260  case 29 : return -3112; // Sigma+ (PB)*/
261  case 6 : return 13; // mu-
262  case 30 : return -3322; // Xi0
263  case 7 : return 111; // pi0
264  case 31 : return -3312; // Xi+
265  case 8 : return 211; // pi+
266  case 32 : return -3334; // Omega+ (PB)
267  case 9 : return -211; // pi-
268  case 33 : return -15; // tau+
269  case 10 : return 130; // K long
270  case 34 : return 15; // tau-
271  case 11 : return 321; // K+
272  case 35 : return 411; // D+
273  case 12 : return -321; // K-
274  case 36 : return -411; // D-
275  case 13 : return 2112; // n
276  case 37 : return 421; // D0
277  case 14 : return 2212; // p
278  case 38 : return -421; // D0
279  case 15 : return -2212; // anti-proton
280  case 39 : return 431; // Ds+
281  case 16 : return 310; // K short
282  case 40 : return -431; // anti Ds-
283  case 17 : return 221; // eta
284  case 41 : return 4122; // Lamba_c+
285  case 18 : return 3122; // Lambda
286  case 42 : return 24; // W+
287  case 19 : return 3222; // Sigma+
288  case 43 : return -24; // W-
289  case 20 : return 3212; // Sigma0
290  case 44 : return 23; // Z
291  case 21 : return 3112; // Sigma-
292  case 45 : return 0; // deuteron
293  case 22 : return 3322; // Xi0
294  case 46 : return 0; // triton
295  case 23 : return 3312; // Xi-
296  case 47 : return 0; // alpha
297  case 24 : return 3334; // Omega- (PB)
298  case 48 : return 0; // G nu ? PDG ID 0 is undefined
299 
300  default : return 0;
301 
302  }
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// Converts pdg code to geant3 id
307 
309 {
310  switch(pdgNumber) {
311 
312  case 22 : return 1; // photon
313  case -2112 : return 25; // anti-neutron
314  case -11 : return 2; // e+
315  case -3122 : return 26; // anti-Lambda
316  case 11 : return 3; // e-
317  case -3222 : return 27; // Sigma-
318  case 12 : return 4; // e-neutrino (NB: flavour undefined by Geant)
319  case -3212 : return 28; // Sigma0
320  case -13 : return 5; // mu+
321  case -3112 : return 29; // Sigma+ (PB)*/
322  case 13 : return 6; // mu-
323  case -3322 : return 30; // Xi0
324  case 111 : return 7; // pi0
325  case -3312 : return 31; // Xi+
326  case 211 : return 8; // pi+
327  case -3334 : return 32; // Omega+ (PB)
328  case -211 : return 9; // pi-
329  case -15 : return 33; // tau+
330  case 130 : return 10; // K long
331  case 15 : return 34; // tau-
332  case 321 : return 11; // K+
333  case 411 : return 35; // D+
334  case -321 : return 12; // K-
335  case -411 : return 36; // D-
336  case 2112 : return 13; // n
337  case 421 : return 37; // D0
338  case 2212 : return 14; // p
339  case -421 : return 38; // D0
340  case -2212 : return 15; // anti-proton
341  case 431 : return 39; // Ds+
342  case 310 : return 16; // K short
343  case -431 : return 40; // anti Ds-
344  case 221 : return 17; // eta
345  case 4122 : return 41; // Lamba_c+
346  case 3122 : return 18; // Lambda
347  case 24 : return 42; // W+
348  case 3222 : return 19; // Sigma+
349  case -24 : return 43; // W-
350  case 3212 : return 20; // Sigma0
351  case 23 : return 44; // Z
352  case 3112 : return 21; // Sigma-
353  case 3322 : return 22; // Xi0
354  case 3312 : return 23; // Xi-
355  case 3334 : return 24; // Omega- (PB)
356 
357  default : return 0;
358 
359  }
360 }
361 
362 ////////////////////////////////////////////////////////////////////////////////
363 ///
364 /// Converts the ISAJET Particle number into the PDG MC number
365 ///
366 
368 {
369  switch (isaNumber) {
370  case 1 : return 2; // UP .30000E+00 .67
371  case -1 : return -2; // UB .30000E+00 -.67
372  case 2 : return 1; // DN .30000E+00 -.33
373  case -2 : return -1; // DB .30000E+00 .33
374  case 3 : return 3; // ST .50000E+00 -.33
375  case -3 : return -3; // SB .50000E+00 .33
376  case 4 : return 4; // CH .16000E+01 .67
377  case -4 : return -4; // CB .16000E+01 -.67
378  case 5 : return 5; // BT .49000E+01 -.33
379  case -5 : return -5; // BB .49000E+01 .33
380  case 6 : return 6; // TP .17500E+03 .67
381  case -6 : return -6; // TB .17500E+03 -.67
382  case 9 : return 21; // GL 0. 0.00
383  case 80 : return 24; // W+ SIN2W=.23 1.00
384  case -80 : return -24; // W- SIN2W=.23 -1.00
385  case 90 : return 23; // Z0 SIN2W=.23 0.00
386  case 230 : return 311; // K0 .49767E+00 0.00
387  case -230 : return -311; // AK0 .49767E+00 0.00
388  case 330 : return 331; // ETAP .95760E+00 0.00
389  case 340 : return 0; // F- .20300E+01 -1.00
390  case -340 : return 0; // F+ .20300E+01 1.00
391  case 440 : return 441; // ETAC .29760E+01 0.00
392  case 111 : return 113; // RHO0 .77000E+00 0.00
393  case 121 : return 213; // RHO+ .77000E+00 1.00
394  case -121 : return -213; // RHO- .77000E+00 -1.00
395  case 221 : return 223; // OMEG .78260E+00 0.00
396  case 131 : return 323; // K*+ .88810E+00 1.00
397  case -131 : return -323; // K*- .88810E+00 -1.00
398  case 231 : return 313; // K*0 .89220E+00 0.00
399  case -231 : return -313; // AK*0 .89220E+00 0.00
400  case 331 : return 333; // PHI .10196E+01 0.00
401  case -140 : return 421; // D0
402  case 140 : return -421; // D0 bar
403  case 141 : return -423; // AD*0 .20060E+01 0.00
404  case -141 : return 423; // D*0 .20060E+01 0.00
405  case -240 : return -411; // D+
406  case 240 : return 411; // D-
407  case 241 : return -413; // D*- .20086E+01 -1.00
408  case -241 : return 413; // D*+ .20086E+01 1.00
409  case 341 : return 0; // F*- .21400E+01 -1.00
410  case -341 : return 0; // F*+ .21400E+01 1.00
411  case 441 : return 443; // JPSI .30970E+01 0.00
412 
413  // B-mesons, Bc still missing
414  case 250 : return 511; // B0
415  case -250 : return -511; // B0 bar
416  case 150 : return 521; // B+
417  case -150 : return -521; // B-
418  case 350 : return 531; // Bs 0
419  case -350 : return -531; // Bs bar
420  case 351 : return 533; // Bs* 0
421  case -351 : return -533; // Bs* bar
422  case 450 : return 541; // Bc +
423  case -450 : return -541; // Bc bar
424 
425  case 1140 : return 4222; // SC++ .24300E+01 2.00
426  case -1140 : return -4222; // ASC-- .24300E+01 -2.00
427  case 1240 : return 4212; // SC+ .24300E+01 1.00
428  case -1240 : return -4212; // ASC- .24300E+01 -1.00
429  case 2140 : return 4122; // LC+ .22600E+01 1.00
430  case -2140 : return -4122; // ALC- .22600E+01 -1.00
431  case 2240 : return 4112; // SC0 .24300E+01 0.00
432  case -2240 : return -4112; // ASC0 .24300E+01 0.00
433  case 1340 : return 0; // USC. .25000E+01 1.00
434  case -1340 : return 0; // AUSC. .25000E+01 -1.00
435  case 3140 : return 0; // SUC. .24000E+01 1.00
436  case -3140 : return 0; // ASUC. .24000E+01 -1.00
437  case 2340 : return 0; // DSC. .25000E+01 0.00
438  case -2340 : return 0; // ADSC. .25000E+01 0.00
439  case 3240 : return 0; // SDC. .24000E+01 0.00
440  case -3240 : return 0; // ASDC. .24000E+01 0.00
441  case 3340 : return 0; // SSC. .26000E+01 0.00
442  case -3340 : return 0; // ASSC. .26000E+01 0.00
443  case 1440 : return 0; // UCC. .35500E+01 2.00
444  case -1440 : return 0; // AUCC. .35500E+01 -2.00
445  case 2440 : return 0; // DCC. .35500E+01 1.00
446  case -2440 : return 0; // ADCC. .35500E+01 -1.00
447  case 3440 : return 0; // SCC. .37000E+01 1.00
448  case -3440 : return 0; // ASCC. .37000E+01 -1.00
449  case 1111 : return 2224; // DL++ .12320E+01 2.00
450  case -1111 : return -2224; // ADL-- .12320E+01 -2.00
451  case 1121 : return 2214; // DL+ .12320E+01 1.00
452  case -1121 : return -2214; // ADL- .12320E+01 -1.00
453  case 1221 : return 2114; // DL0 .12320E+01 0.00
454  case -1221 : return -2114; // ADL0 .12320E+01 0.00
455  case 2221 : return 1114; // DL- .12320E+01 -1.00
456  case -2221 : return -1114; // ADL+ .12320E+01 1.00
457  case 1131 : return 3224; // S*+ .13823E+01 1.00
458  case -1131 : return -3224; // AS*- .13823E+01 -1.00
459  case 1231 : return 3214; // S*0 .13820E+01 0.00
460  case -1231 : return -3214; // AS*0 .13820E+01 0.00
461  case 2231 : return 3114; // S*- .13875E+01 -1.00
462  case -2231 : return -3114; // AS*+ .13875E+01 1.00
463  case 1331 : return 3324; // XI*0 .15318E+01 0.00
464  case -1331 : return -3324; // AXI*0 .15318E+01 0.00
465  case 2331 : return 3314; // XI*- .15350E+01 -1.00
466  case -2331 : return -3314; // AXI*+ .15350E+01 1.00
467  case 3331 : return 3334; // OM- .16722E+01 -1.00
468  case -3331 : return -3334; // AOM+ .16722E+01 1.00
469  case 1141 : return 0; // UUC* .26300E+01 2.00
470  case -1141 : return 0; // AUUC* .26300E+01 -2.00
471  case 1241 : return 0; // UDC* .26300E+01 1.00
472  case -1241 : return 0; // AUDC* .26300E+01 -1.00
473  case 2241 : return 0; // DDC* .26300E+01 0.00
474  case -2241 : return 0; // ADDC* .26300E+01 0.00
475  case 1341 : return 0; // USC* .27000E+01 1.00
476  case -1341 : return 0; // AUSC* .27000E+01 -1.00
477  case 2341 : return 0; // DSC* .27000E+01 0.00
478  case -2341 : return 0; // ADSC* .27000E+01 0.00
479  case 3341 : return 0; // SSC* .28000E+01 0.00
480  case -3341 : return 0; // ASSC* .28000E+01 0.00
481  case 1441 : return 0; // UCC* .37500E+01 2.00
482  case -1441 : return 0; // AUCC* .37500E+01 -2.00
483  case 2441 : return 0; // DCC* .37500E+01 1.00
484  case -2441 : return 0; // ADCC* .37500E+01 -1.00
485  case 3441 : return 0; // SCC* .39000E+01 1.00
486  case -3441 : return 0; // ASCC* .39000E+01 -1.00
487  case 4441 : return 0; // CCC* .48000E+01 2.00
488  case -4441 : return 0; // ACCC* .48000E+01 -2.00
489  case 10 : return 22; // Photon
490  case 12 : return 11; // Electron
491  case -12 : return -11; // Positron
492  case 14 : return 13; // Muon-
493  case -14 : return -13; // Muon+
494  case 16 : return 15; // Tau-
495  case -16 : return -15; // Tau+
496  case 11 : return 12; // Neutrino e
497  case -11 : return -12; // Anti Neutrino e
498  case 13 : return 14; // Neutrino Muon
499  case -13 : return -14; // Anti Neutrino Muon
500  case 15 : return 16; // Neutrino Tau
501  case -15 : return -16; // Anti Neutrino Tau
502  case 110 : return 111; // Pion0
503  case 120 : return 211; // Pion+
504  case -120 : return -211; // Pion-
505  case 220 : return 221; // Eta
506  case 130 : return 321; // Kaon+
507  case -130 : return -321; // Kaon-
508  case -20 : return 130; // Kaon Long
509  case 20 : return 310; // Kaon Short
510 
511  // baryons
512  case 1120 : return 2212; // Proton
513  case -1120 : return -2212; // Anti Proton
514  case 1220 : return 2112; // Neutron
515  case -1220 : return -2112; // Anti Neutron
516  case 2130 : return 3122; // Lambda
517  case -2130 : return -3122; // Lambda bar
518  case 1130 : return 3222; // Sigma+
519  case -1130 : return -3222; // Sigma bar -
520  case 1230 : return 3212; // Sigma0
521  case -1230 : return -3212; // Sigma bar 0
522  case 2230 : return 3112; // Sigma-
523  case -2230 : return -3112; // Sigma bar +
524  case 1330 : return 3322; // Xi0
525  case -1330 : return -3322; // Xi bar 0
526  case 2330 : return 3312; // Xi-
527  case -2330 : return -3312; // Xi bar +
528  default : return 0; // isajet or pdg number does not exist
529  }
530 }
531 
532 ////////////////////////////////////////////////////////////////////////////////
533 /// read list of particles from a file
534 /// if the particle list does not exist, it is created, otherwise
535 /// particles are added to the existing list
536 /// See $ROOTSYS/etc/pdg_table.txt to see the file format
537 
538 void TDatabasePDG::ReadPDGTable(const char *FileName)
539 {
540  if (fParticleList == 0) {
541  fParticleList = new THashList;
543  }
544 
545  TString default_name;
546  const char *fn;
547 
548  if (!FileName[0]) {
549  default_name = "pdg_table.txt";
550  gSystem->PrependPathName(TROOT::GetEtcDir(), default_name);
551  fn = gEnv->GetValue("Root.DatabasePDG", default_name.Data());
552  } else {
553  fn = FileName;
554  }
555 
556  FILE* file = fopen(fn,"r");
557  if (file == 0) {
558  Error("ReadPDGTable","Could not open PDG particle file %s",fn);
559  return;
560  }
561 
562  char c[512];
563  Int_t class_number, anti, isospin, i3, spin, tracking_code;
564  Int_t ich, kf, nch, charge;
565  char name[30], class_name[30];
566  Double_t mass, width, branching_ratio;
567  Int_t dau[20];
568 
569  Int_t idecay, decay_type, flavor, ndau, stable;
570 
571  Int_t input;
572  while ( (input=getc(file)) != EOF) {
573  c[0] = input;
574  if (c[0] != '#') {
575  ungetc(c[0],file);
576  // read channel number
577  // coverity [secure_coding : FALSE]
578  if (fscanf(file,"%i",&ich)) {;}
579  // coverity [secure_coding : FALSE]
580  if (fscanf(file,"%s",name )) {;}
581  // coverity [secure_coding : FALSE]
582  if (fscanf(file,"%i",&kf )) {;}
583  // coverity [secure_coding : FALSE]
584  if (fscanf(file,"%i",&anti )) {;}
585 
586  if (kf < 0) {
587  AddAntiParticle(name,kf);
588  // nothing more on this line
589  if (fgets(c,200,file)) {;}
590  } else {
591  // coverity [secure_coding : FALSE]
592  if (fscanf(file,"%i",&class_number)) {;}
593  // coverity [secure_coding : FALSE]
594  if (fscanf(file,"%s",class_name)) {;}
595  // coverity [secure_coding : FALSE]
596  if (fscanf(file,"%i",&charge)) {;}
597  // coverity [secure_coding : FALSE]
598  if (fscanf(file,"%le",&mass)) {;}
599  // coverity [secure_coding : FALSE]
600  if (fscanf(file,"%le",&width)) {;}
601  // coverity [secure_coding : FALSE]
602  if (fscanf(file,"%i",&isospin)) {;}
603  // coverity [secure_coding : FALSE]
604  if (fscanf(file,"%i",&i3)) {;}
605  // coverity [secure_coding : FALSE]
606  if (fscanf(file,"%i",&spin)) {;}
607  // coverity [secure_coding : FALSE]
608  if (fscanf(file,"%i",&flavor)) {;}
609  // coverity [secure_coding : FALSE]
610  if (fscanf(file,"%i",&tracking_code)) {;}
611  // coverity [secure_coding : FALSE]
612  if (fscanf(file,"%i",&nch)) {;}
613  // nothing more on this line
614  if (fgets(c,200,file)) {;}
615  if (width > 1e-10) stable = 0;
616  else stable = 1;
617 
618  // create particle
619 
620  TParticlePDG* part = AddParticle(name,
621  name,
622  mass,
623  stable,
624  width,
625  charge,
626  class_name,
627  kf,
628  anti,
629  tracking_code);
630 
631  if (nch) {
632  // read in decay channels
633  ich = 0;
634  Int_t c_input = 0;
635  while ( ((c_input=getc(file)) != EOF) && (ich <nch)) {
636  c[0] = c_input;
637  if (c[0] != '#') {
638  ungetc(c[0],file);
639 
640  // coverity [secure_coding : FALSE]
641  if (fscanf(file,"%i",&idecay)) {;}
642  // coverity [secure_coding : FALSE]
643  if (fscanf(file,"%i",&decay_type)) {;}
644  // coverity [secure_coding : FALSE]
645  if (fscanf(file,"%le",&branching_ratio)) {;}
646  // coverity [secure_coding : FALSE]
647  if (fscanf(file,"%i",&ndau)) {;}
648  for (int idau=0; idau<ndau; idau++) {
649  // coverity [secure_coding : FALSE]
650  if (fscanf(file,"%i",&dau[idau])) {;}
651  }
652  // add decay channel
653 
654  if (part) part->AddDecayChannel(decay_type,branching_ratio,ndau,dau);
655  ich++;
656  }
657  // skip end of line
658  if (fgets(c,200,file)) {;}
659  }
660  }
661  }
662  } else {
663  // skip end of line
664  if (fgets(c,200,file)) {;}
665  }
666  }
667  // in the end loop over the antiparticles and
668  // define their decay lists
669  TIter it(fParticleList);
670 
671  Int_t code[20];
672  TParticlePDG *ap, *p, *daughter;
673  TDecayChannel *dc;
674 
675  while ((p = (TParticlePDG*) it.Next())) {
676 
677  // define decay channels for antiparticles
678  if (p->PdgCode() < 0) {
679  ap = GetParticle(-p->PdgCode());
680  if (!ap) continue;
681  nch = ap->NDecayChannels();
682  for (ich=0; ich<nch; ich++) {
683  dc = ap->DecayChannel(ich);
684  if (!dc) continue;
685  ndau = dc->NDaughters();
686  for (int i=0; i<ndau; i++) {
687  // conserve CPT
688 
689  code[i] = dc->DaughterPdgCode(i);
690  daughter = GetParticle(code[i]);
691  if (daughter && daughter->AntiParticle()) {
692  // this particle does have an
693  // antiparticle
694  code[i] = -code[i];
695  }
696  }
698  dc->BranchingRatio(),
699  dc->NDaughters(),
700  code);
701  }
702  p->SetAntiParticle(ap);
703  ap->SetAntiParticle(p);
704  }
705  }
706 
707  fclose(file);
708  return;
709 }
710 
711 
712 ////////////////////////////////////////////////////////////////////////////////
713 ///browse data base
714 
716 {
718 }
719 
720 
721 ////////////////////////////////////////////////////////////////////////////////
722 /// write contents of the particle DB into a file
723 
725 {
726  if (fParticleList == 0) {
727  Error("WritePDGTable","Do not have a valid PDG particle list;"
728  " consider loading it with ReadPDGTable first.");
729  return -1;
730  }
731 
732  FILE *file = fopen(filename,"w");
733  if (file == 0) {
734  Error("WritePDGTable","Could not open PDG particle file %s",filename);
735  return -1;
736  }
737 
738  fprintf(file,"#--------------------------------------------------------------------\n");
739  fprintf(file,"# i NAME............. KF AP CLASS Q MASS WIDTH 2*I+1 I3 2*S+1 FLVR TrkCod N(dec)\n");
740  fprintf(file,"#--------------------------------------------------------------------\n");
741 
742  Int_t nparts=fParticleList->GetEntries();
743  for(Int_t i=0;i<nparts;++i) {
744  TParticlePDG *p = dynamic_cast<TParticlePDG*>(fParticleList->At(i));
745  if(!p) continue;
746 
747  Int_t ich=i+1;
748  Int_t kf=p->PdgCode();
749  fprintf(file,"%5i %-20s %- 6i ", ich, p->GetName(), kf);
750 
751  Int_t anti=p->AntiParticle() ? 1:0;
752  if(kf<0) {
753  for(Int_t j=0;j<nparts;++j) {
754  TParticlePDG *dummy = dynamic_cast<TParticlePDG*>(fParticleList->At(j));
755  if(dummy==p->AntiParticle()) {
756  anti=j+1;
757  break;
758  }
759  }
760  fprintf(file,"%i 0\n",anti);
761  continue;
762  }
763 
764  fprintf(file,"%i ",anti);
765  fprintf(file,"%i ",100);
766  fprintf(file,"%s ",p->ParticleClass());
767  fprintf(file,"% i ",(Int_t)p->Charge());
768  fprintf(file,"%.5le ",p->Mass());
769  fprintf(file,"%.5le ",p->Width());
770  fprintf(file,"%i ",(Int_t)p->Isospin());
771  fprintf(file,"%i ",(Int_t)p->I3());
772  fprintf(file,"%i ",(Int_t)p->Spin());
773  fprintf(file,"%i ",-1);
774  fprintf(file,"%i ",p->TrackingCode());
775  Int_t nch=p->NDecayChannels();
776  fprintf(file,"%i\n",nch);
777  if(nch==0) {
778  continue;
779  }
780  fprintf(file,"#----------------------------------------------------------------------\n");
781  fprintf(file,"# decay type(PY6) BR Nd daughters(codes, then names)\n");
782  fprintf(file,"#----------------------------------------------------------------------\n");
783  for(Int_t j=0;j<nch;++j) {
784  TDecayChannel *dc=p->DecayChannel(j);
785  if (!dc) continue;
786  fprintf(file,"%9i ",dc->Number()+1);
787  fprintf(file,"%3i ",dc->MatrixElementCode());
788  fprintf(file,"%.5le ",dc->BranchingRatio());
789  Int_t ndau=dc->NDaughters();
790  fprintf(file,"%3i ",ndau);
791  for (int idau=0; idau<ndau; idau++) {
792  fprintf(file,"%- 6i ",dc->DaughterPdgCode(idau));
793  }
794  for (int idau=0; idau<ndau; idau++) {
796  if(dummy)
797  fprintf(file,"%-10s ",dummy->GetName());
798  else
799  fprintf(file,"%-10s ","???");
800  }
801  fprintf(file,"\n");
802  }
803  }
804  fclose(file);
805  return nparts;
806 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
void Add(ULong64_t hash, Long64_t key, Long64_t value)
Add an (key,value) pair to the table. The key should be unique.
Definition: TExMap.cxx:85
An array of TObjects.
Definition: TObjArray.h:39
Int_t PdgCode() const
Definition: TParticlePDG.h:70
Int_t NDaughters()
Definition: TDecayChannel.h:48
Double_t Width() const
Definition: TParticlePDG.h:74
const char Option_t
Definition: RtypesCore.h:62
TDatabasePDG()
Create PDG database.
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:328
const char * ParticleClass() const
Definition: TParticlePDG.h:86
virtual ~TDatabasePDG()
Cleanup the PDG database.
virtual Int_t GetEntries() const
Definition: TCollection.h:92
Int_t NDecayChannels() const
Definition: TParticlePDG.h:90
static const char * filename()
void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: THashList.cxx:183
#define gROOT
Definition: TROOT.h:352
Basic string class.
Definition: TString.h:137
Double_t I3() const
Definition: TParticlePDG.h:78
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TExMap * fPdgMap
Definition: TDatabasePDG.h:31
TObject * FindObject(const char *name) const
Find object using its name.
Definition: THashList.cxx:212
TParticleClassPDG * GetParticleClass(const char *name)
Definition: TDatabasePDG.h:72
TDecayChannel * DecayChannel(Int_t i)
return pointer to decay channel object at index i
const char * Name
Definition: TXMLSetup.cxx:67
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)
Particle definition normal constructor.
static TDatabasePDG * Instance()
static function
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition: THashList.h:36
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
virtual TParticlePDG * AddAntiParticle(const char *Name, Int_t PdgCode)
assuming particle has already been defined
TParticlePDG * AntiParticle()
Definition: TParticlePDG.h:98
virtual void Print(Option_t *opt="") const
Print the entire information of this kind of particle.
Int_t DaughterPdgCode(Int_t i)
Definition: TDecayChannel.h:50
static Vc_ALWAYS_INLINE Vector< T > abs(const Vector< T > &x)
Definition: vector.h:450
virtual const char * PrependPathName(const char *dir, TString &name)
Concatenate a directory and a file name.
Definition: TSystem.cxx:1054
Int_t MatrixElementCode()
Definition: TDecayChannel.h:47
virtual void Browse(TBrowser *b)
browse data base
virtual Int_t ConvertIsajetToPdg(Int_t isaNumber) const
Converts the ISAJET Particle number into the PDG MC number.
void SetAntiParticle(TParticlePDG *ap)
Definition: TParticlePDG.h:103
TObjArray * fListOfClasses
Definition: TDatabasePDG.h:30
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:41
Long64_t GetValue(ULong64_t hash, Long64_t key)
Return the value belonging to specified key and hash value.
Definition: TExMap.cxx:171
R__EXTERN TSystem * gSystem
Definition: TSystem.h:549
THashList * fParticleList
Definition: TDatabasePDG.h:29
virtual Int_t GetValue(const char *name, Int_t dflt)
Returns the integer value for a resource.
Definition: TEnv.cxx:480
Double_t Spin() const
Definition: TParticlePDG.h:76
TObject * Next()
Definition: TCollection.h:158
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
Int_t TrackingCode() const
Definition: TParticlePDG.h:94
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:310
Double_t Charge() const
Definition: TParticlePDG.h:72
void Warning(const char *location, const char *msgfmt,...)
Int_t Number()
Definition: TDecayChannel.h:46
void AddParticle(TObject *p)
long Long_t
Definition: RtypesCore.h:50
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
virtual 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
R__EXTERN TEnv * gEnv
Definition: TEnv.h:174
static const TString & GetEtcDir()
Get the sysconfig directory in the installation. Static utility function.
Definition: TROOT.cxx:2601
static RooMathCoreReg dummy
static TDatabasePDG * fgInstance
Definition: TDatabasePDG.h:28
void Browse(TBrowser *b)
Browse this collection (called by TBrowser).
virtual Int_t ConvertGeant3ToPdg(Int_t Geant3Number) const
Converts Geant3 particle codes to PDG convention.
TParticlePDG * GetParticle(Int_t pdgCode) const
Get a pointer to the particle object according to the MC code number.
virtual Int_t WritePDGTable(const char *filename)
write contents of the particle DB into a file
#define name(a, b)
Definition: linkTestLib0.cpp:5
void BuildPdgMap() const
Build fPdgMap mapping pdg-code to particle.
virtual void Add(TObject *obj)
Definition: TList.h:81
Double_t Isospin() const
Definition: TParticlePDG.h:77
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
virtual void Print(Option_t *opt="") const
Print contents of PDG database.
#define NULL
Definition: Rtypes.h:82
Double_t BranchingRatio()
Definition: TDecayChannel.h:49
void Add(TObject *obj)
Definition: TObjArray.h:75
Double_t Mass() const
Definition: TParticlePDG.h:71
virtual Int_t ConvertPdgToGeant3(Int_t pdgNumber) const
Converts pdg code to geant3 id.
This class stores a (key,value) pair using an external hash.
Definition: TExMap.h:35
Int_t AddDecayChannel(Int_t Type, Double_t BranchingRatio, Int_t NDaughters, Int_t *DaughterPdgCode)
add new decay channel, Particle owns those...
const char * Data() const
Definition: TString.h:349