Logo ROOT   6.14/05
Reference Guide
TGeoElement.cxx
Go to the documentation of this file.
1 // @(#)root/geom:$Id$
2 // Author: Andrei Gheata 17/06/04
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 /** \class TGeoElement
13 \ingroup Geometry_classes
14 Base class for chemical elements
15 */
16 
17 /** \class TGeoElementRN
18 \ingroup Geometry_classes
19 Class representing a radionuclide
20 */
21 
22 /** \class TGeoElemIter
23 \ingroup Geometry_classes
24 Iterator for decay branches
25 */
26 
27 /** \class TGeoDecayChannel
28 \ingroup Geometry_classes
29 A decay channel for a radionuclide
30 */
31 
32 /** \class TGeoElementTable
33 \ingroup Geometry_classes
34 Table of elements
35 */
36 
37 #include "RConfigure.h"
38 
39 #include "Riostream.h"
40 
41 #include "TSystem.h"
42 #include "TROOT.h"
43 #include "TObjArray.h"
44 #include "TVirtualGeoPainter.h"
45 #include "TGeoManager.h"
46 #include "TGeoElement.h"
47 #include "TMath.h"
48 #include "TGeoPhysicalConstants.h"
49 
50 // statics and globals
51 static const Int_t gMaxElem = 110;
52 static const Int_t gMaxLevel = 8;
53 static const Int_t gMaxDecay = 15;
54 
55 static const char gElName[gMaxElem][3] = {
56  "H ","He","Li","Be","B ","C ","N ","O ","F ","Ne","Na","Mg",
57  "Al","Si","P ","S ","Cl","Ar","K ","Ca","Sc","Ti","V ","Cr",
58  "Mn","Fe","Co","Ni","Cu","Zn","Ga","Ge","As","Se","Br","Kr",
59  "Rb","Sr","Y ","Zr","Nb","Mo","Tc","Ru","Rh","Pd","Ag","Cd",
60  "In","Sn","Sb","Te","I ","Xe","Cs","Ba","La","Ce","Pr","Nd",
61  "Pm","Sm","Eu","Gd","Tb","Dy","Ho","Er","Tm","Yb","Lu","Hf",
62  "Ta","W ","Re","Os","Ir","Pt","Au","Hg","Tl","Pb","Bi","Po",
63  "At","Rn","Fr","Ra","Ac","Th","Pa","U ","Np","Pu","Am","Cm",
64  "Bk","Cf","Es","Fm","Md","No","Lr","Rf","Db","Sg","Bh","Hs",
65  "Mt","Ds" };
66 
67 static const char *gDecayName[gMaxDecay+1] = {
68  "2BetaMinus", "BetaMinus", "NeutronEm", "ProtonEm", "Alpha", "ECF",
69  "ElecCapt", "IsoTrans", "I", "SpontFiss", "2ProtonEm", "2NeutronEm",
70  "2Alpha", "Carbon12", "Carbon14", "Stable" };
71 
72 static const Int_t gDecayDeltaA[gMaxDecay] = {
73  0, 0, -1, -1, -4,
74  -99, 0, 0, -99, -99,
75  -2, -2, -8, -12, -14 };
76 
77 static const Int_t gDecayDeltaZ[gMaxDecay] = {
78  2, 1, 0, -1, -2,
79  -99, -1, 0, -99, -99,
80  -2, 0, -4, -6, -6 };
81 static const char gLevName[gMaxLevel]=" mnpqrs";
82 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Default constructor
87 
89 {
91  SetUsed(kFALSE);
92  fZ = 0;
93  fN = 0;
94  fNisotopes = 0;
95  fA = 0.0;
96  fIsotopes = NULL;
97  fAbundances = NULL;
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Obsolete constructor
102 
103 TGeoElement::TGeoElement(const char *name, const char *title, Int_t z, Double_t a)
104  :TNamed(name, title)
105 {
107  SetUsed(kFALSE);
108  fZ = z;
109  fN = Int_t(a);
110  fNisotopes = 0;
111  fA = a;
112  fIsotopes = NULL;
113  fAbundances = NULL;
115 }
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// Element having isotopes.
119 
120 TGeoElement::TGeoElement(const char *name, const char *title, Int_t nisotopes)
121  :TNamed(name, title)
122 {
124  SetUsed(kFALSE);
125  fZ = 0;
126  fN = 0;
127  fNisotopes = nisotopes;
128  fA = 0.0;
129  fIsotopes = new TObjArray(nisotopes);
130  fAbundances = new Double_t[nisotopes];
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// Constructor
135 
136 TGeoElement::TGeoElement(const char *name, const char *title, Int_t z, Int_t n, Double_t a)
137  :TNamed(name, title)
138 {
140  SetUsed(kFALSE);
141  fZ = z;
142  fN = n;
143  fNisotopes = 0;
144  fA = a;
145  fIsotopes = NULL;
146  fAbundances = NULL;
148 }
149 ////////////////////////////////////////////////////////////////////////////////
150 /// Calculate properties for an atomic number
151 
153 {
154  // Radiation Length
157 }
158 ////////////////////////////////////////////////////////////////////////////////
159 /// Compute Coulomb correction factor (Phys Rev. D50 3-1 (1994) page 1254)
160 
162 {
163  static const Double_t k1 = 0.0083 , k2 = 0.20206 ,k3 = 0.0020 , k4 = 0.0369 ;
164 
166  Double_t az4 = az2 * az2;
167 
168  fCoulomb = (k1*az4 + k2 + 1./(1.+az2))*az2 - (k3*az4 + k4)*az4;
169 }
170 ////////////////////////////////////////////////////////////////////////////////
171 /// Compute Tsai's Expression for the Radiation Length (Phys Rev. D50 3-1 (1994) page 1254)
172 
174 {
175  static const Double_t Lrad_light[] = {5.31 , 4.79 , 4.74 , 4.71} ;
176  static const Double_t Lprad_light[] = {6.144 , 5.621 , 5.805 , 5.924} ;
177 
178  fRadTsai = 0.0;
179  if (fZ == 0) return;
180  const Double_t logZ3 = TMath::Log(fZ)/3.;
181 
182  Double_t Lrad, Lprad;
183  Int_t iz = static_cast<Int_t>(fZ+0.5) - 1 ; // The static cast comes from G4lrint
184  static const Double_t log184 = TMath::Log(184.15);
185  static const Double_t log1194 = TMath::Log(1194.);
186  if (iz <= 3) { Lrad = Lrad_light[iz] ; Lprad = Lprad_light[iz] ; }
187  else { Lrad = log184 - logZ3 ; Lprad = log1194 - 2*logZ3;}
188 
189  fRadTsai = 4*TGeoUnit::alpha_rcl2*fZ*(fZ*(Lrad-fCoulomb) + Lprad);
190 }
191 ////////////////////////////////////////////////////////////////////////////////
192 /// Print this isotope
193 
194 void TGeoElement::Print(Option_t *option) const
195 {
196  printf("Element: %s Z=%d N=%f A=%f [g/mole]\n", GetName(), fZ,Neff(),fA);
197  if (HasIsotopes()) {
198  for (Int_t i=0; i<fNisotopes; i++) {
199  TGeoIsotope *iso = GetIsotope(i);
200  printf("=>Isotope %s, abundance=%f :\n", iso->GetName(), fAbundances[i]);
201  iso->Print(option);
202  }
203  }
204 }
205 
206 ////////////////////////////////////////////////////////////////////////////////
207 /// Returns pointer to the table.
208 
210 {
211  if (!gGeoManager) {
212  ::Error("TGeoElementTable::GetElementTable", "Create a geometry manager first");
213  return NULL;
214  }
215  return gGeoManager->GetElementTable();
216 }
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 /// Add an isotope for this element. All isotopes have to be isotopes of the same element.
220 
221 void TGeoElement::AddIsotope(TGeoIsotope *isotope, Double_t relativeAbundance)
222 {
223  if (!fIsotopes) {
224  Fatal("AddIsotope", "Cannot add isotopes to normal elements. Use constructor with number of isotopes.");
225  return;
226  }
227  Int_t ncurrent = 0;
228  TGeoIsotope *isocrt;
229  for (ncurrent=0; ncurrent<fNisotopes; ncurrent++)
230  if (!fIsotopes->At(ncurrent)) break;
231  if (ncurrent==fNisotopes) {
232  Error("AddIsotope", "All %d isotopes of element %s already defined", fNisotopes, GetName());
233  return;
234  }
235  // Check Z of the new isotope
236  if ((fZ!=0) && (isotope->GetZ()!=fZ)) {
237  Fatal("AddIsotope", "Trying to add isotope %s with different Z to the same element %s",
238  isotope->GetName(), GetName());
239  return;
240  } else {
241  fZ = isotope->GetZ();
242  }
243  fIsotopes->Add(isotope);
244  fAbundances[ncurrent] = relativeAbundance;
245  if (ncurrent==fNisotopes-1) {
246  Double_t weight = 0.0;
247  Double_t aeff = 0.0;
248  Double_t neff = 0.0;
249  for (Int_t i=0; i<fNisotopes; i++) {
250  isocrt = (TGeoIsotope*)fIsotopes->At(i);
251  aeff += fAbundances[i]*isocrt->GetA();
252  neff += fAbundances[i]*isocrt->GetN();
253  weight += fAbundances[i];
254  }
255  aeff /= weight;
256  neff /= weight;
257  fN = (Int_t)neff;
258  fA = aeff;
259  }
261 }
262 
263 ////////////////////////////////////////////////////////////////////////////////
264 /// Returns effective number of nucleons.
265 
267 {
268  if (!fNisotopes) return fN;
269  TGeoIsotope *isocrt;
270  Double_t weight = 0.0;
271  Double_t neff = 0.0;
272  for (Int_t i=0; i<fNisotopes; i++) {
273  isocrt = (TGeoIsotope*)fIsotopes->At(i);
274  neff += fAbundances[i]*isocrt->GetN();
275  weight += fAbundances[i];
276  }
277  neff /= weight;
278  return neff;
279 }
280 
281 ////////////////////////////////////////////////////////////////////////////////
282 /// Return i-th isotope in the element.
283 
285 {
286  if (i>=0 && i<fNisotopes) {
287  return (TGeoIsotope*)fIsotopes->At(i);
288  }
289  return NULL;
290 }
291 
292 ////////////////////////////////////////////////////////////////////////////////
293 /// Return relative abundance of i-th isotope in this element.
294 
296 {
297  if (i>=0 && i<fNisotopes) return fAbundances[i];
298  return 0.0;
299 }
300 
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// Dummy I/O constructor
305 
307  :TNamed(),
308  fZ(0),
309  fN(0),
310  fA(0)
311 {
312 }
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Constructor
316 
318  :TNamed(name,""),
319  fZ(z),
320  fN(n),
321  fA(a)
322 {
323  if (z<1) Fatal("ctor", "Not allowed Z=%d (<1) for isotope: %s", z,name);
324  if (n<z) Fatal("ctor", "Not allowed Z=%d < N=%d for isotope: %s", z,n,name);
326 }
327 
328 ////////////////////////////////////////////////////////////////////////////////
329 /// Find existing isotope by name.
330 
332 {
334  if (!elTable) return 0;
335  return elTable->FindIsotope(name);
336 }
337 
338 ////////////////////////////////////////////////////////////////////////////////
339 /// Print this isotope
340 
342 {
343  printf("Isotope: %s Z=%d N=%d A=%f [g/mole]\n", GetName(), fZ,fN,fA);
344 }
345 
347 
348 ////////////////////////////////////////////////////////////////////////////////
349 /// Default constructor
350 
352 {
353  TObject::SetBit(kElementChecked,kFALSE);
354  fENDFcode = 0;
355  fIso = 0;
356  fLevel = 0;
357  fDeltaM = 0;
358  fHalfLife = 0;
359  fNatAbun = 0;
360  fTH_F = 0;
361  fTG_F = 0;
362  fTH_S = 0;
363  fTG_S = 0;
364  fStatus = 0;
365  fRatio = 0;
366  fDecays = 0;
367 }
368 
369 ////////////////////////////////////////////////////////////////////////////////
370 /// Constructor.
371 
373  Double_t deltaM, Double_t halfLife, const char* JP,
374  Double_t natAbun, Double_t th_f, Double_t tg_f, Double_t th_s,
375  Double_t tg_s, Int_t status)
376  :TGeoElement("", JP, zz, aa)
377 {
379  fENDFcode = ENDF(aa,zz,iso);
380  fIso = iso;
381  fLevel = level;
382  fDeltaM = deltaM;
383  fHalfLife = halfLife;
384  fTitle = JP;
385  if (!fTitle.Length()) fTitle = "?";
386  fNatAbun = natAbun;
387  fTH_F = th_f;
388  fTG_F = tg_f;
389  fTH_S = th_s;
390  fTG_S = tg_s;
391  fStatus = status;
392  fDecays = 0;
393  fRatio = 0;
394  MakeName(aa,zz,iso);
395  if ((TMath::Abs(fHalfLife)<1.e-30) || fHalfLife<-1) Warning("ctor","Element %s has T1/2=%g [s]", fName.Data(), fHalfLife);
396 }
397 
398 ////////////////////////////////////////////////////////////////////////////////
399 /// Destructor.
400 
402 {
403  if (fDecays) {
404  fDecays->Delete();
405  delete fDecays;
406  }
407  if (fRatio) delete fRatio;
408 }
409 
410 ////////////////////////////////////////////////////////////////////////////////
411 /// Adds a decay mode for this element.
412 
413 void TGeoElementRN::AddDecay(Int_t decay, Int_t diso, Double_t branchingRatio, Double_t qValue)
414 {
415  if (branchingRatio<1E-20) {
416  TString decayName;
417  TGeoDecayChannel::DecayName(decay, decayName);
418  Warning("AddDecay", "Decay %s of %s has BR=0. Not added.", decayName.Data(),fName.Data());
419  return;
420  }
421  TGeoDecayChannel *dc = new TGeoDecayChannel(decay,diso,branchingRatio, qValue);
422  dc->SetParent(this);
423  if (!fDecays) fDecays = new TObjArray(5);
424  fDecays->Add(dc);
425 }
426 
427 ////////////////////////////////////////////////////////////////////////////////
428 /// Adds a decay channel to the list of decays.
429 
431 {
432  dc->SetParent(this);
433  if (!fDecays) fDecays = new TObjArray(5);
434  fDecays->Add(dc);
435 }
436 
437 ////////////////////////////////////////////////////////////////////////////////
438 /// Get number of decay channels of this element.
439 
441 {
442  if (!fDecays) return 0;
443  return fDecays->GetEntriesFast();
444 }
445 
446 ////////////////////////////////////////////////////////////////////////////////
447 /// Get the activity in Bq of a gram of material made from this element.
448 
450 {
451  static const Double_t ln2 = TMath::Log(2.);
452  Double_t sa = (fHalfLife>0 && fA>0)?(ln2*TMath::Na()/fHalfLife/fA):0.;
453  return sa;
454 }
455 
456 ////////////////////////////////////////////////////////////////////////////////
457 /// Check if all decay chain of the element is OK.
458 
460 {
462  TObject *oelem = (TObject*)this;
463  TGeoDecayChannel *dc;
464  TGeoElementRN *elem;
466  TString decayName;
467  if (!table) {
468  Error("CheckDecays", "Element table not present");
469  return kFALSE;
470  }
471  Bool_t resultOK = kTRUE;
472  if (!fDecays) {
473  oelem->SetBit(kElementChecked,kTRUE);
474  return resultOK;
475  }
476  Double_t br = 0.;
477  Int_t decayResult = 0;
478  TIter next(fDecays);
479  while ((dc=(TGeoDecayChannel*)next())) {
480  br += dc->BranchingRatio();
481  decayResult = DecayResult(dc);
482  if (decayResult) {
483  elem = table->GetElementRN(decayResult);
484  if (!elem) {
485  TGeoDecayChannel::DecayName(dc->Decay(),decayName);
486  Error("CheckDecays", "Element after decay %s of %s not found in DB", decayName.Data(),fName.Data());
487  return kFALSE;
488  }
489  dc->SetDaughter(elem);
490  resultOK = elem->CheckDecays();
491  }
492  }
493  if (TMath::Abs(br-100) > 1.E-3) {
494  Warning("CheckDecays", "BR for decays of element %s sum-up = %f", fName.Data(), br);
495  resultOK = kFALSE;
496  }
497  oelem->SetBit(kElementChecked,kTRUE);
498  return resultOK;
499 }
500 
501 ////////////////////////////////////////////////////////////////////////////////
502 /// Returns ENDF code of decay result.
503 
505 {
506  Int_t da, dz, diso;
507  dc->DecayShift(da, dz, diso);
508  if (da == -99 || dz == -99) return 0;
509  return ENDF(Int_t(fA)+da,fZ+dz,fIso+diso);
510 }
511 
512 ////////////////////////////////////////////////////////////////////////////////
513 /// Fills the input array with the set of RN elements resulting from the decay of
514 /// this one. All element in the list will contain the time evolution of their
515 /// proportion by number with respect to this element. The proportion can be
516 /// retrieved via the method TGeoElementRN::Ratio().
517 /// The precision represent the minimum cumulative branching ratio for
518 /// which decay products are still taken into account.
519 
520 void TGeoElementRN::FillPopulation(TObjArray *population, Double_t precision, Double_t factor)
521 {
522  TGeoElementRN *elem;
523  TGeoElemIter next(this, precision);
524  TGeoBatemanSol s(this);
525  s.Normalize(factor);
526  AddRatio(s);
527  if (!population->FindObject(this)) population->Add(this);
528  while ((elem=next())) {
529  TGeoBatemanSol ratio(next.GetBranch());
530  ratio.Normalize(factor);
531  elem->AddRatio(ratio);
532  if (!population->FindObject(elem)) population->Add(elem);
533  }
534 }
535 
536 ////////////////////////////////////////////////////////////////////////////////
537 /// Generate a default name for the element.
538 
540 {
541  fName = "";
542  if (z==0 && a==1) {
543  fName = "neutron";
544  return;
545  }
546  if (z>=1 && z<= gMaxElem) fName += TString::Format("%3d-%s-",z,gElName[z-1]);
547  else fName = "?? -?? -";
548  if (a>=1 && a<=999) fName += TString::Format("%3.3d",a);
549  else fName += "??";
550  if (iso>0 && iso<gMaxLevel) fName += TString::Format("%c", gLevName[iso]);
551  fName.ReplaceAll(" ","");
552 }
553 
554 ////////////////////////////////////////////////////////////////////////////////
555 /// Print info about the element;
556 
557 void TGeoElementRN::Print(Option_t *option) const
558 {
559  printf("\n%-12s ",fName.Data());
560  printf("ENDF=%d; ",fENDFcode);
561  printf("A=%d; ",(Int_t)fA);
562  printf("Z=%d; ",fZ);
563  printf("Iso=%d; ",fIso);
564  printf("Level=%g[MeV]; ",fLevel);
565  printf("Dmass=%g[MeV]; ",fDeltaM);
566  if (fHalfLife>0) printf("Hlife=%g[s]\n",fHalfLife);
567  else printf("Hlife=INF\n");
568  printf("%13s"," ");
569  printf("J/P=%s; ",fTitle.Data());
570  printf("Abund=%g; ",fNatAbun);
571  printf("Htox=%g; ",fTH_F);
572  printf("Itox=%g; ",fTG_F);
573  printf("Stat=%d\n",fStatus);
574  if(!fDecays) return;
575  printf("Decay modes:\n");
576  TIter next(fDecays);
577  TGeoDecayChannel *dc;
578  while ((dc=(TGeoDecayChannel*)next())) dc->Print(option);
579 }
580 
581 ////////////////////////////////////////////////////////////////////////////////
582 /// Create element from line record.
583 
585 {
586  Int_t a,z,iso,status;
587  Double_t level, deltaM, halfLife, natAbun, th_f, tg_f, th_s, tg_s;
588  char name[20],jp[20];
589  sscanf(&line[0], "%s%d%d%d%lg%lg%lg%s%lg%lg%lg%lg%lg%d%d", name,&a,&z,&iso,&level,&deltaM,
590  &halfLife,jp,&natAbun,&th_f,&tg_f,&th_s,&tg_s,&status,&ndecays);
591  TGeoElementRN *elem = new TGeoElementRN(a,z,iso,level,deltaM,halfLife,
592  jp,natAbun,th_f,tg_f,th_s,tg_s,status);
593  return elem;
594 }
595 
596 ////////////////////////////////////////////////////////////////////////////////
597 /// Save primitive for RN elements.
598 
599 void TGeoElementRN::SavePrimitive(std::ostream &out, Option_t *option)
600 {
601  if (!strcmp(option,"h")) {
602  // print a header if requested
603  out << "#====================================================================================================================================" << std::endl;
604  out << "# Name A Z ISO LEV[MeV] DM[MeV] T1/2[s] J/P ABUND[%] HTOX ITOX HTOX ITOX STAT NDCY" << std::endl;
605  out << "#====================================================================================================================================" << std::endl;
606  }
607  out << std::setw(11) << fName.Data();
608  out << std::setw(5) << (Int_t)fA;
609  out << std::setw(5) << fZ;
610  out << std::setw(5) << fIso;
611  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fLevel;
612  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fDeltaM;
613  out << std::setw(10) << std::setiosflags(std::ios::scientific) << std::setprecision(3) << fHalfLife;
614  out << std::setw(13) << fTitle.Data();
615  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fNatAbun;
616  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fTH_F;
617  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fTG_F;
618  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fTH_S;
619  out << std::setw(10) << std::setiosflags(std::ios::fixed) << std::setprecision(5) << fTG_S;
620  out << std::setw(5) << fStatus;
621  Int_t ndecays = 0;
622  if (fDecays) ndecays = fDecays->GetEntries();
623  out << std::setw(5) << ndecays;
624  out << std::endl;
625  if (fDecays) {
626  TIter next(fDecays);
627  TGeoDecayChannel *dc;
628  while ((dc=(TGeoDecayChannel*)next())) dc->SavePrimitive(out);
629  }
630 }
631 
632 ////////////////////////////////////////////////////////////////////////////////
633 /// Adds a proportion ratio to the existing one.
634 
636 {
637  if (!fRatio) fRatio = new TGeoBatemanSol(ratio);
638  else *fRatio += ratio;
639 }
640 
641 ////////////////////////////////////////////////////////////////////////////////
642 /// Clears the existing ratio.
643 
645 {
646  if (fRatio) {
647  delete fRatio;
648  fRatio = 0;
649  }
650 }
651 
653 
654 ////////////////////////////////////////////////////////////////////////////////
655 /// Assignment.
656 ///assignment operator
657 
659 {
660  if(this!=&dc) {
661  TObject::operator=(dc);
662  fDecay = dc.fDecay;
663  fDiso = dc.fDiso;
664  fBranchingRatio = dc.fBranchingRatio;
665  fQvalue = dc.fQvalue;
666  fParent = dc.fParent;
667  fDaughter = dc.fDaughter;
668  }
669  return *this;
670 }
671 
672 ////////////////////////////////////////////////////////////////////////////////
673 /// Returns name of decay.
674 
675 const char *TGeoDecayChannel::GetName() const
676 {
677  static TString name = "";
678  name = "";
679  if (!fDecay) return gDecayName[gMaxDecay];
680  for (Int_t i=0; i<gMaxDecay; i++) {
681  if (1<<i & fDecay) {
682  if (name.Length()) name += "+";
683  name += gDecayName[i];
684  }
685  }
686  return name.Data();
687 }
688 
689 ////////////////////////////////////////////////////////////////////////////////
690 /// Returns decay name.
691 
693 {
694  if (!decay) {
695  name = gDecayName[gMaxDecay];
696  return;
697  }
698  name = "";
699  for (Int_t i=0; i<gMaxDecay; i++) {
700  if (1<<i & decay) {
701  if (name.Length()) name += "+";
702  name += gDecayName[i];
703  }
704  }
705 }
706 
707 ////////////////////////////////////////////////////////////////////////////////
708 /// Get index of this channel in the list of decays of the parent nuclide.
709 
711 {
712  return fParent->Decays()->IndexOf(this);
713 }
714 
715 ////////////////////////////////////////////////////////////////////////////////
716 /// Prints decay info.
717 
719 {
720  TString name;
721  DecayName(fDecay, name);
722  printf("%-20s Diso: %3d BR: %9.3f%% Qval: %g\n", name.Data(),fDiso,fBranchingRatio,fQvalue);
723 }
724 
725 ////////////////////////////////////////////////////////////////////////////////
726 /// Create element from line record.
727 
729 {
730  char name[80];
731  Int_t decay,diso;
732  Double_t branchingRatio, qValue;
733  sscanf(&line[0], "%s%d%d%lg%lg", name,&decay,&diso,&branchingRatio,&qValue);
734  TGeoDecayChannel *dc = new TGeoDecayChannel(decay,diso,branchingRatio,qValue);
735  return dc;
736 }
737 
738 ////////////////////////////////////////////////////////////////////////////////
739 /// Save primitive for decays.
740 
741 void TGeoDecayChannel::SavePrimitive(std::ostream &out, Option_t *)
742 {
743  TString decayName;
744  DecayName(fDecay, decayName);
745  out << std::setw(50) << decayName.Data();
746  out << std::setw(10) << fDecay;
747  out << std::setw(10) << fDiso;
748  out << std::setw(12) << std::setiosflags(std::ios::fixed) << std::setprecision(6) << fBranchingRatio;
749  out << std::setw(12) << std::setiosflags(std::ios::fixed) << std::setprecision(6) << fQvalue;
750  out << std::endl;
751 }
752 
753 ////////////////////////////////////////////////////////////////////////////////
754 /// Returns variation in A, Z and Iso after decay.
755 
757 {
758  dA=dZ=0;
759  dI=fDiso;
760  for(Int_t i=0; i<gMaxDecay; ++i) {
761  if(1<<i & fDecay) {
762  if(gDecayDeltaA[i] == -99 || gDecayDeltaZ[i] == -99 ) {
763  dA=dZ=-99;
764  return;
765  }
766  dA += gDecayDeltaA[i];
767  dZ += gDecayDeltaZ[i];
768  }
769  }
770 }
771 
773 
774 ////////////////////////////////////////////////////////////////////////////////
775 /// Default constructor.
776 
778  : fTop(top), fElem(top), fBranch(0), fLevel(0), fLimitRatio(limit), fRatio(1.)
779 {
780  fBranch = new TObjArray(10);
781 }
782 
783 ////////////////////////////////////////////////////////////////////////////////
784 /// Copy ctor.
785 
787  :fTop(iter.fTop),
788  fElem(iter.fElem),
789  fBranch(0),
790  fLevel(iter.fLevel),
791  fLimitRatio(iter.fLimitRatio),
792  fRatio(iter.fRatio)
793 {
794  if (iter.fBranch) {
795  fBranch = new TObjArray(10);
796  for (Int_t i=0; i<fLevel; i++) fBranch->Add(iter.fBranch->At(i));
797  }
798 }
799 
800 ////////////////////////////////////////////////////////////////////////////////
801 /// Destructor.
802 
804 {
805  if (fBranch) delete fBranch;
806 }
807 
808 ////////////////////////////////////////////////////////////////////////////////
809 /// Assignment.
810 
812 {
813  if (&iter == this) return *this;
814  fTop = iter.fTop;
815  fElem = iter.fElem;
816  fLevel = iter.fLevel;
817  if (iter.fBranch) {
818  fBranch = new TObjArray(10);
819  for (Int_t i=0; i<fLevel; i++) fBranch->Add(iter.fBranch->At(i));
820  }
821  fLimitRatio = iter.fLimitRatio;
822  fRatio = iter.fRatio;
823  return *this;
824 }
825 
826 ////////////////////////////////////////////////////////////////////////////////
827 /// () operator.
828 
830 {
831  return Next();
832 }
833 
834 ////////////////////////////////////////////////////////////////////////////////
835 /// Go upwards from the current location until the next branching, then down.
836 
838 {
839  TGeoDecayChannel *dc;
840  Int_t ind, nd;
841  while (fLevel) {
842  // Current decay channel
843  dc = (TGeoDecayChannel*)fBranch->At(fLevel-1);
844  ind = dc->GetIndex();
845  nd = dc->Parent()->GetNdecays();
846  fRatio /= 0.01*dc->BranchingRatio();
847  fElem = dc->Parent();
848  fBranch->RemoveAt(--fLevel);
849  ind++;
850  while (ind<nd) {
851  if (Down(ind++)) return (TGeoElementRN*)fElem;
852  }
853  }
854  fElem = NULL;
855  return NULL;
856 }
857 
858 ////////////////////////////////////////////////////////////////////////////////
859 /// Go downwards from current level via ibranch as low in the tree as possible.
860 /// Return value flags if the operation was successful.
861 
863 {
864  TGeoDecayChannel *dc = (TGeoDecayChannel*)fElem->Decays()->At(ibranch);
865  if (!dc->Daughter()) return NULL;
866  Double_t br = 0.01*fRatio*dc->BranchingRatio();
867  if (br < fLimitRatio) return NULL;
868  fLevel++;
869  fRatio = br;
870  fBranch->Add(dc);
871  fElem = dc->Daughter();
872  return (TGeoElementRN*)fElem;
873 }
874 
875 ////////////////////////////////////////////////////////////////////////////////
876 /// Return next element.
877 
879 {
880  if (!fElem) return NULL;
881  // Check if this is the first iteration.
882  Int_t nd = fElem->GetNdecays();
883  for (Int_t i=0; i<nd; i++) if (Down(i)) return (TGeoElementRN*)fElem;
884  return Up();
885 }
886 
887 ////////////////////////////////////////////////////////////////////////////////
888 /// Print info about the current decay branch.
889 
890 void TGeoElemIter::Print(Option_t * /*option*/) const
891 {
892  TGeoElementRN *elem;
893  TGeoDecayChannel *dc;
894  TString indent = "";
895  printf("=== Chain with %g %%\n", 100*fRatio);
896  for (Int_t i=0; i<fLevel; i++) {
897  dc = (TGeoDecayChannel*)fBranch->At(i);
898  elem = dc->Parent();
899  printf("%s%s (%g%% %s) T1/2=%g\n", indent.Data(), elem->GetName(),dc->BranchingRatio(),dc->GetName(),elem->HalfLife());
900  indent += " ";
901  if (i==fLevel-1) {
902  elem = dc->Daughter();
903  printf("%s%s\n", indent.Data(), elem->GetName());
904  }
905  }
906 }
907 
909 
910 ////////////////////////////////////////////////////////////////////////////////
911 /// default constructor
912 
914 {
915  fNelements = 0;
916  fNelementsRN = 0;
917  fNisotopes = 0;
918  fList = 0;
919  fListRN = 0;
920  fIsotopes = 0;
921 }
922 
923 ////////////////////////////////////////////////////////////////////////////////
924 /// constructor
925 
927 {
928  fNelements = 0;
929  fNelementsRN = 0;
930  fNisotopes = 0;
931  fList = new TObjArray(128);
932  fListRN = 0;
933  fIsotopes = 0;
934  BuildDefaultElements();
935 // BuildElementsRN();
936 }
937 
938 ////////////////////////////////////////////////////////////////////////////////
939 ///copy constructor
940 
942  TObject(get),
943  fNelements(get.fNelements),
944  fNelementsRN(get.fNelementsRN),
945  fNisotopes(get.fNisotopes),
946  fList(get.fList),
947  fListRN(get.fListRN),
948  fIsotopes(0)
949 {
950 }
951 
952 ////////////////////////////////////////////////////////////////////////////////
953 ///assignment operator
954 
956 {
957  if(this!=&get) {
958  TObject::operator=(get);
959  fNelements=get.fNelements;
960  fNelementsRN=get.fNelementsRN;
961  fNisotopes=get.fNisotopes;
962  fList=get.fList;
963  fListRN=get.fListRN;
964  fIsotopes = 0;
965  }
966  return *this;
967 }
968 
969 ////////////////////////////////////////////////////////////////////////////////
970 /// destructor
971 
973 {
974  if (fList) {
975  fList->Delete();
976  delete fList;
977  }
978  if (fListRN) {
979  fListRN->Delete();
980  delete fListRN;
981  }
982  if (fIsotopes) {
983  fIsotopes->Delete();
984  delete fIsotopes;
985  }
986 }
987 
988 ////////////////////////////////////////////////////////////////////////////////
989 /// Add an element to the table. Obsolete.
990 
991 void TGeoElementTable::AddElement(const char *name, const char *title, Int_t z, Double_t a)
992 {
993  if (!fList) fList = new TObjArray(128);
994  fList->AddAtAndExpand(new TGeoElement(name,title,z,a), fNelements++);
995 }
996 
997 ////////////////////////////////////////////////////////////////////////////////
998 /// Add an element to the table.
999 
1000 void TGeoElementTable::AddElement(const char *name, const char *title, Int_t z, Int_t n, Double_t a)
1001 {
1002  if (!fList) fList = new TObjArray(128);
1003  fList->AddAtAndExpand(new TGeoElement(name,title,z,n,a), fNelements++);
1004 }
1005 
1006 ////////////////////////////////////////////////////////////////////////////////
1007 /// Add a custom element to the table.
1008 
1010 {
1011  if (!fList) fList = new TObjArray(128);
1012  TGeoElement *orig = FindElement(elem->GetName());
1013  if (orig) {
1014  Error("AddElement", "Found element with same name: %s (%s). Cannot add to table.",
1015  orig->GetName(), orig->GetTitle());
1016  return;
1017  }
1018  fList->AddAtAndExpand(elem, fNelements++);
1019 }
1020 
1021 ////////////////////////////////////////////////////////////////////////////////
1022 /// Add a radionuclide to the table and map it.
1023 
1025 {
1026  if (!fListRN) fListRN = new TObjArray(3600);
1027  if (HasRNElements() && GetElementRN(elem->ENDFCode())) return;
1028 // elem->Print();
1029  fListRN->Add(elem);
1030  fNelementsRN++;
1031  fElementsRN.insert(ElementRNMap_t::value_type(elem->ENDFCode(), elem));
1032 }
1033 
1034 ////////////////////////////////////////////////////////////////////////////////
1035 /// Add isotope to the table.
1036 
1038 {
1039  if (FindIsotope(isotope->GetName())) {
1040  Error("AddIsotope", "Isotope with the same name: %s already in table. Not adding.",isotope->GetName());
1041  return;
1042  }
1043  if (!fIsotopes) fIsotopes = new TObjArray();
1044  fIsotopes->Add(isotope);
1045 }
1046 
1047 ////////////////////////////////////////////////////////////////////////////////
1048 /// Creates the default element table
1049 
1051 {
1052  if (HasDefaultElements()) return;
1053  AddElement("VACUUM","VACUUM" ,0, 0, 0.0);
1054  AddElement("H" ,"HYDROGEN" ,1, 1, 1.00794);
1055  AddElement("HE" ,"HELIUM" ,2, 4, 4.002602);
1056  AddElement("LI" ,"LITHIUM" ,3, 7, 6.941);
1057  AddElement("BE" ,"BERYLLIUM" ,4, 9, 9.01218);
1058  AddElement("B" ,"BORON" ,5, 11, 10.811);
1059  AddElement("C" ,"CARBON" ,6, 12, 12.0107);
1060  AddElement("N" ,"NITROGEN" ,7, 14, 14.00674);
1061  AddElement("O" ,"OXYGEN" ,8, 16, 15.9994);
1062  AddElement("F" ,"FLUORINE" ,9, 19, 18.9984032);
1063  AddElement("NE" ,"NEON" ,10, 20, 20.1797);
1064  AddElement("NA" ,"SODIUM" ,11, 23, 22.989770);
1065  AddElement("MG" ,"MAGNESIUM" ,12, 24, 24.3050);
1066  AddElement("AL" ,"ALUMINIUM" ,13, 27, 26.981538);
1067  AddElement("SI" ,"SILICON" ,14, 28, 28.0855);
1068  AddElement("P" ,"PHOSPHORUS" ,15, 31, 30.973761);
1069  AddElement("S" ,"SULFUR" ,16, 32, 32.066);
1070  AddElement("CL" ,"CHLORINE" ,17, 35, 35.4527);
1071  AddElement("AR" ,"ARGON" ,18, 40, 39.948);
1072  AddElement("K" ,"POTASSIUM" ,19, 39, 39.0983);
1073  AddElement("CA" ,"CALCIUM" ,20, 40, 40.078);
1074  AddElement("SC" ,"SCANDIUM" ,21, 45, 44.955910);
1075  AddElement("TI" ,"TITANIUM" ,22, 48, 47.867);
1076  AddElement("V" ,"VANADIUM" ,23, 51, 50.9415);
1077  AddElement("CR" ,"CHROMIUM" ,24, 52, 51.9961);
1078  AddElement("MN" ,"MANGANESE" ,25, 55, 54.938049);
1079  AddElement("FE" ,"IRON" ,26, 56, 55.845);
1080  AddElement("CO" ,"COBALT" ,27, 59, 58.933200);
1081  AddElement("NI" ,"NICKEL" ,28, 59, 58.6934);
1082  AddElement("CU" ,"COPPER" ,29, 64, 63.546);
1083  AddElement("ZN" ,"ZINC" ,30, 65, 65.39);
1084  AddElement("GA" ,"GALLIUM" ,31, 70, 69.723);
1085  AddElement("GE" ,"GERMANIUM" ,32, 73, 72.61);
1086  AddElement("AS" ,"ARSENIC" ,33, 75, 74.92160);
1087  AddElement("SE" ,"SELENIUM" ,34, 79, 78.96);
1088  AddElement("BR" ,"BROMINE" ,35, 80, 79.904);
1089  AddElement("KR" ,"KRYPTON" ,36, 84, 83.80);
1090  AddElement("RB" ,"RUBIDIUM" ,37, 85, 85.4678);
1091  AddElement("SR" ,"STRONTIUM" ,38, 88, 87.62);
1092  AddElement("Y" ,"YTTRIUM" ,39, 89, 88.90585);
1093  AddElement("ZR" ,"ZIRCONIUM" ,40, 91, 91.224);
1094  AddElement("NB" ,"NIOBIUM" ,41, 93, 92.90638);
1095  AddElement("MO" ,"MOLYBDENUM" ,42, 96, 95.94);
1096  AddElement("TC" ,"TECHNETIUM" ,43, 98, 98.0);
1097  AddElement("RU" ,"RUTHENIUM" ,44, 101, 101.07);
1098  AddElement("RH" ,"RHODIUM" ,45, 103, 102.90550);
1099  AddElement("PD" ,"PALLADIUM" ,46, 106, 106.42);
1100  AddElement("AG" ,"SILVER" ,47, 108, 107.8682);
1101  AddElement("CD" ,"CADMIUM" ,48, 112, 112.411);
1102  AddElement("IN" ,"INDIUM" ,49, 115, 114.818);
1103  AddElement("SN" ,"TIN" ,50, 119, 118.710);
1104  AddElement("SB" ,"ANTIMONY" ,51, 122, 121.760);
1105  AddElement("TE" ,"TELLURIUM" ,52, 128, 127.60);
1106  AddElement("I" ,"IODINE" ,53, 127, 126.90447);
1107  AddElement("XE" ,"XENON" ,54, 131, 131.29);
1108  AddElement("CS" ,"CESIUM" ,55, 133, 132.90545);
1109  AddElement("BA" ,"BARIUM" ,56, 137, 137.327);
1110  AddElement("LA" ,"LANTHANUM" ,57, 139, 138.9055);
1111  AddElement("CE" ,"CERIUM" ,58, 140, 140.116);
1112  AddElement("PR" ,"PRASEODYMIUM" ,59, 141, 140.90765);
1113  AddElement("ND" ,"NEODYMIUM" ,60, 144, 144.24);
1114  AddElement("PM" ,"PROMETHIUM" ,61, 145, 145.0);
1115  AddElement("SM" ,"SAMARIUM" ,62, 150, 150.36);
1116  AddElement("EU" ,"EUROPIUM" ,63, 152, 151.964);
1117  AddElement("GD" ,"GADOLINIUM" ,64, 157, 157.25);
1118  AddElement("TB" ,"TERBIUM" ,65, 159, 158.92534);
1119  AddElement("DY" ,"DYSPROSIUM" ,66, 162, 162.50);
1120  AddElement("HO" ,"HOLMIUM" ,67, 165, 164.93032);
1121  AddElement("ER" ,"ERBIUM" ,68, 167, 167.26);
1122  AddElement("TM" ,"THULIUM" ,69, 169, 168.93421);
1123  AddElement("YB" ,"YTTERBIUM" ,70, 173, 173.04);
1124  AddElement("LU" ,"LUTETIUM" ,71, 175, 174.967);
1125  AddElement("HF" ,"HAFNIUM" ,72, 178, 178.49);
1126  AddElement("TA" ,"TANTALUM" ,73, 181, 180.9479);
1127  AddElement("W" ,"TUNGSTEN" ,74, 184, 183.84);
1128  AddElement("RE" ,"RHENIUM" ,75, 186, 186.207);
1129  AddElement("OS" ,"OSMIUM" ,76, 190, 190.23);
1130  AddElement("IR" ,"IRIDIUM" ,77, 192, 192.217);
1131  AddElement("PT" ,"PLATINUM" ,78, 195, 195.078);
1132  AddElement("AU" ,"GOLD" ,79, 197, 196.96655);
1133  AddElement("HG" ,"MERCURY" ,80, 200, 200.59);
1134  AddElement("TL" ,"THALLIUM" ,81, 204, 204.3833);
1135  AddElement("PB" ,"LEAD" ,82, 207, 207.2);
1136  AddElement("BI" ,"BISMUTH" ,83, 209, 208.98038);
1137  AddElement("PO" ,"POLONIUM" ,84, 209, 209.0);
1138  AddElement("AT" ,"ASTATINE" ,85, 210, 210.0);
1139  AddElement("RN" ,"RADON" ,86, 222, 222.0);
1140  AddElement("FR" ,"FRANCIUM" ,87, 223, 223.0);
1141  AddElement("RA" ,"RADIUM" ,88, 226, 226.0);
1142  AddElement("AC" ,"ACTINIUM" ,89, 227, 227.0);
1143  AddElement("TH" ,"THORIUM" ,90, 232, 232.0381);
1144  AddElement("PA" ,"PROTACTINIUM" ,91, 231, 231.03588);
1145  AddElement("U" ,"URANIUM" ,92, 238, 238.0289);
1146  AddElement("NP" ,"NEPTUNIUM" ,93, 237, 237.0);
1147  AddElement("PU" ,"PLUTONIUM" ,94, 244, 244.0);
1148  AddElement("AM" ,"AMERICIUM" ,95, 243, 243.0);
1149  AddElement("CM" ,"CURIUM" ,96, 247, 247.0);
1150  AddElement("BK" ,"BERKELIUM" ,97, 247, 247.0);
1151  AddElement("CF" ,"CALIFORNIUM",98, 251, 251.0);
1152  AddElement("ES" ,"EINSTEINIUM",99, 252, 252.0);
1153  AddElement("FM" ,"FERMIUM" ,100, 257, 257.0);
1154  AddElement("MD" ,"MENDELEVIUM",101, 258, 258.0);
1155  AddElement("NO" ,"NOBELIUM" ,102, 259, 259.0);
1156  AddElement("LR" ,"LAWRENCIUM" ,103, 262, 262.0);
1157  AddElement("RF" ,"RUTHERFORDIUM",104, 261, 261.0);
1158  AddElement("DB" ,"DUBNIUM" ,105, 262, 262.0);
1159  AddElement("SG" ,"SEABORGIUM" ,106, 263, 263.0);
1160  AddElement("BH" ,"BOHRIUM" ,107, 262, 262.0);
1161  AddElement("HS" ,"HASSIUM" ,108, 265, 265.0);
1162  AddElement("MT" ,"MEITNERIUM" ,109, 266, 266.0);
1163  AddElement("UUN" ,"UNUNNILIUM" ,110, 269, 269.0);
1164  AddElement("UUU" ,"UNUNUNIUM" ,111, 272, 272.0);
1165  AddElement("UUB" ,"UNUNBIUM" ,112, 277, 277.0);
1166 
1168 }
1169 
1170 ////////////////////////////////////////////////////////////////////////////////
1171 /// Creates the list of radionuclides.
1172 
1174 {
1175  if (HasRNElements()) return;
1176  TGeoElementRN *elem;
1177  TString rnf = "RadioNuclides.txt";
1179  FILE *fp = fopen(rnf, "r");
1180  if (!fp) {
1181  Error("ImportElementsRN","File RadioNuclides.txt not found");
1182  return;
1183  }
1184  char line[150];
1185  Int_t ndecays = 0;
1186  Int_t i;
1187  while (fgets(&line[0],140,fp)) {
1188  if (line[0]=='#') continue;
1189  elem = TGeoElementRN::ReadElementRN(line, ndecays);
1190  for (i=0; i<ndecays; i++) {
1191  if (!fgets(&line[0],140,fp)) {
1192  Error("ImportElementsRN", "Error parsing RadioNuclides.txt file");
1193  fclose(fp);
1194  return;
1195  }
1197  elem->AddDecay(dc);
1198  }
1199  AddElementRN(elem);
1200 // elem->Print();
1201  }
1203  CheckTable();
1204  fclose(fp);
1205 }
1206 
1207 ////////////////////////////////////////////////////////////////////////////////
1208 /// Checks status of element table.
1209 
1211 {
1212  if (!HasRNElements()) return HasDefaultElements();
1213  TGeoElementRN *elem;
1214  Bool_t result = kTRUE;
1215  TIter next(fListRN);
1216  while ((elem=(TGeoElementRN*)next())) {
1217  if (!elem->CheckDecays()) result = kFALSE;
1218  }
1219  return result;
1220 }
1221 
1222 ////////////////////////////////////////////////////////////////////////////////
1223 /// Export radionuclides in a file.
1224 
1225 void TGeoElementTable::ExportElementsRN(const char *filename)
1226 {
1227  if (!HasRNElements()) return;
1228  TString sname = filename;
1229  if (!sname.Length()) sname = "RadioNuclides.txt";
1230  std::ofstream out;
1231  out.open(sname.Data(), std::ios::out);
1232  if (!out.good()) {
1233  Error("ExportElementsRN", "Cannot open file %s", sname.Data());
1234  return;
1235  }
1236 
1237  TGeoElementRN *elem;
1238  TIter next(fListRN);
1239  Int_t i=0;
1240  while ((elem=(TGeoElementRN*)next())) {
1241  if ((i%48)==0) elem->SavePrimitive(out,"h");
1242  else elem->SavePrimitive(out);
1243  i++;
1244  }
1245  out.close();
1246 }
1247 
1248 ////////////////////////////////////////////////////////////////////////////////
1249 /// Search an element by symbol or full name
1250 /// Exact matching
1251 
1253 {
1254  TGeoElement *elem;
1255  elem = (TGeoElement*)fList->FindObject(name);
1256  if (elem) return elem;
1257  // Search case insensitive by element name
1258  TString s(name);
1259  s.ToUpper();
1260  elem = (TGeoElement*)fList->FindObject(s.Data());
1261  if (elem) return elem;
1262  // Search by full name
1263  TIter next(fList);
1264  while ((elem=(TGeoElement*)next())) {
1265  if (s == elem->GetTitle()) return elem;
1266  }
1267  return 0;
1268 }
1269 
1270 ////////////////////////////////////////////////////////////////////////////////
1271 /// Find existing isotope by name. Not optimized for a big number of isotopes.
1272 
1274 {
1275  if (!fIsotopes) return NULL;
1276  return (TGeoIsotope*)fIsotopes->FindObject(name);
1277 }
1278 
1279 ////////////////////////////////////////////////////////////////////////////////
1280 /// Retrieve a radionuclide by ENDF code.
1281 
1283 {
1284  if (!HasRNElements()) {
1285  TGeoElementTable *table = (TGeoElementTable*)this;
1286  table->ImportElementsRN();
1287  if (!fListRN) return 0;
1288  }
1289  ElementRNMap_t::const_iterator it = fElementsRN.find(ENDFcode);
1290  if (it != fElementsRN.end()) return it->second;
1291  return 0;
1292 }
1293 
1294 ////////////////////////////////////////////////////////////////////////////////
1295 /// Retrieve a radionuclide by a, z, and isomeric state.
1296 
1298 {
1299  return GetElementRN(TGeoElementRN::ENDF(a,z,iso));
1300 }
1301 
1302 ////////////////////////////////////////////////////////////////////////////////
1303 /// Print table of elements. The accepted options are:
1304 /// "" - prints everything by default
1305 /// "D" - prints default elements only
1306 /// "I" - prints isotopes
1307 /// "R" - prints radio-nuclides only if imported
1308 /// "U" - prints user-defined elements only
1309 
1311 {
1312  TString opt(option);
1313  opt.ToUpper();
1314  Int_t induser = HasDefaultElements() ? 113 : 0;
1315  // Default elements
1316  if (opt=="" || opt=="D") {
1317  if (induser) printf("================\nDefault elements\n================\n");
1318  for (Int_t iel=0; iel<induser; ++iel) fList->At(iel)->Print();
1319  }
1320  // Isotopes
1321  if (opt=="" || opt=="I") {
1322  if (fIsotopes) {
1323  printf("================\nIsotopes\n================\n");
1324  fIsotopes->Print();
1325  }
1326  }
1327  // Radio-nuclides
1328  if (opt=="" || opt=="R") {
1329  if (HasRNElements()) {
1330  printf("================\nRadio-nuclides\n================\n");
1331  fListRN->Print();
1332  }
1333  }
1334  // User-defined elements
1335  if (opt=="" || opt=="U") {
1336  if (fNelements>induser) printf("================\nUser elements\n================\n");
1337  for (Int_t iel=induser; iel<fNelements; ++iel) fList->At(iel)->Print();
1338  }
1339 }
1340 
1342 
1343 ////////////////////////////////////////////////////////////////////////////////
1344 /// Default ctor.
1345 
1347  :TObject(), TAttLine(), TAttFill(), TAttMarker(),
1348  fElem(elem),
1349  fElemTop(elem),
1350  fCsize(10),
1351  fNcoeff(0),
1352  fFactor(1.),
1353  fTmin(0.),
1354  fTmax(0.),
1355  fCoeff(NULL)
1356 {
1357  fCoeff = new BtCoef_t[fCsize];
1358  fNcoeff = 1;
1359  fCoeff[0].cn = 1.;
1360  Double_t t12 = elem->HalfLife();
1361  if (t12 == 0.) t12 = 1.e-30;
1362  if (elem->Stable()) fCoeff[0].lambda = 0.;
1363  else fCoeff[0].lambda = TMath::Log(2.)/t12;
1364 }
1365 
1366 ////////////////////////////////////////////////////////////////////////////////
1367 /// Default ctor.
1368 
1370  :TObject(), TAttLine(), TAttFill(), TAttMarker(),
1371  fElem(NULL),
1372  fElemTop(NULL),
1373  fCsize(0),
1374  fNcoeff(0),
1375  fFactor(1.),
1376  fTmin(0.),
1377  fTmax(0.),
1378  fCoeff(NULL)
1379 {
1380  TGeoDecayChannel *dc = (TGeoDecayChannel*)chain->At(0);
1381  if (dc) fElemTop = dc->Parent();
1382  dc = (TGeoDecayChannel*)chain->At(chain->GetEntriesFast()-1);
1383  if (dc) {
1384  fElem = dc->Daughter();
1385  fCsize = chain->GetEntriesFast()+1;
1386  fCoeff = new BtCoef_t[fCsize];
1387  FindSolution(chain);
1388  }
1389 }
1390 
1391 ////////////////////////////////////////////////////////////////////////////////
1392 /// Copy constructor.
1393 
1395  :TObject(other), TAttLine(other), TAttFill(other), TAttMarker(other),
1396  fElem(other.fElem),
1397  fElemTop(other.fElemTop),
1398  fCsize(other.fCsize),
1399  fNcoeff(other.fNcoeff),
1400  fFactor(other.fFactor),
1401  fTmin(other.fTmin),
1402  fTmax(other.fTmax),
1403  fCoeff(NULL)
1404 {
1405  if (fCsize) {
1406  fCoeff = new BtCoef_t[fCsize];
1407  for (Int_t i=0; i<fNcoeff; i++) {
1408  fCoeff[i].cn = other.fCoeff[i].cn;
1409  fCoeff[i].lambda = other.fCoeff[i].lambda;
1410  }
1411  }
1412 }
1413 
1414 ////////////////////////////////////////////////////////////////////////////////
1415 /// Destructor.
1416 
1418 {
1419  if (fCoeff) delete [] fCoeff;
1420 }
1421 
1422 ////////////////////////////////////////////////////////////////////////////////
1423 /// Assignment.
1424 
1426 {
1427  if (this == &other) return *this;
1428  TObject::operator=(other);
1429  TAttLine::operator=(other);
1430  TAttFill::operator=(other);
1431  TAttMarker::operator=(other);
1432  fElem = other.fElem;
1433  fElemTop = other.fElemTop;
1434  if (fCoeff) delete [] fCoeff;
1435  fCoeff = 0;
1436  fCsize = other.fCsize;
1437  fNcoeff = other.fNcoeff;
1438  fFactor = other.fFactor;
1439  fTmin = other.fTmin;
1440  fTmax = other.fTmax;
1441  if (fCsize) {
1442  fCoeff = new BtCoef_t[fCsize];
1443  for (Int_t i=0; i<fNcoeff; i++) {
1444  fCoeff[i].cn = other.fCoeff[i].cn;
1445  fCoeff[i].lambda = other.fCoeff[i].lambda;
1446  }
1447  }
1448  return *this;
1449 }
1450 
1451 ////////////////////////////////////////////////////////////////////////////////
1452 /// Addition of other solution.
1453 
1455 {
1456  if (other.GetElement() != fElem) {
1457  Error("operator+=", "Cannot add 2 solutions for different elements");
1458  return *this;
1459  }
1460  Int_t i,j;
1461  BtCoef_t *coeff = fCoeff;
1462  Int_t ncoeff = fNcoeff + other.fNcoeff;
1463  if (ncoeff > fCsize) {
1464  fCsize = ncoeff;
1465  coeff = new BtCoef_t[ncoeff];
1466  for (i=0; i<fNcoeff; i++) {
1467  coeff[i].cn = fCoeff[i].cn;
1468  coeff[i].lambda = fCoeff[i].lambda;
1469  }
1470  delete [] fCoeff;
1471  fCoeff = coeff;
1472  }
1473  ncoeff = fNcoeff;
1474  for (j=0; j<other.fNcoeff; j++) {
1475  for (i=0; i<fNcoeff; i++) {
1476  if (coeff[i].lambda == other.fCoeff[j].lambda) {
1477  coeff[i].cn += fFactor * other.fCoeff[j].cn;
1478  break;
1479  }
1480  }
1481  if (i == fNcoeff) {
1482  coeff[ncoeff].cn = fFactor * other.fCoeff[j].cn;
1483  coeff[ncoeff].lambda = other.fCoeff[j].lambda;
1484  ncoeff++;
1485  }
1486  }
1487  fNcoeff = ncoeff;
1488  return *this;
1489 }
1490 ////////////////////////////////////////////////////////////////////////////////
1491 /// Find concentration of the element at a given time.
1492 
1494 {
1495  Double_t conc = 0.;
1496  for (Int_t i=0; i<fNcoeff; i++)
1497  conc += fCoeff[i].cn * TMath::Exp(-fCoeff[i].lambda * time);
1498  return conc;
1499 }
1500 
1501 ////////////////////////////////////////////////////////////////////////////////
1502 /// Draw the solution of Bateman equation versus time.
1503 
1505 {
1506  if (!gGeoManager) return;
1507  gGeoManager->GetGeomPainter()->DrawBatemanSol(this, option);
1508 }
1509 
1510 ////////////////////////////////////////////////////////////////////////////////
1511 /// Find the solution for the Bateman equations corresponding to the decay
1512 /// chain described by an array ending with element X.
1513 /// A->B->...->X
1514 /// Cn = SUM [Ain * exp(-LMBDi*t)];
1515 /// Cn - concentration Nx/Na
1516 /// n - order of X in chain (A->B->X => n=3)
1517 /// LMBDi - decay constant for element of order i in the chain
1518 /// Ain = LMBD1*...*LMBD(n-1) * br1*...*br(n-1)/(LMBD1-LMBDi)...(LMBDn-LMBDi)
1519 /// bri - branching ratio for decay Ei->Ei+1
1520 
1522 {
1523  fNcoeff = 0;
1524  if (!array || !array->GetEntriesFast()) return;
1525  Int_t n = array->GetEntriesFast();
1526  TGeoDecayChannel *dc = (TGeoDecayChannel*)array->At(n-1);
1527  TGeoElementRN *elem = dc->Daughter();
1528  if (elem != fElem) {
1529  Error("FindSolution", "Last element in the list must be %s\n", fElem->GetName());
1530  return;
1531  }
1532  Int_t i,j;
1533  Int_t order = n+1;
1534  if (!fCoeff) {
1535  fCsize = order;
1536  fCoeff = new BtCoef_t[fCsize];
1537  }
1538  if (fCsize < order) {
1539  delete [] fCoeff;
1540  fCsize = order;
1541  fCoeff = new BtCoef_t[fCsize];
1542  }
1543 
1544  Double_t *lambda = new Double_t[order];
1545  Double_t *br = new Double_t[n];
1546  Double_t halflife;
1547  for (i=0; i<n; i++) {
1548  dc = (TGeoDecayChannel*)array->At(i);
1549  elem = dc->Parent();
1550  br[i] = 0.01 * dc->BranchingRatio();
1551  halflife = elem->HalfLife();
1552  if (halflife==0.) halflife = 1.e-30;
1553  if (elem->Stable()) lambda[i] = 0.;
1554  else lambda[i] = TMath::Log(2.)/halflife;
1555  if (i==n-1) {
1556  elem = dc->Daughter();
1557  halflife = elem->HalfLife();
1558  if (halflife==0.) halflife = 1.e-30;
1559  if (elem->Stable()) lambda[n] = 0.;
1560  else lambda[n] = TMath::Log(2.)/halflife;
1561  }
1562  }
1563  // Check if we have equal lambdas
1564  for (i=0; i<order-1; i++) {
1565  for (j=i+1; j<order; j++) {
1566  if (lambda[j] == lambda[i]) lambda[j] += 0.001*lambda[j];
1567  }
1568  }
1569  Double_t ain;
1570  Double_t pdlambda, plambdabr=1.;
1571  for (j=0; j<n; j++) plambdabr *= lambda[j]*br[j];
1572  for (i=0; i<order; i++) {
1573  pdlambda = 1.;
1574  for (j=0; j<n+1; j++) {
1575  if (j == i) continue;
1576  pdlambda *= lambda[j] - lambda[i];
1577  }
1578  if (pdlambda == 0.) {
1579  Error("FindSolution", "pdlambda=0 !!!");
1580  delete [] lambda;
1581  delete [] br;
1582  return;
1583  }
1584  ain = plambdabr/pdlambda;
1585  fCoeff[i].cn = ain;
1586  fCoeff[i].lambda = lambda[i];
1587  }
1588  fNcoeff = order;
1589  Normalize(fFactor);
1590  delete [] lambda;
1591  delete [] br;
1592 }
1593 
1594 ////////////////////////////////////////////////////////////////////////////////
1595 /// Normalize all coefficients with a given factor.
1596 
1598 {
1599  for (Int_t i=0; i<fNcoeff; i++) fCoeff[i].cn *= factor;
1600 }
1601 
1602 ////////////////////////////////////////////////////////////////////////////////
1603 /// Print concentration evolution.
1604 
1605 void TGeoBatemanSol::Print(Option_t * /*option*/) const
1606 {
1607  TString formula;
1608  formula.Form("N[%s]/N[%s] = ", fElem->GetName(), fElemTop->GetName());
1609  for (Int_t i=0; i<fNcoeff; i++) {
1610  if (i == fNcoeff-1) formula += TString::Format("%g*exp(-%g*t)", fCoeff[i].cn, fCoeff[i].lambda);
1611  else formula += TString::Format("%g*exp(-%g*t) + ", fCoeff[i].cn, fCoeff[i].lambda);
1612  }
1613  printf("%s\n", formula.Data());
1614 }
1615 
static void DecayName(UInt_t decay, TString &name)
Returns decay name.
Int_t GetZ() const
Definition: TGeoElement.h:120
TString fTitle
Definition: TNamed.h:33
TObjArray * GetBranch() const
Definition: TGeoElement.h:352
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
void FillPopulation(TObjArray *population, Double_t precision=0.001, Double_t factor=1.)
Fills the input array with the set of RN elements resulting from the decay of this one...
static TGeoElementRN * ReadElementRN(const char *record, Int_t &ndecays)
Create element from line record.
void Normalize(Double_t factor)
Normalize all coefficients with a given factor.
An array of TObjects.
Definition: TObjArray.h:37
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive for RN elements.
static TGeoDecayChannel * ReadDecay(const char *record)
Create element from line record.
Table of elements.
Definition: TGeoElement.h:369
Double_t fHalfLife
Definition: TGeoElement.h:145
void AddDecay(Int_t decay, Int_t diso, Double_t branchingRatio, Double_t qValue)
Adds a decay mode for this element.
Double_t Log(Double_t x)
Definition: TMath.h:759
Double_t fCoulomb
Definition: TGeoElement.h:51
Double_t fTmin
Definition: TGeoElement.h:293
TLine * line
virtual void Draw(Option_t *option="")
Draw the solution of Bateman equation versus time.
~TGeoBatemanSol()
Destructor.
TGeoElementRN * Daughter() const
Definition: TGeoElement.h:260
const char Option_t
Definition: RtypesCore.h:62
virtual void Print(Option_t *opt=" ") const
Prints decay info.
TObjArray * fListRN
Definition: TGeoElement.h:377
virtual void Print(Option_t *option="") const
Print table of elements.
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:355
TGeoElement()
Default constructor.
Definition: TGeoElement.cxx:88
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
TGeoElementRN * operator()()
() operator.
static const Int_t gDecayDeltaA[gMaxDecay]
Definition: TGeoElement.cxx:72
TGeoElementRN * fParent
Definition: TGeoElement.h:223
TGeoElementTable * GetElementTable()
Returns material table. Creates it if not existing.
TGeoElementRN * GetElement() const
Definition: TGeoElement.h:310
static const char gElName[gMaxElem][3]
Definition: TGeoElement.cxx:55
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1113
Double_t fFactor
Definition: TGeoElement.h:292
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
void AddRatio(TGeoBatemanSol &ratio)
Adds a proportion ratio to the existing one.
UInt_t Decay() const
Definition: TGeoElement.h:256
Basic string class.
Definition: TString.h:131
Double_t GetA() const
Definition: TGeoElement.h:122
int Int_t
Definition: RtypesCore.h:41
TGeoElementTable & operator=(const TGeoElementTable &)
assignment operator
bool Bool_t
Definition: RtypesCore.h:59
Int_t GetNdecays() const
Get number of decay channels of this element.
void ComputeDerivedQuantities()
Calculate properties for an atomic number.
virtual ~TGeoElementRN()
Destructor.
void AddIsotope(TGeoIsotope *isotope)
Add isotope to the table.
Bool_t CheckDecays() const
Check if all decay chain of the element is OK.
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TObject.cxx:550
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
TGeoIsotope * GetIsotope(Int_t i) const
Return i-th isotope in the element.
TGeoIsotope * FindIsotope(const char *name) const
Find existing isotope by name. Not optimized for a big number of isotopes.
Double_t fTH_F
Definition: TGeoElement.h:148
virtual void Print(Option_t *option="") const
Print info about the current decay branch.
virtual void DecayShift(Int_t &dA, Int_t &dZ, Int_t &dI) const
Returns variation in A, Z and Iso after decay.
TObjArray * fIsotopes
Definition: TGeoElement.h:378
TObjArray * fDecays
Definition: TGeoElement.h:155
Double_t * fAbundances
Definition: TGeoElement.h:50
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
TGeoElementRN * fElemTop
Definition: TGeoElement.h:289
const TGeoElementRN * fTop
Definition: TGeoElement.h:331
TObjArray * fList
Definition: TGeoElement.h:376
virtual Int_t ENDFCode() const
Definition: TGeoElement.h:178
Marker Attributes class.
Definition: TAttMarker.h:19
TGeoElementRN * Parent() const
Definition: TGeoElement.h:261
static TGeoIsotope * FindIsotope(const char *name)
Find existing isotope by name.
Bool_t HasIsotopes() const
Definition: TGeoElement.h:85
virtual ~TGeoElementTable()
destructor
TGeoBatemanSol & operator=(const TGeoBatemanSol &other)
Assignment.
Fill Area Attributes class.
Definition: TAttFill.h:19
TObjArray * Decays() const
Definition: TGeoElement.h:195
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2286
void MakeName(Int_t a, Int_t z, Int_t iso)
Generate a default name for the element.
Double_t fA
Definition: TGeoElement.h:48
Int_t GetN() const
Definition: TGeoElement.h:121
Double_t Neff() const
Returns effective number of nucleons.
Double_t fRadTsai
Definition: TGeoElement.h:52
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive for decays.
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual TObject * FindObject(const char *name) const
Find an object in this collection using its name.
Definition: TObjArray.cxx:414
TObjArray * fBranch
Definition: TGeoElement.h:333
Double_t Concentration(Double_t time) const
Find concentration of the element at a given time.
void AddIsotope(TGeoIsotope *isotope, Double_t relativeAbundance)
Add an isotope for this element. All isotopes have to be isotopes of the same element.
Base class for chemical elements.
Definition: TGeoElement.h:36
static TGeoElementTable * GetElementTable()
Returns pointer to the table.
virtual const char * PrependPathName(const char *dir, TString &name)
Concatenate a directory and a file name.
Definition: TSystem.cxx:1062
static const char gLevName[gMaxLevel]
Definition: TGeoElement.cxx:81
Double_t fTmax
Definition: TGeoElement.h:294
virtual void DrawBatemanSol(TGeoBatemanSol *sol, Option_t *option="")=0
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.h:271
TGeoElementRN()
Default constructor.
TGeoElementRN * Down(Int_t ibranch)
Go downwards from current level via ibranch as low in the tree as possible.
BtCoef_t * fCoeff
Definition: TGeoElement.h:295
Double_t HalfLife() const
Definition: TGeoElement.h:186
Bool_t HasRNElements() const
Definition: TGeoElement.h:415
static const Int_t gDecayDeltaZ[gMaxDecay]
Definition: TGeoElement.cxx:77
void SetDaughter(TGeoElementRN *daughter)
Definition: TGeoElement.h:265
Double_t fLimitRatio
Definition: TGeoElement.h:335
TObjArray * fIsotopes
Definition: TGeoElement.h:49
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:234
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
auto * a
Definition: textangle.C:12
static const Int_t gMaxDecay
Definition: TGeoElement.cxx:53
virtual const char * GetName() const
Returns name of decay.
void FindSolution(const TObjArray *array)
Find the solution for the Bateman equations corresponding to the decay chain described by an array en...
virtual TObject * RemoveAt(Int_t idx)
Remove object at index idx.
Definition: TObjArray.cxx:678
const TGeoElementRN * fElem
Definition: TGeoElement.h:332
Double_t fLevel
Definition: TGeoElement.h:143
virtual void Print(Option_t *option="") const
Print info about the element;.
Double_t fTG_F
Definition: TGeoElement.h:149
Double_t fBranchingRatio
Definition: TGeoElement.h:221
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2264
unsigned int UInt_t
Definition: RtypesCore.h:42
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
void AddElementRN(TGeoElementRN *elem)
Add a radionuclide to the table and map it.
Ssiz_t Length() const
Definition: TString.h:405
TGeoDecayChannel & operator=(const TGeoDecayChannel &dc)
Assignment.
static const Int_t gMaxLevel
Definition: TGeoElement.cxx:52
Iterator for decay branches.
Definition: TGeoElement.h:328
TGeoElement * FindElement(const char *name) const
Search an element by symbol or full name Exact matching.
TGeoElementRN * GetElementRN(Int_t ENDFcode) const
Retrieve a radionuclide by ENDF code.
void AddElement(const char *name, const char *title, Int_t z, Double_t a)
Add an element to the table. Obsolete.
TString fName
Definition: TNamed.h:32
constexpr Double_t Na()
Avogadro constant (Avogadro&#39;s Number) in .
Definition: TMath.h:283
TGeoBatemanSol & operator+=(const TGeoBatemanSol &other)
Addition of other solution.
Double_t fDeltaM
Definition: TGeoElement.h:144
static const char * gDecayName[gMaxDecay+1]
Definition: TGeoElement.cxx:67
constexpr Double_t E()
Base of natural log: .
Definition: TMath.h:97
const Bool_t kFALSE
Definition: RtypesCore.h:88
Double_t fNatAbun
Definition: TGeoElement.h:146
Double_t Exp(Double_t x)
Definition: TMath.h:726
void SetUsed(Bool_t flag=kTRUE)
Definition: TGeoElement.h:91
Bool_t HasDefaultElements() const
Definition: TGeoElement.h:414
#define ClassImp(name)
Definition: Rtypes.h:359
void ComputeLradTsaiFactor()
Compute Tsai&#39;s Expression for the Radiation Length (Phys Rev. D50 3-1 (1994) page 1254) ...
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:562
ElementRNMap_t fElementsRN
Definition: TGeoElement.h:382
Double_t fQvalue
Definition: TGeoElement.h:222
Int_t DecayResult(TGeoDecayChannel *dc) const
Returns ENDF code of decay result.
double Double_t
Definition: RtypesCore.h:55
TGeoElementRN * Next()
Return next element.
Double_t fTH_S
Definition: TGeoElement.h:150
Bool_t CheckTable() const
Checks status of element table.
virtual void Print(Option_t *option="") const
Print this isotope.
static const TString & GetEtcDir()
Get the sysconfig directory in the installation. Static utility function.
Definition: TROOT.cxx:2991
TGeoElementRN * fElem
Definition: TGeoElement.h:288
virtual ~TGeoElemIter()
Destructor.
static constexpr double alpha_rcl2
void SetParent(TGeoElementRN *parent)
Definition: TGeoElement.h:264
void ComputeCoulombFactor()
Compute Coulomb correction factor (Phys Rev. D50 3-1 (1994) page 1254)
static constexpr double s
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
void SetDefined(Bool_t flag=kTRUE)
Definition: TGeoElement.h:90
Int_t fNisotopes
Definition: TGeoElement.h:47
TGeoIsotope()
Dummy I/O constructor.
Binding & operator=(OUT(*fun)(void))
A decay channel for a radionuclide.
Definition: TGeoElement.h:216
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Double_t GetSpecificActivity() const
Get the activity in Bq of a gram of material made from this element.
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Double_t fRatio
Definition: TGeoElement.h:336
static constexpr double fine_structure_const
virtual void Print(Option_t *option="") const
Print concentration evolution.
TGeoElementRN * Up()
Go upwards from the current location until the next branching, then down.
void ImportElementsRN()
Creates the list of radionuclides.
virtual void Print(Option_t *option="") const
Print this isotope.
Int_t GetIndex() const
Get index of this channel in the list of decays of the parent nuclide.
Class representing a radionuclide.
Definition: TGeoElement.h:138
TGeoElementRN * fDaughter
Definition: TGeoElement.h:224
void BuildDefaultElements()
Creates the default element table.
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
TGeoElementTable()
default constructor
static Int_t ENDF(Int_t a, Int_t z, Int_t iso)
Definition: TGeoElement.h:175
TGeoBatemanSol * fRatio
Definition: TGeoElement.h:153
void Add(TObject *obj)
Definition: TObjArray.h:73
Double_t fTG_S
Definition: TGeoElement.h:151
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:908
Double_t GetRelativeAbundance(Int_t i) const
Return relative abundance of i-th isotope in this element.
static const Int_t gMaxElem
Definition: TGeoElement.cxx:51
virtual void Print(Option_t *option="") const
Default print for collections, calls Print(option, 1).
Bool_t Stable() const
Definition: TGeoElement.h:194
const Bool_t kTRUE
Definition: RtypesCore.h:87
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
const Int_t n
Definition: legend1.C:16
Line Attributes class.
Definition: TAttLine.h:18
Double_t BranchingRatio() const
Definition: TGeoElement.h:257
char name[80]
Definition: TGX11.cxx:109
void ResetRatio()
Clears the existing ratio.
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
Double_t fA
Definition: TGeoElement.h:113
const char * Data() const
Definition: TString.h:364
void ExportElementsRN(const char *filename="")
Export radionuclides in a file.
TGeoElemIter & operator=(const TGeoElemIter &iter)
Assignment.