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