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