Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
advanced/microbeam/plot.C
이 파일의 문서화 페이지로 가기
1 // -------------------------------------------------------------------
2 // $Id: plot.C 69157 2013-04-19 14:05:11Z gcosmo $
3 // -------------------------------------------------------------------
4 //
5 // *********************************************************************
6 // To execute this macro under ROOT,
7 // 1 - launch ROOT (usually type 'root' at your machine's prompt)
8 // 2 - type '.X plot.C' at the ROOT session prompt
9 // This macro needs five files : dose.txt, stoppingPower.txt, range.txt,
10 // 3DDose.txt and beamPosition.txt
11 // written by S. Incerti and O. Boissonnade, 10/04/2006
12 // *********************************************************************
13 {
14 
15 gROOT->Reset();
16 
17 //****************
18 // DOSE IN NUCLEUS
19 //****************
20 
21 gStyle->SetOptStat(0000);
22 gStyle->SetOptFit();
23 gStyle->SetPalette(1);
24 gROOT->SetStyle("Plain");
26 
27 
28 c1 = new TCanvas ("c1","",20,20,1200,900);
29 c1->Divide(4,3);
30 
31 //*********************
32 // INTENSITY HISTOGRAMS
33 //*********************
34 
35 FILE * fp = fopen("phantom.dat","r");
37 
39 Float_t vox = 0, mat = 0;
41 
42 TH1F *h1 = new TH1F("h1","Nucleus marker intensity",100,1,300);
43 TH1F *h11 = new TH1F("h11 ","",100,1,300);
44 
45 TH1F *h2 = new TH1F("h2","Cytoplasm marker intensity",100,1,300);
46 TH1F *h20 = new TH1F("h20 ","",100,1,300);
47 
48 TNtuple *ntupleYXN = new TNtuple("NUCLEUS","ntuple","Y:X:vox");
49 TNtuple *ntupleZX = new TNtuple("CYTOPLASM","ntuple","Z:X:vox");
50 TNtuple *ntupleYX = new TNtuple("CYTOPLASM","ntuple","Y:X:vox");
51 
54 
55 while (1)
56  {
57  if ( nlines == 0 ) ncols = fscanf(fp,"%f %f %f",&tmp,&tmp,&tmp);
58  if ( nlines == 1 ) ncols = fscanf(fp,"%f %f %f",&voxelSizeX,&voxelSizeY,&voxelSizeZ);
59  if ( nlines == 2 ) ncols = fscanf(fp,"%f %f %f",&tmp,&tmp,&tmp);
60  if ( nlines >= 3 ) ncols = fscanf(fp,"%f %f %f %f %f %f", &X, &Y, &Z, &mat, &den, &vox);
61  if (ncols < 0) break;
62 
63  X= X*voxelSizeX;
64  Y= Y*voxelSizeY;
65  Z= Z*voxelSizeZ;
66 
67  if ( mat == 2 ) // noyau
68  {
69  if (den==1) h1->Fill( vox );
70  if (den==2) h11->Fill( vox );
71  ntupleYXN->Fill(Y,X,vox);
72  }
73  if ( mat == 1 ) // cytoplasm
74  {
75  if (den==1) h2->Fill( vox );
76  if (den==2) h20->Fill( vox );
77  ntupleZX->Fill(Z,X,vox);
78  ntupleYX->Fill(Y,X,vox);
79  }
80  nlines++;
81 
82  }
83 fclose(fp);
84 
85 // HISTO NUCLEUS
86 
87 c1->cd(1);
88  h1->Draw();
89  h1->GetXaxis()->SetLabelSize(0.025);
90  h1->GetYaxis()->SetLabelSize(0.025);
91  h1->GetXaxis()->SetTitleSize(0.035);
92  h1->GetYaxis()->SetTitleSize(0.035);
93  h1->GetXaxis()->SetTitleOffset(1.4);
94  h1->GetYaxis()->SetTitleOffset(1.4);
95  h1->GetXaxis()->SetTitle("Voxel intensity (0-255)");
96  h1->GetYaxis()->SetTitle("Number of events");
97  h1->SetLineColor(3);
98  h1->SetFillColor(3); // green
99 
100  h11->SetLineColor(8);
101  h11->SetFillColor(8); // dark green
102  h11->Draw("same");
103 
104 // HISTO CYTOPLASM
105 
106 c1->cd(5);
107  h2->Draw();
108  h2->GetXaxis()->SetLabelSize(0.025);
109  h2->GetYaxis()->SetLabelSize(0.025);
110  h2->GetXaxis()->SetTitleSize(0.035);
111  h2->GetYaxis()->SetTitleSize(0.035);
112  h2->GetXaxis()->SetTitleOffset(1.4);
113  h2->GetYaxis()->SetTitleOffset(1.4);
114  h2->GetXaxis()->SetTitle("Voxel intensity (0-255)");
115  h2->GetYaxis()->SetTitle("Number of events");
116  h2->SetLineColor(2);
117  h2->SetFillColor(2); // red
118 
119  h20->SetLineColor(5);
120  h20->SetFillColor(5); // yellow (nucleoli)
121  h20->Draw("same");
122 
123 //*************************
124 // CUMULATED CELL INTENSITY
125 //*************************
126 
127 gStyle->SetOptStat(0000);
128 gStyle->SetOptFit();
129 gStyle->SetPalette(1);
130 gROOT->SetStyle("Plain");
131 
132 //CYTOPLASM
133 
134 c1->cd(7); // axe YX
135  TH2F *hist = new TH2F("hist","hist",50,-20,20,50,-20,20);
136  ntupleYX->Draw("Y:X>>hist","vox","colz");
137  gPad->SetLogz();
138  hist->Draw("colz");
139  hist->GetXaxis()->SetLabelSize(0.025);
140  hist->GetYaxis()->SetLabelSize(0.025);
141  hist->GetZaxis()->SetLabelSize(0.025);
142  hist->GetXaxis()->SetTitleSize(0.035);
143  hist->GetYaxis()->SetTitleSize(0.035);
144  hist->GetXaxis()->SetTitleOffset(1.4);
145  hist->GetYaxis()->SetTitleOffset(1.4);
146  hist->GetXaxis()->SetTitle("Y (um)");
147  hist->GetYaxis()->SetTitle("X (um)");
148  hist->SetTitle("Cytoplasm intensity on transverse section");
149 
150 //NUCLEUS
151 
152 c1->cd(3); // axe YX
153  TH2F *hist2 = new TH2F("hist2","hist2",50,-20,20,50,-20,20);
154  ntupleYXN->Draw("Y:X>>hist2","vox","colz");
155  gPad->SetLogz();
156  hist2->Draw("colz");
157  hist2->GetXaxis()->SetLabelSize(0.025);
158  hist2->GetYaxis()->SetLabelSize(0.025);
159  hist2->GetZaxis()->SetLabelSize(0.025);
160  hist2->GetXaxis()->SetTitleSize(0.035);
161  hist2->GetYaxis()->SetTitleSize(0.035);
162  hist2->GetXaxis()->SetTitleOffset(1.4);
163  hist2->GetYaxis()->SetTitleOffset(1.4);
164  hist2->GetXaxis()->SetTitle("Y (um)");
165  hist2->GetYaxis()->SetTitle("X (um)");
166  hist2->SetTitle("Nucleus intensity on transverse section");
167 
168 //
169 
170 system ("rm -rf microbeam.root");
171 system ("hadd -O microbeam.root microbeam_*.root");
172 
173 TFile f("microbeam.root");
174 
175 TNtuple* ntuple0;
176 TNtuple* ntuple1;
177 TNtuple* ntuple2;
178 TNtuple* ntuple3;
179 TNtuple* ntuple4;
180 
181 ntuple0 = (TNtuple*)f.Get("ntuple0");
182 ntuple1 = (TNtuple*)f.Get("ntuple1");
183 ntuple2 = (TNtuple*)f.Get("ntuple2");
184 ntuple3 = (TNtuple*)f.Get("ntuple3");
185 ntuple4 = (TNtuple*)f.Get("ntuple4");
186 
187 TH1F *h1bis = new TH1F("h1bis","Dose distribution in Nucleus",100,0.001,1.);
188 TH1F *h10 = new TH1F("h10bis","Dose distribution in Cytoplasm",100,0.001,.2);
189 
190 c1->cd(2);
191 
192  ntuple3->Project("h1bis","doseN");
193  scale = 1/h1bis->Integral();
194  h1bis->Scale(scale);
195  h1bis->Draw();
196  h1bis->GetXaxis()->SetLabelSize(0.025);
197  h1bis->GetYaxis()->SetLabelSize(0.025);
198  h1bis->GetXaxis()->SetTitleSize(0.035);
199  h1bis->GetYaxis()->SetTitleSize(0.035);
200  h1bis->GetXaxis()->SetTitleOffset(1.4);
201  h1bis->GetYaxis()->SetTitleOffset(1.4);
202  h1bis->GetXaxis()->SetTitle("Absorbed dose (Gy)");
203  h1bis->GetYaxis()->SetTitle("Fraction of events");
204  h1bis->SetLineColor(3);
205  h1bis->SetFillColor(3);
206 
207 //*****************
208 // DOSE IN CYTOPLASM
209 //*****************
210 
211 c1->cd(6);
212  ntuple3->Project("h10bis","doseC");
213  scale = 1/h10->Integral();
214  h10->Scale(scale);
215  h10->Draw();
216  h10->GetXaxis()->SetLabelSize(0.025);
217  h10->GetYaxis()->SetLabelSize(0.025);
218  h10->GetXaxis()->SetTitleSize(0.035);
219  h10->GetYaxis()->SetTitleSize(0.035);
220  h10->GetXaxis()->SetTitleOffset(1.4);
221  h10->GetYaxis()->SetTitleOffset(1.4);
222  h10->GetXaxis()->SetTitle("Absorbed dose (Gy)");
223  h10->GetYaxis()->SetTitle("Fraction of events");
224  h10->SetLineColor(2);
225  h10->SetFillColor(2);
226 
227 //********************************
228 // STOPPING POWER AT CELL ENTRANCE
229 //********************************
230 
231 gStyle->SetOptStat(0000);
232 gStyle->SetOptFit();
233 gStyle->SetPalette(1);
234 gROOT->SetStyle("Plain");
235 
237 
238 TH1F *h2bis = new TH1F("h2bis","Beam stopping power at cell entrance",200,0,300);
239 
240 c1->cd(9);
241  ntuple0->Project("h2bis","sp");
242  scale = 1/h2bis->Integral();
243  h2bis->Scale(scale);
244  h2bis->Draw();
245  h2bis->GetXaxis()->SetLabelSize(0.025);
246  h2bis->GetYaxis()->SetLabelSize(0.025);
247  h2bis->GetXaxis()->SetTitleSize(0.035);
248  h2bis->GetYaxis()->SetTitleSize(0.035);
249  h2bis->GetXaxis()->SetTitleOffset(1.4);
250  h2bis->GetYaxis()->SetTitleOffset(1.4);
251  h2bis->GetXaxis()->SetTitle("dE/dx (keV/um)");
252  h2bis->GetYaxis()->SetTitle("Fraction of events");
253  h2bis->SetTitle("dE/dx at cell entrance");
254  h2bis->SetFillColor(4);
255  h2bis->SetLineColor(4);
256  h2bis->Fit("gaus");
257  gaus->SetLineColor(6);
258  h2bis->Fit("gaus");
259 
260 //**************
261 // RANGE IN CELL
262 //**************
263 
265 
266 // X position of target in World
267 Xc = -1295.59e3 - 955e3*sin(10*TMath::Pi()/180);
268 
269 // Z position of target in World
270 Zc = -1327e3 + 955e3*cos(10*TMath::Pi()/180);
271 
272 // Line alignment (cf MicrobeamEMField.cc)
273 Xc = Xc + 5.24*cos(10*TMath::Pi()/180);
274 Zc = Zc + 5.24*sin(10*TMath::Pi()/180);
275 
276 TNtuple *ntupleR = new TNtuple("Rmax","ntuple","Z2:Y2:X2");
278 ntuple2->SetBranchAddress("x",&x);
279 ntuple2->SetBranchAddress("y",&y);
280 ntuple2->SetBranchAddress("z",&z);
281 Int_t nentries = (Int_t)ntuple2->GetEntries();
282 for (Int_t i=0;i<nentries;i++)
283 {
284  ntuple2->GetEntry(i);
285  X1=x;
286  Y1=y;
287  Z1=z;
288  xx = X1-Xc;
289  zz = Z1-Zc;
290  Z2 = zz*cos(10*TMath::Pi()/180)-xx*sin(10*TMath::Pi()/180);
291  X2 = zz*sin(10*TMath::Pi()/180)+xx*cos(10*TMath::Pi()/180);
292  Y2 = Y1;
293  ntupleR->Fill(Z2,Y2,X2);
294 }
295 
296 c1->cd(10);
297  ntupleR->Draw("X2:Z2","abs(X2)<50","surf3");
298  gPad->SetLogz();
299 
300 //****************
301 // ENERGY DEPOSITS
302 //****************
303 
304 gStyle->SetOptStat(0000);
305 gStyle->SetOptFit();
306 gStyle->SetPalette(1);
307 gROOT->SetStyle("Plain");
308 
309 c1->cd(11);
310  TH2F *histbis = new TH2F("histbis","histbis",50,-20,20,50,-20,20);
311  ntuple4->Draw("y*0.359060:x*0.359060>>histbis","doseV","contz");
312  gPad->SetLogz();
313  histbis->Draw("contz");
314  histbis->GetXaxis()->SetLabelSize(0.025);
315  histbis->GetYaxis()->SetLabelSize(0.025);
316  histbis->GetZaxis()->SetLabelSize(0.025);
317  histbis->GetXaxis()->SetTitleSize(0.035);
318  histbis->GetYaxis()->SetTitleSize(0.035);
319  histbis->GetXaxis()->SetTitleOffset(1.4);
320  histbis->GetYaxis()->SetTitleOffset(1.4);
321  histbis->GetXaxis()->SetTitle("Y (um)");
322  histbis->GetYaxis()->SetTitle("X (um)");
323  histbis->SetTitle("Energy deposit -transverse- (z axis in eV)");
324 
325 c1->cd(12);
326  TH2F *histter = new TH2F("histter","histter",50,-20,20,50,-20,20);
327  ntuple4->Draw("x*0.359060:(z+1500/0.162810+21)*0.162810>>histter","doseV","contz");
328  gPad->SetLogz();
329  histter->Draw("contz");
330  histter->GetXaxis()->SetLabelSize(0.025);
331  histter->GetYaxis()->SetLabelSize(0.025);
332  histter->GetZaxis()->SetLabelSize(0.025);
333  histter->GetXaxis()->SetTitleSize(0.035);
334  histter->GetYaxis()->SetTitleSize(0.035);
335  histter->GetXaxis()->SetTitleOffset(1.4);
336  histter->GetYaxis()->SetTitleOffset(1.4);
337  histter->GetXaxis()->SetTitle("Z (um)");
338  histter->GetYaxis()->SetTitle("X (um)");
339  histter->SetTitle("Energy deposit -longitudinal- (z axis in eV)");
340 
341 //*******************************
342 // BEAM POSITION AT CELL ENTRANCE
343 //*******************************
344 
345 gStyle->SetOptStat(0000);
346 gStyle->SetOptFit();
347 gStyle->SetPalette(1);
348 gROOT->SetStyle("Plain");
349 
350 TH1F *h77 = new TH1F("hx","h1",200,-10,10);
351 TH1F *h88 = new TH1F("hy","h1",200,-10,10);
352 
353 c1->cd(4);
354  ntuple1->Project("hx","x");
355  scale = 1/h77->Integral();
356  h77->Scale(scale);
357  h77->Draw();
358  h77->GetXaxis()->SetLabelSize(0.025);
359  h77->GetYaxis()->SetLabelSize(0.025);
360  h77->GetXaxis()->SetTitleSize(0.035);
361  h77->GetYaxis()->SetTitleSize(0.035);
362  h77->GetXaxis()->SetTitleOffset(1.4);
363  h77->GetYaxis()->SetTitleOffset(1.4);
364  h77->GetXaxis()->SetTitle("Position (um)");
365  h77->GetYaxis()->SetTitle("Fraction of events");
366  h77->SetTitle("Beam X position on cell");
367  h77->SetFillColor(4);
368  h77->SetLineColor(4);
369  h77->Fit("gaus");
370 
371 c1->cd(8);
372  ntuple1->Project("hy","y");
373  scale = 1/h88->Integral();
374  h88->Scale(scale);
375  h88->Draw();
376  h88->GetXaxis()->SetLabelSize(0.025);
377  h88->GetYaxis()->SetLabelSize(0.025);
378  h88->GetXaxis()->SetTitleSize(0.035);
379  h88->GetYaxis()->SetTitleSize(0.035);
380  h88->GetXaxis()->SetTitleOffset(1.4);
381  h88->GetYaxis()->SetTitleOffset(1.4);
382  h88->GetXaxis()->SetTitle("Position (um)");
383  h88->GetYaxis()->SetTitle("Fraction of events");
384  h88->SetTitle("Beam Y position on cell");
385  h88->SetFillColor(4);
386  h88->SetLineColor(4);
387  h88->Fit("gaus");
388 }
Float_t x
Definition: compare.C:6
G4double Double_t
TNtuple * ntuple4
Float_t yVox
Double_t xx
system("rm -rf microbeam.root")
TH2F * histbis
Float_t vox
Float_t dose
TNtuple * ntupleYX
Float_t den
Float_t y
Definition: compare.C:6
Double_t z
Float_t Y
Float_t xVox
Double_t zz
Float_t tmp
TH1F * h10
Float_t zVox
Double_t X2
TNtuple * ntuple2
Float_t Z
nlines
TNtuple * ntupleYXN
G4float Float_t
fclose(fg1)
Float_t voxelSizeY
Double_t scale
TNtuple * ntupleZX
G4int Int_t
TNtuple * ntuple0
Float_t voxelSizeZ
Double_t X1
Double_t Y1
TNtuple * ntuple1
TH2F * histter
Float_t d
TNtuple * ntuple3
TH1F * h2bis
Float_t mat
FILE * fp
TH1F * h77
TH1F * h2
Double_t Z2
Float_t voxelSizeX
Double_t Z1
Int_t nentries
Definition: comparison.C:29
TH1F * h11
TH2F * hist
Double_t Zc
TH1F * h20
TH1F * h1
TH1F * h1bis
Float_t X
Double_t Xc
Int_t ncols
TH1F * h88
TNtuple * ntupleR
Double_t Y2