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