Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HadrontherapyAnalysisManager.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // This is the *BASIC* version of Hadrontherapy, a Geant4-based application
27 // See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28 //
29 // Visit the Hadrontherapy web site (http://www.lns.infn.it/link/Hadrontherapy) to request
30 // the *COMPLETE* version of this program, together with its documentation;
31 // Hadrontherapy (both basic and full version) are supported by the Italian INFN
32 // Institute in the framework of the MC-INFN Group
33 //
34 
36 #include "HadrontherapyMatrix.hh"
38 #include "G4SystemOfUnits.hh"
39 #include <time.h>
40 
42 
43 HadrontherapyAnalysisManager::HadrontherapyAnalysisManager()
44 #ifdef G4ANALYSIS_USE_ROOT
45 :
46 analysisFileName("DoseDistribution.root"),theTFile(0), histo1(0), histo2(0), histo3(0),
47 histo4(0), histo5(0), histo6(0), histo7(0), histo8(0), histo9(0), histo10(0), histo11(0), histo12(0), histo13(0), histo14(0), histo15(0), histo16(0),
48 kinFragNtuple(0),
49 kineticEnergyPrimaryNtuple(0),
50 doseFragNtuple(0),
51 fluenceFragNtuple(0),
52 letFragNtuple(0),
53 theROOTNtuple(0),
54 theROOTIonTuple(0),
55 fragmentNtuple(0),
56 metaData(0),
57 eventCounter(0)
58 #endif
59 {
60  fMess = new HadrontherapyAnalysisFileMessenger(this);
61 }
63 
65 {
66  delete fMess;
67 #ifdef G4ANALYSIS_USE_ROOT
68  Clear();
69 #endif
70 }
71 
73 {
75  return instance;
76 }
77 #ifdef G4ANALYSIS_USE_ROOT
79 {
80  if (theTFile)
81  {
82  delete metaData;
83  metaData = 0;
84 
85  delete fragmentNtuple;
86  fragmentNtuple = 0;
87 
88  delete theROOTIonTuple;
89  theROOTIonTuple = 0;
90 
91  delete theROOTNtuple;
92  theROOTNtuple = 0;
93 
94  delete histo16;
95  histo16 = 0;
96 
97  delete histo15;
98  histo15 = 0;
99 
100  delete histo14;
101  histo14 = 0;
102 
103  delete histo13;
104  histo13 = 0;
105 
106  delete histo12;
107  histo12 = 0;
108 
109  delete histo11;
110  histo11 = 0;
111 
112  delete histo10;
113  histo10 = 0;
114 
115  delete histo9;
116  histo9 = 0;
117 
118  delete histo8;
119  histo8 = 0;
120 
121  delete histo7;
122  histo7 = 0;
123 
124  delete histo6;
125  histo6 = 0;
126 
127  delete histo5;
128  histo5 = 0;
129 
130  delete histo4;
131  histo4 = 0;
132 
133  delete histo3;
134  histo3 = 0;
135 
136  delete histo2;
137  histo2 = 0;
138 
139  delete histo1;
140  histo1 = 0;
141  }
142 }
144 
145 void HadrontherapyAnalysisManager::SetAnalysisFileName(G4String aFileName)
146 {
147  this->analysisFileName = aFileName;
148 }
149 
151 G4bool HadrontherapyAnalysisManager::IsTheTFile()
152 {
153  return (theTFile) ? true:false;
154 }
155 void HadrontherapyAnalysisManager::book()
156 {
157  delete theTFile; // this is similar to theTFile->Close() => delete all associated variables created via new, moreover it delete itself.
158 
159  theTFile = new TFile(analysisFileName, "RECREATE");
160 
161  // Create the histograms with the energy deposit along the X axis
162  histo1 = createHistogram1D("braggPeak","slice, energy", 400, 0., 80); //<different waterthicknesses are accoutned for in ROOT-analysis stage
163  histo2 = createHistogram1D("h20","Secondary protons - slice, energy", 400, 0., 400.);
164  histo3 = createHistogram1D("h30","Secondary neutrons - slice, energy", 400, 0., 400.);
165  histo4 = createHistogram1D("h40","Secondary alpha - slice, energy", 400, 0., 400.);
166  histo5 = createHistogram1D("h50","Secondary gamma - slice, energy", 400, 0., 400.);
167  histo6 = createHistogram1D("h60","Secondary electron - slice, energy", 400, 0., 400.);
168  histo7 = createHistogram1D("h70","Secondary triton - slice, energy", 400, 0., 400.);
169  histo8 = createHistogram1D("h80","Secondary deuteron - slice, energy", 400, 0., 400.);
170  histo9 = createHistogram1D("h90","Secondary pion - slice, energy", 400, 0., 400.);
171  histo10 = createHistogram1D("h100","Energy distribution of secondary electrons", 70, 0., 70.);
172  histo11 = createHistogram1D("h110","Energy distribution of secondary photons", 70, 0., 70.);
173  histo12 = createHistogram1D("h120","Energy distribution of secondary deuterons", 70, 0., 70.);
174  histo13 = createHistogram1D("h130","Energy distribution of secondary tritons", 70, 0., 70.);
175  histo14 = createHistogram1D("h140","Energy distribution of secondary alpha particles", 70, 0., 70.);
176  histo15 = createHistogram1D("heliumEnergyAfterPhantom","Energy distribution of secondary helium fragments after the phantom",
177  70, 0., 500.);
178  histo16 = createHistogram1D("hydrogenEnergyAfterPhantom","Energy distribution of secondary helium fragments after the phantom",
179  70, 0., 500.);
180 
181  kinFragNtuple = new TNtuple("kinFragNtuple",
182  "Kinetic energy by voxel & fragment",
183  "i:j:k:A:Z:kineticEnergy");
184  kineticEnergyPrimaryNtuple= new TNtuple("kineticEnergyPrimaryNtuple",
185  "Kinetic energy by voxel of primary",
186  "i:j:k:kineticEnergy");
187  doseFragNtuple = new TNtuple("doseFragNtuple",
188  "Energy deposit by voxel & fragment",
189  "i:j:k:A:Z:energy");
190 
191  fluenceFragNtuple = new TNtuple("fluenceFragNtuple",
192  "Fluence by voxel & fragment",
193  "i:j:k:A:Z:fluence");
194 
195  letFragNtuple = new TNtuple("letFragNtuple",
196  "Let by voxel & fragment",
197  "i:j:k:A:Z:letT:letD");
198 
199  theROOTNtuple = new TNtuple("theROOTNtuple",
200  "Energy deposit by slice",
201  "i:j:k:energy");
202 
203  theROOTIonTuple = new TNtuple("theROOTIonTuple",
204  "Generic ion information",
205  "a:z:occupancy:energy");
206 
207  fragmentNtuple = new TNtuple("fragmentNtuple",
208  "Fragments",
209  "A:Z:energy:posX:posY:posZ");
210 
211  metaData = new TNtuple("metaData",
212  "Metadata",
213  "events:detectorDistance:waterThickness:beamEnergy:energyError:phantomCenterDistance");
214 }
215 
217 void HadrontherapyAnalysisManager::FillEnergyDeposit(G4int i,
218  G4int j,
219  G4int k,
221 {
222  if (theROOTNtuple)
223  {
224  theROOTNtuple->Fill(i, j, k, energy);
225  }
226 }
227 
229 void HadrontherapyAnalysisManager::BraggPeak(G4int slice, G4double energy)
230 {
231  histo1->SetBinContent(slice, energy); //This uses setbincontent instead of fill to get labels correct
232 }
233 
235 void HadrontherapyAnalysisManager::SecondaryProtonEnergyDeposit(G4int slice, G4double energy)
236 {
237  histo2->Fill(slice, energy);
238 }
239 
241 void HadrontherapyAnalysisManager::SecondaryNeutronEnergyDeposit(G4int slice, G4double energy)
242 {
243  histo3->Fill(slice, energy);
244 }
245 
247 void HadrontherapyAnalysisManager::SecondaryAlphaEnergyDeposit(G4int slice, G4double energy)
248 {
249  histo4->Fill(slice, energy);
250 }
251 
253 void HadrontherapyAnalysisManager::SecondaryGammaEnergyDeposit(G4int slice, G4double energy)
254 {
255  histo5->Fill(slice, energy);
256 }
257 
259 void HadrontherapyAnalysisManager::SecondaryElectronEnergyDeposit(G4int slice, G4double energy)
260 {
261  histo6->Fill(slice, energy);
262 }
263 
265 void HadrontherapyAnalysisManager::SecondaryTritonEnergyDeposit(G4int slice, G4double energy)
266 {
267  histo7->Fill(slice, energy);
268 }
269 
271 void HadrontherapyAnalysisManager::SecondaryDeuteronEnergyDeposit(G4int slice, G4double energy)
272 {
273  histo8->Fill(slice, energy);
274 }
275 
277 void HadrontherapyAnalysisManager::SecondaryPionEnergyDeposit(G4int slice, G4double energy)
278 {
279  histo9->Fill(slice, energy);
280 }
281 
283 void HadrontherapyAnalysisManager::electronEnergyDistribution(G4double energy)
284 {
285  histo10->Fill(energy);
286 }
287 
289 void HadrontherapyAnalysisManager::gammaEnergyDistribution(G4double energy)
290 {
291  histo11->Fill(energy);
292 }
293 
295 void HadrontherapyAnalysisManager::deuteronEnergyDistribution(G4double energy)
296 {
297  histo12->Fill(energy);
298 }
299 
301 void HadrontherapyAnalysisManager::tritonEnergyDistribution(G4double energy)
302 {
303  histo13->Fill(energy);
304 }
305 
307 void HadrontherapyAnalysisManager::alphaEnergyDistribution(G4double energy)
308 {
309  histo14->Fill(energy);
310 }
312 void HadrontherapyAnalysisManager::heliumEnergy(G4double secondaryParticleKineticEnergy)
313 {
314  histo15->Fill(secondaryParticleKineticEnergy);
315 }
316 
318 void HadrontherapyAnalysisManager::hydrogenEnergy(G4double secondaryParticleKineticEnergy)
319 {
320  histo16->Fill(secondaryParticleKineticEnergy);
321 }
322 
324  // FillKineticFragmentTuple create an ntuple where the voxel indexs, the atomic number and mass and the kinetic
325  // energy of all the particles interacting with the phantom, are stored
326 void HadrontherapyAnalysisManager::FillKineticFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double kinEnergy)
327 {
328  kinFragNtuple -> Fill(i, j, k, A, Z, kinEnergy);
329 }
330 
332  // FillKineticEnergyPrimaryNTuple creates a ntuple where the voxel indexs and the kinetic
333  // energies of ONLY primary particles interacting with the phantom, are stored
334 void HadrontherapyAnalysisManager::FillKineticEnergyPrimaryNTuple(G4int i, G4int j, G4int k, G4double kinEnergy)
335 {
336  kineticEnergyPrimaryNtuple -> Fill(i, j, k, kinEnergy);
337 }
338 
340  // This function is called only if ROOT is activated.
341  // It is called by the HadrontherapyMatric.cc class file and it is used to create two ntuples containing
342  // the total energy deposited and the fluence values, in each voxel and per any particle (primary
343  // and secondary particles beam)
344 void HadrontherapyAnalysisManager::FillVoxelFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double energy, G4double fluence)
345 {
346  // Fill the ntuple containing the voxel, mass and atomic number and the energy deposited
347  doseFragNtuple -> Fill( i, j, k, A, Z, energy );
348 
349  // Fill the ntuple containing the voxel, mass and atomic number and the fluence
350  if (i==1 && Z==1) {
351  fluenceFragNtuple -> Fill( i, j, k, A, Z, fluence );
352 
353  }
354 }
355 
356 void HadrontherapyAnalysisManager::FillLetFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double letT, G4double letD)
357 {
358  letFragNtuple -> Fill( i, j, k, A, Z, letT, letD);
359 
360 }
362 void HadrontherapyAnalysisManager::FillFragmentTuple(G4int A, G4double Z, G4double energy, G4double posX, G4double posY, G4double posZ)
363 {
364  fragmentNtuple->Fill(A, Z, energy, posX, posY, posZ);
365 }
366 
368 void HadrontherapyAnalysisManager::genericIonInformation(G4int a,
369  G4double z,
370  G4int electronOccupancy,
371  G4double energy)
372 {
373  if (theROOTIonTuple) {
374  theROOTIonTuple->Fill(a, z, electronOccupancy, energy);
375  }
376 }
377 
379 void HadrontherapyAnalysisManager::startNewEvent()
380 {
381  eventCounter++;
382 }
384 void HadrontherapyAnalysisManager::setGeometryMetaData(G4double endDetectorPosition, G4double waterThickness, G4double phantomCenter)
385 {
386  this->detectorDistance = endDetectorPosition;
387  this->phantomDepth = waterThickness;
388  this->phantomCenterDistance = phantomCenter;
389 }
390 void HadrontherapyAnalysisManager::setBeamMetaData(G4double meanKineticEnergy,G4double sigmaEnergy)
391 {
392  this->beamEnergy = meanKineticEnergy;
393  this->energyError = sigmaEnergy;
394 }
396 // Flush data & close the file
397 void HadrontherapyAnalysisManager::flush()
398 {
399  if (theTFile)
400  {
401  theTFile -> Write();
402  theTFile -> Close();
403  }
404  theTFile = 0;
405  eventCounter = 0;
406 }
407 
408 #endif