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