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