ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RadioNuclides.C
Go to the documentation of this file.
1 void DrawPopulation(TObjArray *vect, TCanvas *can, Double_t tmin=0.,
2  Double_t tmax=0., Bool_t logx=kFALSE);
3 
5 {
6 // Macro that demonstrates usage of radioactive elements/materials/mixtures
7 // with TGeo package.
8 //
9 // A radionuclide (TGeoElementRN) derives from the class TGeoElement and
10 // provides additional information related to its radioactive properties and
11 // decay modes.
12 //
13 // The radionuclides table is loaded on demand by any call:
14 // TGeoElementRN *TGeoElementTable::GetElementRN(Int_t atomic_number,
15 // Int_t atomic_charge,
16 // Int_t isomeric_number)
17 // The isomeric number is optional and the default value is 0.
18 //
19 // To create a radioactive material based on a radionuclide, one should use the
20 // constructor:
21 // TGeoMaterial(const char *name, TGeoElement *elem, Double_t density)
22 // To create a radioactive mixture, one can use radionuclides as well as stable
23 // elements:
24 // TGeoMixture(const char *name, Int_t nelements, Double_t density);
25 // TGeoMixture::AddElement(TGeoElement *elem, Double_t weight_fraction);
26 // Once defined, one can retrieve the time evolution for the radioactive
27 // materials/mixtures by using one of the 2 methods:
28 //
29 // void TGeoMaterial::FillMaterialEvolution(TObjArray *population,
30 // Double_t precision=0.001)
31 // To use this method, one has to provide an empty TObjArray object that will
32 // be filled with all elements coming from the decay chain of the initial
33 // radionuclides contained by the material/mixture. The precision represent the
34 // cumulative branching ratio for which decay products are still considered.
35 // The POPULATION list may contain stable elements as well as radionuclides,
36 // depending on the initial elements. To test if an element is a radionuclide:
37 // Bool_t TGeoElement::IsRadioNuclide() const
38 // All radionuclides in the output population list have attached objects that
39 // represent the time evolution of their fraction of nuclei with respect to the
40 // top radionuclide in the decay chain. These objects (Bateman solutions) can be
41 // retrieved and drawn:
42 // TGeoBatemanSol *TGeoElementRN::Ratio();
43 // void TGeoBatemanSol::Draw();
44 //
45 // Another method allows to create the evolution of a given radioactive
46 // material/mixture at a given moment in time:
47 // TGeoMaterial::DecayMaterial(Double_t time, Double_t precision=0.001)
48 // The method will create the mixture that result from the decay of a initial
49 // material/mixture at TIME, while all resulting elements having a fractional
50 // weight less than PRECISION are excluded.
51 //Author: Mihaela Gheata
52 
53  TGeoManager *geom = new TGeoManager("","");
55  TGeoElementRN *c14 = table->GetElementRN(14,6);
56  TGeoElementRN *el1 = table->GetElementRN(53,20);
57  TGeoElementRN *el2 = table->GetElementRN(78,38);
58  // Radioactive material
59  TGeoMaterial *mat = new TGeoMaterial("C14", c14, 1.3);
60  printf("___________________________________________________________\n");
61  printf("Radioactive material:\n");
62  mat->Print();
63  Double_t time = 1.5e11; // seconds
64  TGeoMaterial *decaymat = mat->DecayMaterial(time);
65  printf("Radioactive material evolution after %g years:\n", time/3.1536e7);
66  decaymat->Print();
67  //Radioactive mixture
68  TGeoMixture *mix = new TGeoMixture("mix", 2, 7.3);
69  mix->AddElement(el1, 0.35);
70  mix->AddElement(el2, 0.65);
71  printf("___________________________________________________________\n");
72  printf("Radioactive mixture:\n");
73  mix->Print();
74  time = 1000.;
75  decaymat = mix->DecayMaterial(time);
76  printf("Radioactive mixture evolution after %g seconds:\n", time);
77  decaymat->Print();
78  TObjArray *vect = new TObjArray();
79  TCanvas *c1 = new TCanvas("c1","C14 decay", 800,600);
80  c1->SetGrid();
81  mat->FillMaterialEvolution(vect);
82  DrawPopulation(vect, c1, 0, 1.4e12);
83  TLatex *tex = new TLatex(8.35e11,0.564871,"C_{N^{14}_{7}}");
84  tex->SetTextSize(0.0388601);
85  tex->SetLineWidth(2);
86  tex->Draw();
87  tex = new TLatex(3.33e11,0.0620678,"C_{C^{14}_{6}}");
88  tex->SetTextSize(0.0388601);
89  tex->SetLineWidth(2);
90  tex->Draw();
91  tex = new TLatex(9.4e11,0.098,"C_{X}=#frac{N_{X}(t)}{N_{0}(t=0)}=\
92  #sum_{j}#alpha_{j}e^{-#lambda_{j}t}");
93  tex->SetTextSize(0.0388601);
94  tex->SetLineWidth(2);
95  tex->Draw();
96  TPaveText *pt = new TPaveText(2.6903e+11,0.0042727,1.11791e+12,0.0325138,"br");
97  pt->SetFillColor(5);
98  pt->SetTextAlign(12);
99  pt->SetTextColor(4);
100  text = pt->AddText("Time evolution of a population of radionuclides.");
101  text = pt->AddText("The concentration of a nuclide X represent the ");
102  text = pt->AddText("ratio between the number of X nuclei and the ");
103  text = pt->AddText("number of nuclei of the top element of the decay");
104  text = pt->AddText("from which X derives from at T=0. ");
105  pt->Draw();
106  c1->Modified();
107  vect->Clear();
108  TCanvas *c2 = new TCanvas("c2","Mixture decay", 1000,800);
109  c2->SetGrid();
110  mix->FillMaterialEvolution(vect);
111  DrawPopulation(vect, c2, 0.01, 1000., kTRUE);
112  tex = new TLatex(0.019,0.861,"C_{Ca^{53}_{20}}");
113  tex->SetTextSize(0.0388601);
114  tex->SetTextColor(1);
115  tex->Draw();
116  tex = new TLatex(0.0311,0.078064,"C_{Sc^{52}_{21}}");
117  tex->SetTextSize(0.0388601);
118  tex->SetTextColor(2);
119  tex->Draw();
120  tex = new TLatex(0.1337,0.010208,"C_{Ti^{52}_{22}}");
121  tex->SetTextSize(0.0388601);
122  tex->SetTextColor(3);
123  tex->Draw();
124  tex = new TLatex(1.54158,0.00229644,"C_{V^{52}_{23}}");
125  tex->SetTextSize(0.0388601);
126  tex->SetTextColor(4);
127  tex->Draw();
128  tex = new TLatex(25.0522,0.00135315,"C_{Cr^{52}_{24}}");
129  tex->SetTextSize(0.0388601);
130  tex->SetTextColor(5);
131  tex->Draw();
132  tex = new TLatex(0.1056,0.5429,"C_{Sc^{53}_{21}}");
133  tex->SetTextSize(0.0388601);
134  tex->SetTextColor(6);
135  tex->Draw();
136  tex = new TLatex(0.411,0.1044,"C_{Ti^{53}_{22}}");
137  tex->SetTextSize(0.0388601);
138  tex->SetTextColor(7);
139  tex->Draw();
140  tex = new TLatex(2.93358,0.0139452,"C_{V^{53}_{23}}");
141  tex->SetTextSize(0.0388601);
142  tex->SetTextColor(8);
143  tex->Draw();
144  tex = new TLatex(10.6235,0.00440327,"C_{Cr^{53}_{24}}");
145  tex->SetTextSize(0.0388601);
146  tex->SetTextColor(9);
147  tex->Draw();
148  tex = new TLatex(15.6288,0.782976,"C_{Sr^{78}_{38}}");
149  tex->SetTextSize(0.0388601);
150  tex->SetTextColor(1);
151  tex->Draw();
152  tex = new TLatex(20.2162,0.141779,"C_{Rb^{78}_{37}}");
153  tex->SetTextSize(0.0388601);
154  tex->SetTextColor(2);
155  tex->Draw();
156  tex = new TLatex(32.4055,0.0302101,"C_{Kr^{78}_{36}}");
157  tex->SetTextSize(0.0388601);
158  tex->SetTextColor(3);
159  tex->Draw();
160  tex = new TLatex(117.,1.52,"C_{X}=#frac{N_{X}(t)}{N_{0}(t=0)}=#sum_{j}\
161  #alpha_{j}e^{-#lambda_{j}t}");
162  tex->SetTextSize(0.03);
163  tex->SetLineWidth(2);
164  tex->Draw();
165  TArrow *arrow = new TArrow(0.0235313,0.74106,0.0385371,0.115648,0.02,">");
166  arrow->SetFillColor(1);
167  arrow->SetFillStyle(1001);
168  arrow->SetLineWidth(2);
169  arrow->SetAngle(30);
170  arrow->Draw();
171  arrow = new TArrow(0.0543138,0.0586338,0.136594,0.0146596,0.02,">");
172  arrow->SetFillColor(1);
173  arrow->SetFillStyle(1001);
174  arrow->SetLineWidth(2);
175  arrow->SetAngle(30);
176  arrow->Draw();
177  arrow = new TArrow(0.31528,0.00722919,1.29852,0.00306079,0.02,">");
178  arrow->SetFillColor(1);
179  arrow->SetFillStyle(1001);
180  arrow->SetLineWidth(2);
181  arrow->SetAngle(30);
182  arrow->Draw();
183  arrow = new TArrow(4.13457,0.00201942,22.5047,0.00155182,0.02,">");
184  arrow->SetFillColor(1);
185  arrow->SetFillStyle(1001);
186  arrow->SetLineWidth(2);
187  arrow->SetAngle(30);
188  arrow->Draw();
189  arrow = new TArrow(0.0543138,0.761893,0.0928479,0.67253,0.02,">");
190  arrow->SetFillColor(1);
191  arrow->SetFillStyle(1001);
192  arrow->SetLineWidth(2);
193  arrow->SetAngle(30);
194  arrow->Draw();
195  arrow = new TArrow(0.238566,0.375717,0.416662,0.154727,0.02,">");
196  arrow->SetFillColor(1);
197  arrow->SetFillStyle(1001);
198  arrow->SetLineWidth(2);
199  arrow->SetAngle(30);
200  arrow->Draw();
201  arrow = new TArrow(0.653714,0.074215,2.41863,0.0213142,0.02,">");
202  arrow->SetFillColor(1);
203  arrow->SetFillStyle(1001);
204  arrow->SetLineWidth(2);
205  arrow->SetAngle(30);
206  arrow->Draw();
207  arrow = new TArrow(5.58256,0.00953882,10.6235,0.00629343,0.02,">");
208  arrow->SetFillColor(1);
209  arrow->SetFillStyle(1001);
210  arrow->SetLineWidth(2);
211  arrow->SetAngle(30);
212  arrow->Draw();
213  arrow = new TArrow(22.0271,0.601935,22.9926,0.218812,0.02,">");
214  arrow->SetFillColor(1);
215  arrow->SetFillStyle(1001);
216  arrow->SetLineWidth(2);
217  arrow->SetAngle(30);
218  arrow->Draw();
219  arrow = new TArrow(27.2962,0.102084,36.8557,0.045686,0.02,">");
220  arrow->SetFillColor(1);
221  arrow->SetFillStyle(1001);
222  arrow->SetLineWidth(2);
223  arrow->SetAngle(30);
224  arrow->Draw();
225 }
226 
227 void DrawPopulation(TObjArray *vect, TCanvas *can, Double_t tmin,
228  Double_t tmax, Bool_t logx)
229 {
230  Int_t n = vect->GetEntriesFast();
231  TGeoElementRN *elem;
232  TGeoBatemanSol *sol;
233  can->SetLogy();
234 
235  if (logx) can->SetLogx();
236 
237 
238  for (Int_t i=0; i<n; i++) {
239  TGeoElement *el = (TGeoElement*)vect->At(i);
240  if (!el->IsRadioNuclide()) continue;
241  TGeoElementRN *elem = (TGeoElementRN *)el;
242  TGeoBatemanSol *sol = elem->Ratio();
243  if (sol) {
244  sol->SetLineColor(1+(i%9));
245  sol->SetLineWidth(2);
246  if (tmax>0.) sol->SetRange(tmin,tmax);
247  if (i==0) {
248  sol->Draw();
249  TF1 *func = (TF1*)can->FindObject(
250  Form("conc%s",sol->GetElement()->GetName()));
251  if (func) {
252  if (!strcmp(can->GetName(),"c1")) func->SetTitle(
253  "Concentration of C14 derived elements;time[s];Ni/N0(C14)");
254  else func->SetTitle(
255  "Concentration of elements derived from mixture Ca53+Sr78;\
256  time[s];Ni/N0(Ca53)");
257  }
258  }
259  else sol->Draw("SAME");
260  }
261  }
262 }
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
An array of TObjects.
Definition: TObjArray.h:39
virtual void Draw(Option_t *option="")
Draw this pavetext with its current attributes.
Definition: TPaveText.cxx:211
virtual void Draw(Option_t *option="")
Draw the solution of Bateman equation versus time.
virtual void Clear(Option_t *option="")
Remove all objects from the array.
Definition: TObjArray.cxx:298
TCanvas * c1
Definition: legend1.C:2
TGeoElementTable * GetElementTable()
Returns material table. Creates it if not existing.
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
Definition: TPaveText.cxx:160
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition: TObject.cxx:254
virtual void SetFillStyle(Style_t fstyle)
Definition: TAttFill.h:52
Int_t GetEntriesFast() const
Definition: TObjArray.h:66
TGeoElementRN * GetElementRN(Int_t ENDFcode) const
Retreive a radionuclide by ENDF code.
void RadioNuclides()
Definition: RadioNuclides.C:4
virtual void FillMaterialEvolution(TObjArray *population, Double_t precision=0.001)
Fills a user array with all the elements deriving from the possible decay of the top elements composi...
TGeoBatemanSol * Ratio() const
Definition: TGeoElement.h:195
TArrow * arrow
virtual void Draw(Option_t *option="")
Draw this arrow with its current attributes.
Definition: TArrow.cxx:121
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)
Definition: TPad.h:326
To draw Mathematical Formula.
Definition: TLatex.h:33
virtual void SetLogx(Int_t value=1)
Set Lin/Log scale for X.
Definition: TPad.cxx:5301
virtual void FillMaterialEvolution(TObjArray *population, Double_t precision=0.001)
Fills a user array with all the elements deriving from the possible decay of the top element composin...
virtual void SetTextAlign(Short_t align=11)
Definition: TAttText.h:55
const char * GetName() const
Returns name of object.
Definition: TPad.h:255
void AddElement(Double_t a, Double_t z, Double_t weight)
add an element to the mixture using fraction by weight Check if the element is already defined ...
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual TObject * FindObject(const char *name) const
Search if object named name is inside this pad or in pads inside this pad.
Definition: TPad.cxx:2491
TPaveText * pt
virtual TGeoMaterial * DecayMaterial(Double_t time, Double_t precision=0.001)
Create the material representing the decay product of this material at a given time.
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition: TF1.cxx:3227
The Canvas class.
Definition: TCanvas.h:48
return c2
Definition: legend2.C:14
A Pave (see TPave) with text, lines or/and boxes inside.
Definition: TPaveText.h:35
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:556
double Double_t
Definition: RtypesCore.h:55
TText * text
virtual Bool_t IsRadioNuclide() const
Definition: TGeoElement.h:88
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
double func(double *x, double *p)
Definition: stressTF1.cxx:213
virtual void Print(const Option_t *option="") const
print characteristics of this material
TGeoElementRN * GetElement() const
Definition: TGeoElement.h:308
1-Dim function class
Definition: TF1.h:149
virtual void SetTextColor(Color_t tcolor=1)
Definition: TAttText.h:57
virtual TGeoMaterial * DecayMaterial(Double_t time, Double_t precision=0.001)
Create the mixture representing the decay product of this material at a given time.
virtual void SetTextSize(Float_t tsize=1)
Definition: TAttText.h:60
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
Draw all kinds of Arrows.
Definition: TArrow.h:35
void DrawPopulation(TObjArray *vect, TCanvas *can, Double_t tmin=0., Double_t tmax=0., Bool_t logx=kFALSE)
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void Print(const Option_t *option="") const
print characteristics of this material
void SetRange(Double_t tmin=0., Double_t tmax=0.)
Definition: TGeoElement.h:312
const Int_t n
Definition: legend1.C:16
void Modified(Bool_t flag=1)
Definition: TPad.h:407
virtual void SetAngle(Float_t angle=60)
Definition: TArrow.h:64
virtual void SetLogy(Int_t value=1)
Set Lin/Log scale for Y.
Definition: TPad.cxx:5314