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