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 * TMath::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,
158 Double_t width, Double_t charge,
159 const char* ParticleClass,
160 Int_t PDGcode,
161 Int_t Anti,
162 Int_t TrackingCode)
163{
164 TParticlePDG* old = GetParticle(PDGcode);
165
166 if (old) {
167 Warning("AddParticle", "Particle with PDGcode=%d already defined", PDGcode);
168 return 0;
169 }
170
171 TParticlePDG* p = new TParticlePDG(name, title, mass, stable, width,
172 charge, ParticleClass, PDGcode, Anti,
173 TrackingCode);
175 if (fPdgMap)
176 fPdgMap->Add((Long64_t)PDGcode, (Long64_t)p);
177
178 TParticleClassPDG* pclass = GetParticleClass(ParticleClass);
179
180 if (!pclass) {
181 pclass = new TParticleClassPDG(ParticleClass);
182 fListOfClasses->Add(pclass);
183 }
184
185 pclass->AddParticle(p);
186
187 return p;
188}
189
190////////////////////////////////////////////////////////////////////////////////
191/// assuming particle has already been defined
192
194{
195 TParticlePDG* old = GetParticle(PdgCode);
196
197 if (old) {
198 Warning("AddAntiParticle", "Can't redefine parameters");
199 return NULL;
200 }
201
202 Int_t pdg_code = abs(PdgCode);
203 TParticlePDG* p = GetParticle(pdg_code);
204
205 if (!p) {
206 Warning("AddAntiParticle", "Particle with pdg code %d not known", pdg_code);
207 return NULL;
208 }
209
210 TParticlePDG* ap = AddParticle(Name,
211 Name,
212 p->Mass(),
213 1,
214 p->Width(),
215 -p->Charge(),
216 p->ParticleClass(),
217 PdgCode,
218 1,
219 p->TrackingCode());
220 return ap;
221}
222
223
224////////////////////////////////////////////////////////////////////////////////
225///
226/// Get a pointer to the particle object according to the name given
227///
228
230{
231 {
233 if (fParticleList == 0)
234 ((TDatabasePDG *)this)->ReadPDGTableImpl("", false);
235 }
237// if (!def) {
238// Error("GetParticle","No match for %s exists!",name);
239// }
240 return def;
241}
242
243////////////////////////////////////////////////////////////////////////////////
244///
245/// Get a pointer to the particle object according to the MC code number
246///
247
249{
250 {
252 if (fParticleList == 0)
253 ((TDatabasePDG *)this)->ReadPDGTableImpl("", false);
254 if (fPdgMap == 0)
255 BuildPdgMap();
256 }
257
258 return (TParticlePDG *)fPdgMap->GetValue((Long64_t)PDGcode);
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// Print contents of PDG database.
263
265{
266 if (fParticleList == 0) ((TDatabasePDG*)this)->ReadPDGTable();
267
268 TIter next(fParticleList);
270 while ((p = (TParticlePDG *)next())) {
271 p->Print(option);
272 }
273}
274
275////////////////////////////////////////////////////////////////////////////////
276/// Converts Geant3 particle codes to PDG convention. (Geant4 uses
277/// PDG convention already)
278/// Source: BaBar User Guide, Neil I. Geddes,
279///
280/// see <A href="http://www.slac.stanford.edu/BFROOT/www/Computing/Environment/NewUser/htmlbug/node51.html"> Conversion table</A>
281///
282/// with some fixes by PB, marked with (PB) below. Checked against
283/// PDG listings from 2000.
284///
285/// Paul Balm, Nov 19, 2001
286
288{
289
290
291 switch(Geant3number) {
292
293 case 1 : return 22; // photon
294 case 25 : return -2112; // anti-neutron
295 case 2 : return -11; // e+
296 case 26 : return -3122; // anti-Lambda
297 case 3 : return 11; // e-
298 case 27 : return -3222; // Sigma-
299 case 4 : return 12; // e-neutrino (NB: flavour undefined by Geant)
300 case 28 : return -3212; // Sigma0
301 case 5 : return -13; // mu+
302 case 29 : return -3112; // Sigma+ (PB)*/
303 case 6 : return 13; // mu-
304 case 30 : return -3322; // Xi0
305 case 7 : return 111; // pi0
306 case 31 : return -3312; // Xi+
307 case 8 : return 211; // pi+
308 case 32 : return -3334; // Omega+ (PB)
309 case 9 : return -211; // pi-
310 case 33 : return -15; // tau+
311 case 10 : return 130; // K long
312 case 34 : return 15; // tau-
313 case 11 : return 321; // K+
314 case 35 : return 411; // D+
315 case 12 : return -321; // K-
316 case 36 : return -411; // D-
317 case 13 : return 2112; // n
318 case 37 : return 421; // D0
319 case 14 : return 2212; // p
320 case 38 : return -421; // D0
321 case 15 : return -2212; // anti-proton
322 case 39 : return 431; // Ds+
323 case 16 : return 310; // K short
324 case 40 : return -431; // anti Ds-
325 case 17 : return 221; // eta
326 case 41 : return 4122; // Lamba_c+
327 case 18 : return 3122; // Lambda
328 case 42 : return 24; // W+
329 case 19 : return 3222; // Sigma+
330 case 43 : return -24; // W-
331 case 20 : return 3212; // Sigma0
332 case 44 : return 23; // Z
333 case 21 : return 3112; // Sigma-
334 case 45 : return 0; // deuteron
335 case 22 : return 3322; // Xi0
336 case 46 : return 0; // triton
337 case 23 : return 3312; // Xi-
338 case 47 : return 0; // alpha
339 case 24 : return 3334; // Omega- (PB)
340 case 48 : return 0; // G nu ? PDG ID 0 is undefined
341
342 default : return 0;
343
344 }
345}
346
347////////////////////////////////////////////////////////////////////////////////
348/// Converts pdg code to geant3 id
349
351{
352 switch(pdgNumber) {
353
354 case 22 : return 1; // photon
355 case -2112 : return 25; // anti-neutron
356 case -11 : return 2; // e+
357 case -3122 : return 26; // anti-Lambda
358 case 11 : return 3; // e-
359 case -3222 : return 27; // Sigma-
360 case 12 : return 4; // e-neutrino (NB: flavour undefined by Geant)
361 case -3212 : return 28; // Sigma0
362 case -13 : return 5; // mu+
363 case -3112 : return 29; // Sigma+ (PB)*/
364 case 13 : return 6; // mu-
365 case -3322 : return 30; // Xi0
366 case 111 : return 7; // pi0
367 case -3312 : return 31; // Xi+
368 case 211 : return 8; // pi+
369 case -3334 : return 32; // Omega+ (PB)
370 case -211 : return 9; // pi-
371 case -15 : return 33; // tau+
372 case 130 : return 10; // K long
373 case 15 : return 34; // tau-
374 case 321 : return 11; // K+
375 case 411 : return 35; // D+
376 case -321 : return 12; // K-
377 case -411 : return 36; // D-
378 case 2112 : return 13; // n
379 case 421 : return 37; // D0
380 case 2212 : return 14; // p
381 case -421 : return 38; // D0
382 case -2212 : return 15; // anti-proton
383 case 431 : return 39; // Ds+
384 case 310 : return 16; // K short
385 case -431 : return 40; // anti Ds-
386 case 221 : return 17; // eta
387 case 4122 : return 41; // Lamba_c+
388 case 3122 : return 18; // Lambda
389 case 24 : return 42; // W+
390 case 3222 : return 19; // Sigma+
391 case -24 : return 43; // W-
392 case 3212 : return 20; // Sigma0
393 case 23 : return 44; // Z
394 case 3112 : return 21; // Sigma-
395 case 3322 : return 22; // Xi0
396 case 3312 : return 23; // Xi-
397 case 3334 : return 24; // Omega- (PB)
398
399 default : return 0;
400
401 }
402}
403
404////////////////////////////////////////////////////////////////////////////////
405///
406/// Converts the ISAJET Particle number into the PDG MC number
407///
408
410{
411 switch (isaNumber) {
412 case 1 : return 2; // UP .30000E+00 .67
413 case -1 : return -2; // UB .30000E+00 -.67
414 case 2 : return 1; // DN .30000E+00 -.33
415 case -2 : return -1; // DB .30000E+00 .33
416 case 3 : return 3; // ST .50000E+00 -.33
417 case -3 : return -3; // SB .50000E+00 .33
418 case 4 : return 4; // CH .16000E+01 .67
419 case -4 : return -4; // CB .16000E+01 -.67
420 case 5 : return 5; // BT .49000E+01 -.33
421 case -5 : return -5; // BB .49000E+01 .33
422 case 6 : return 6; // TP .17500E+03 .67
423 case -6 : return -6; // TB .17500E+03 -.67
424 case 9 : return 21; // GL 0. 0.00
425 case 80 : return 24; // W+ SIN2W=.23 1.00
426 case -80 : return -24; // W- SIN2W=.23 -1.00
427 case 90 : return 23; // Z0 SIN2W=.23 0.00
428 case 230 : return 311; // K0 .49767E+00 0.00
429 case -230 : return -311; // AK0 .49767E+00 0.00
430 case 330 : return 331; // ETAP .95760E+00 0.00
431 case 340 : return 0; // F- .20300E+01 -1.00
432 case -340 : return 0; // F+ .20300E+01 1.00
433 case 440 : return 441; // ETAC .29760E+01 0.00
434 case 111 : return 113; // RHO0 .77000E+00 0.00
435 case 121 : return 213; // RHO+ .77000E+00 1.00
436 case -121 : return -213; // RHO- .77000E+00 -1.00
437 case 221 : return 223; // OMEG .78260E+00 0.00
438 case 131 : return 323; // K*+ .88810E+00 1.00
439 case -131 : return -323; // K*- .88810E+00 -1.00
440 case 231 : return 313; // K*0 .89220E+00 0.00
441 case -231 : return -313; // AK*0 .89220E+00 0.00
442 case 331 : return 333; // PHI .10196E+01 0.00
443 case -140 : return 421; // D0
444 case 140 : return -421; // D0 bar
445 case 141 : return -423; // AD*0 .20060E+01 0.00
446 case -141 : return 423; // D*0 .20060E+01 0.00
447 case -240 : return -411; // D+
448 case 240 : return 411; // D-
449 case 241 : return -413; // D*- .20086E+01 -1.00
450 case -241 : return 413; // D*+ .20086E+01 1.00
451 case 341 : return 0; // F*- .21400E+01 -1.00
452 case -341 : return 0; // F*+ .21400E+01 1.00
453 case 441 : return 443; // JPSI .30970E+01 0.00
454
455 // B-mesons, Bc still missing
456 case 250 : return 511; // B0
457 case -250 : return -511; // B0 bar
458 case 150 : return 521; // B+
459 case -150 : return -521; // B-
460 case 350 : return 531; // Bs 0
461 case -350 : return -531; // Bs bar
462 case 351 : return 533; // Bs* 0
463 case -351 : return -533; // Bs* bar
464 case 450 : return 541; // Bc +
465 case -450 : return -541; // Bc bar
466
467 case 1140 : return 4222; // SC++ .24300E+01 2.00
468 case -1140 : return -4222; // ASC-- .24300E+01 -2.00
469 case 1240 : return 4212; // SC+ .24300E+01 1.00
470 case -1240 : return -4212; // ASC- .24300E+01 -1.00
471 case 2140 : return 4122; // LC+ .22600E+01 1.00
472 case -2140 : return -4122; // ALC- .22600E+01 -1.00
473 case 2240 : return 4112; // SC0 .24300E+01 0.00
474 case -2240 : return -4112; // ASC0 .24300E+01 0.00
475 case 1340 : return 0; // USC. .25000E+01 1.00
476 case -1340 : return 0; // AUSC. .25000E+01 -1.00
477 case 3140 : return 0; // SUC. .24000E+01 1.00
478 case -3140 : return 0; // ASUC. .24000E+01 -1.00
479 case 2340 : return 0; // DSC. .25000E+01 0.00
480 case -2340 : return 0; // ADSC. .25000E+01 0.00
481 case 3240 : return 0; // SDC. .24000E+01 0.00
482 case -3240 : return 0; // ASDC. .24000E+01 0.00
483 case 3340 : return 0; // SSC. .26000E+01 0.00
484 case -3340 : return 0; // ASSC. .26000E+01 0.00
485 case 1440 : return 0; // UCC. .35500E+01 2.00
486 case -1440 : return 0; // AUCC. .35500E+01 -2.00
487 case 2440 : return 0; // DCC. .35500E+01 1.00
488 case -2440 : return 0; // ADCC. .35500E+01 -1.00
489 case 3440 : return 0; // SCC. .37000E+01 1.00
490 case -3440 : return 0; // ASCC. .37000E+01 -1.00
491 case 1111 : return 2224; // DL++ .12320E+01 2.00
492 case -1111 : return -2224; // ADL-- .12320E+01 -2.00
493 case 1121 : return 2214; // DL+ .12320E+01 1.00
494 case -1121 : return -2214; // ADL- .12320E+01 -1.00
495 case 1221 : return 2114; // DL0 .12320E+01 0.00
496 case -1221 : return -2114; // ADL0 .12320E+01 0.00
497 case 2221 : return 1114; // DL- .12320E+01 -1.00
498 case -2221 : return -1114; // ADL+ .12320E+01 1.00
499 case 1131 : return 3224; // S*+ .13823E+01 1.00
500 case -1131 : return -3224; // AS*- .13823E+01 -1.00
501 case 1231 : return 3214; // S*0 .13820E+01 0.00
502 case -1231 : return -3214; // AS*0 .13820E+01 0.00
503 case 2231 : return 3114; // S*- .13875E+01 -1.00
504 case -2231 : return -3114; // AS*+ .13875E+01 1.00
505 case 1331 : return 3324; // XI*0 .15318E+01 0.00
506 case -1331 : return -3324; // AXI*0 .15318E+01 0.00
507 case 2331 : return 3314; // XI*- .15350E+01 -1.00
508 case -2331 : return -3314; // AXI*+ .15350E+01 1.00
509 case 3331 : return 3334; // OM- .16722E+01 -1.00
510 case -3331 : return -3334; // AOM+ .16722E+01 1.00
511 case 1141 : return 0; // UUC* .26300E+01 2.00
512 case -1141 : return 0; // AUUC* .26300E+01 -2.00
513 case 1241 : return 0; // UDC* .26300E+01 1.00
514 case -1241 : return 0; // AUDC* .26300E+01 -1.00
515 case 2241 : return 0; // DDC* .26300E+01 0.00
516 case -2241 : return 0; // ADDC* .26300E+01 0.00
517 case 1341 : return 0; // USC* .27000E+01 1.00
518 case -1341 : return 0; // AUSC* .27000E+01 -1.00
519 case 2341 : return 0; // DSC* .27000E+01 0.00
520 case -2341 : return 0; // ADSC* .27000E+01 0.00
521 case 3341 : return 0; // SSC* .28000E+01 0.00
522 case -3341 : return 0; // ASSC* .28000E+01 0.00
523 case 1441 : return 0; // UCC* .37500E+01 2.00
524 case -1441 : return 0; // AUCC* .37500E+01 -2.00
525 case 2441 : return 0; // DCC* .37500E+01 1.00
526 case -2441 : return 0; // ADCC* .37500E+01 -1.00
527 case 3441 : return 0; // SCC* .39000E+01 1.00
528 case -3441 : return 0; // ASCC* .39000E+01 -1.00
529 case 4441 : return 0; // CCC* .48000E+01 2.00
530 case -4441 : return 0; // ACCC* .48000E+01 -2.00
531 case 10 : return 22; // Photon
532 case 12 : return 11; // Electron
533 case -12 : return -11; // Positron
534 case 14 : return 13; // Muon-
535 case -14 : return -13; // Muon+
536 case 16 : return 15; // Tau-
537 case -16 : return -15; // Tau+
538 case 11 : return 12; // Neutrino e
539 case -11 : return -12; // Anti Neutrino e
540 case 13 : return 14; // Neutrino Muon
541 case -13 : return -14; // Anti Neutrino Muon
542 case 15 : return 16; // Neutrino Tau
543 case -15 : return -16; // Anti Neutrino Tau
544 case 110 : return 111; // Pion0
545 case 120 : return 211; // Pion+
546 case -120 : return -211; // Pion-
547 case 220 : return 221; // Eta
548 case 130 : return 321; // Kaon+
549 case -130 : return -321; // Kaon-
550 case -20 : return 130; // Kaon Long
551 case 20 : return 310; // Kaon Short
552
553 // baryons
554 case 1120 : return 2212; // Proton
555 case -1120 : return -2212; // Anti Proton
556 case 1220 : return 2112; // Neutron
557 case -1220 : return -2112; // Anti Neutron
558 case 2130 : return 3122; // Lambda
559 case -2130 : return -3122; // Lambda bar
560 case 1130 : return 3222; // Sigma+
561 case -1130 : return -3222; // Sigma bar -
562 case 1230 : return 3212; // Sigma0
563 case -1230 : return -3212; // Sigma bar 0
564 case 2230 : return 3112; // Sigma-
565 case -2230 : return -3112; // Sigma bar +
566 case 1330 : return 3322; // Xi0
567 case -1330 : return -3322; // Xi bar 0
568 case 2330 : return 3312; // Xi-
569 case -2330 : return -3312; // Xi bar +
570 default : return 0; // isajet or pdg number does not exist
571 }
572}
573
574void TDatabasePDG::ReadPDGTableImpl(const char *FileName, bool isParticleListInitializedWhenInvoking)
575{
576
578 // To avoid to rebuild the particle list from multiple threads, we check if at
579 // the moment of invoked ReadPDGTableImpl, the list was initialized. If it was not, and now it is, another
580 // thread filled it. Doing so again would result in an error.
581 // However, if the particle list was not initialized when ReadPDGTableImpl was invoked and it is still empty,
582 // we fill it.
583 // With the protection described above, ReadPDGTable can be invoked multiple times, for example once
584 // automatically by the system lazily and a second time by the user to complement the information
585 // already in memory with a second input file.
586 if (!isParticleListInitializedWhenInvoking && fParticleList) {
587 return;
588 }
589
590 if (fParticleList == 0) {
593 }
594
595 TString default_name;
596 const char *fn;
597
598 if (!FileName[0]) {
599 default_name = "pdg_table.txt";
600 gSystem->PrependPathName(TROOT::GetEtcDir(), default_name);
601 fn = gEnv->GetValue("Root.DatabasePDG", default_name.Data());
602 } else {
603 fn = FileName;
604 }
605
606 FILE* file = fopen(fn,"r");
607 if (file == 0) {
608 Error("ReadPDGTable","Could not open PDG particle file %s",fn);
609 return;
610 }
611
612 char c[512];
613 Int_t class_number, anti, isospin, i3, spin, tracking_code;
614 Int_t ich, kf, nch, charge;
615 char name[30], class_name[30];
616 Double_t mass, width, branching_ratio;
617 Int_t dau[20];
618
619 Int_t idecay, decay_type, flavor, ndau, stable;
620
621 Int_t input;
622 while ( (input=getc(file)) != EOF) {
623 c[0] = input;
624 if (c[0] != '#') {
625 ungetc(c[0],file);
626 // read channel number
627 // coverity [secure_coding : FALSE]
628 if (fscanf(file,"%i",&ich)) {;}
629 // coverity [secure_coding : FALSE]
630 if (fscanf(file,"%s",name )) {;}
631 // coverity [secure_coding : FALSE]
632 if (fscanf(file,"%i",&kf )) {;}
633 // coverity [secure_coding : FALSE]
634 if (fscanf(file,"%i",&anti )) {;}
635
636 if (kf < 0) {
638 // nothing more on this line
639 if (fgets(c,200,file)) {;}
640 } else {
641 // coverity [secure_coding : FALSE]
642 if (fscanf(file,"%i",&class_number)) {;}
643 // coverity [secure_coding : FALSE]
644 if (fscanf(file,"%s",class_name)) {;}
645 // coverity [secure_coding : FALSE]
646 if (fscanf(file,"%i",&charge)) {;}
647 // coverity [secure_coding : FALSE]
648 if (fscanf(file,"%le",&mass)) {;}
649 // coverity [secure_coding : FALSE]
650 if (fscanf(file,"%le",&width)) {;}
651 // coverity [secure_coding : FALSE]
652 if (fscanf(file,"%i",&isospin)) {;}
653 // coverity [secure_coding : FALSE]
654 if (fscanf(file,"%i",&i3)) {;}
655 // coverity [secure_coding : FALSE]
656 if (fscanf(file,"%i",&spin)) {;}
657 // coverity [secure_coding : FALSE]
658 if (fscanf(file,"%i",&flavor)) {;}
659 // coverity [secure_coding : FALSE]
660 if (fscanf(file,"%i",&tracking_code)) {;}
661 // coverity [secure_coding : FALSE]
662 if (fscanf(file,"%i",&nch)) {;}
663 // nothing more on this line
664 if (fgets(c,200,file)) {;}
665 if (width > 1e-10) stable = 0;
666 else stable = 1;
667
668 // create particle
669
671 name,
672 mass,
673 stable,
674 width,
675 charge,
676 class_name,
677 kf,
678 anti,
679 tracking_code);
680
681 if (nch) {
682 // read in decay channels
683 ich = 0;
684 Int_t c_input = 0;
685 while ( ((c_input=getc(file)) != EOF) && (ich <nch)) {
686 c[0] = c_input;
687 if (c[0] != '#') {
688 ungetc(c[0],file);
689
690 // coverity [secure_coding : FALSE]
691 if (fscanf(file,"%i",&idecay)) {;}
692 // coverity [secure_coding : FALSE]
693 if (fscanf(file,"%i",&decay_type)) {;}
694 // coverity [secure_coding : FALSE]
695 if (fscanf(file,"%le",&branching_ratio)) {;}
696 // coverity [secure_coding : FALSE]
697 if (fscanf(file,"%i",&ndau)) {;}
698 for (int idau=0; idau<ndau; idau++) {
699 // coverity [secure_coding : FALSE]
700 if (fscanf(file,"%i",&dau[idau])) {;}
701 }
702 // add decay channel
703
704 if (part) part->AddDecayChannel(decay_type,branching_ratio,ndau,dau);
705 ich++;
706 }
707 // skip end of line
708 if (fgets(c,200,file)) {;}
709 }
710 }
711 }
712 } else {
713 // skip end of line
714 if (fgets(c,200,file)) {;}
715 }
716 }
717 // in the end loop over the antiparticles and
718 // define their decay lists
720
721 Int_t code[20];
722 TParticlePDG *ap, *p, *daughter;
723 TDecayChannel *dc;
724
725 while ((p = (TParticlePDG*) it.Next())) {
726
727 // define decay channels for antiparticles
728 if (p->PdgCode() < 0) {
729 ap = GetParticle(-p->PdgCode());
730 if (!ap) continue;
731 nch = ap->NDecayChannels();
732 for (ich=0; ich<nch; ich++) {
733 dc = ap->DecayChannel(ich);
734 if (!dc) continue;
735 ndau = dc->NDaughters();
736 for (int i=0; i<ndau; i++) {
737 // conserve CPT
738
739 code[i] = dc->DaughterPdgCode(i);
740 daughter = GetParticle(code[i]);
741 if (daughter && daughter->AntiParticle()) {
742 // this particle does have an
743 // antiparticle
744 code[i] = -code[i];
745 }
746 }
747 p->AddDecayChannel(dc->MatrixElementCode(),
748 dc->BranchingRatio(),
749 dc->NDaughters(),
750 code);
751 }
752 p->SetAntiParticle(ap);
753 ap->SetAntiParticle(p);
754 }
755 }
756
757 fclose(file);
758 return;
759}
760
761////////////////////////////////////////////////////////////////////////////////
762/// read list of particles from a file
763/// if the particle list does not exist, it is created, otherwise
764/// particles are added to the existing list
765/// See $ROOTSYS/etc/pdg_table.txt to see the file format
766
767void TDatabasePDG::ReadPDGTable(const char *FileName)
768{
769 ReadPDGTableImpl(FileName, nullptr != fParticleList);
770}
771
772////////////////////////////////////////////////////////////////////////////////
773/// browse data base
774
776{
778}
779
780
781////////////////////////////////////////////////////////////////////////////////
782/// write contents of the particle DB into a file
783
785{
786 if (fParticleList == 0) {
787 Error("WritePDGTable","Do not have a valid PDG particle list;"
788 " consider loading it with ReadPDGTable first.");
789 return -1;
790 }
791
792 FILE *file = fopen(filename,"w");
793 if (file == 0) {
794 Error("WritePDGTable","Could not open PDG particle file %s",filename);
795 return -1;
796 }
797
798 fprintf(file,"#--------------------------------------------------------------------\n");
799 fprintf(file,"# i NAME............. KF AP CLASS Q MASS WIDTH 2*I+1 I3 2*S+1 FLVR TrkCod N(dec)\n");
800 fprintf(file,"#--------------------------------------------------------------------\n");
801
803 for(Int_t i=0;i<nparts;++i) {
804 TParticlePDG *p = dynamic_cast<TParticlePDG*>(fParticleList->At(i));
805 if(!p) continue;
806
807 Int_t ich=i+1;
808 Int_t kf=p->PdgCode();
809 fprintf(file,"%5i %-20s %- 6i ", ich, p->GetName(), kf);
810
811 Int_t anti=p->AntiParticle() ? 1:0;
812 if(kf<0) {
813 for(Int_t j=0;j<nparts;++j) {
814 TParticlePDG *dummy = dynamic_cast<TParticlePDG*>(fParticleList->At(j));
815 if(dummy==p->AntiParticle()) {
816 anti=j+1;
817 break;
818 }
819 }
820 fprintf(file,"%i 0\n",anti);
821 continue;
822 }
823
824 fprintf(file,"%i ",anti);
825 fprintf(file,"%i ",100);
826 fprintf(file,"%s ",p->ParticleClass());
827 fprintf(file,"% i ",(Int_t)p->Charge());
828 fprintf(file,"%.5le ",p->Mass());
829 fprintf(file,"%.5le ",p->Width());
830 fprintf(file,"%i ",(Int_t)p->Isospin());
831 fprintf(file,"%i ",(Int_t)p->I3());
832 fprintf(file,"%i ",(Int_t)p->Spin());
833 fprintf(file,"%i ",-1);
834 fprintf(file,"%i ",p->TrackingCode());
835 Int_t nch=p->NDecayChannels();
836 fprintf(file,"%i\n",nch);
837 if(nch==0) {
838 continue;
839 }
840 fprintf(file,"#----------------------------------------------------------------------\n");
841 fprintf(file,"# decay type(PY6) BR Nd daughters(codes, then names)\n");
842 fprintf(file,"#----------------------------------------------------------------------\n");
843 for(Int_t j=0;j<nch;++j) {
844 TDecayChannel *dc=p->DecayChannel(j);
845 if (!dc) continue;
846 fprintf(file,"%9i ",dc->Number()+1);
847 fprintf(file,"%3i ",dc->MatrixElementCode());
848 fprintf(file,"%.5le ",dc->BranchingRatio());
849 Int_t ndau=dc->NDaughters();
850 fprintf(file,"%3i ",ndau);
851 for (int idau=0; idau<ndau; idau++) {
852 fprintf(file,"%- 6i ",dc->DaughterPdgCode(idau));
853 }
854 for (int idau=0; idau<ndau; idau++) {
855 TParticlePDG *dummy=GetParticle(dc->DaughterPdgCode(idau));
856 if(dummy)
857 fprintf(file,"%-10s ",dummy->GetName());
858 else
859 fprintf(file,"%-10s ","???");
860 }
861 fprintf(file,"\n");
862 }
863 }
864 fclose(file);
865 return nparts;
866}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
long long Long64_t
Definition RtypesCore.h:69
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
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:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
#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.
Int_t MatrixElementCode()
Int_t DaughterPdgCode(Int_t i)
Double_t BranchingRatio()
Int_t NDaughters()
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:47
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 const char * GetName() const
Returns name of object.
Definition TObject.cxx:444
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:979
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition TObject.cxx:642
@ kInvalidObject
if object ctor succeeded but object should not be used
Definition TObject.h:72
Utility class used internally by TDatabasePDG.
void AddParticle(TObject *p)
Description of the static properties of a particle.
TDecayChannel * DecayChannel(Int_t i)
return pointer to decay channel object at index i
Int_t AddDecayChannel(Int_t Type, Double_t BranchingRatio, Int_t NDaughters, Int_t *DaughterPdgCode)
add new decay channel, Particle owns those...
TParticlePDG * AntiParticle()
Int_t NDecayChannels() const
void SetAntiParticle(TParticlePDG *ap)
static const TString & GetEtcDir()
Get the sysconfig directory in the installation. Static utility function.
Definition TROOT.cxx:3055
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
virtual const char * PrependPathName(const char *dir, TString &name)
Concatenate a directory and a file name.
Definition TSystem.cxx:1081
R__EXTERN TVirtualRWMutex * gCoreMutex
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250