Geant4  10.01.p03
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 // Visit the Hadrontherapy web site (http://www.lns.infn.it/link/Hadrontherapy) to request
26 // the *COMPLETE* version of this program, together with its documentation;
27 // Hadrontherapy (both basic and full version) are supported by the Italian INFN
28 // Institute in the framework of the MC-INFN Group
29 //
30 
32 #include "HadrontherapyMatrix.hh"
34 #include "G4SystemOfUnits.hh"
35 #include <time.h>
36 
38 
40 #ifdef G4ANALYSIS_USE_ROOT
41 :
42 analysisFileName("DoseDistribution.root"),theTFile(0), histo1(0), histo2(0), histo3(0),
43 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),
44 kinFragNtuple(0),
45 kineticEnergyPrimaryNtuple(0),
46 doseFragNtuple(0),
47 fluenceFragNtuple(0),
48 letFragNtuple(0),
49 theROOTNtuple(0),
50 theROOTIonTuple(0),
51 fragmentNtuple(0),
52 metaData(0),
53 eventCounter(0)
54 #endif
55 {
56  fMess = new HadrontherapyAnalysisFileMessenger(this);
57 }
59 
61 {
62  delete fMess;
63 #ifdef G4ANALYSIS_USE_ROOT
64  Clear();
65 #endif
66 }
67 
69 {
71  return instance;
72 }
73 #ifdef G4ANALYSIS_USE_ROOT
75 {
76  if (theTFile)
77  {
78  delete metaData;
79  metaData = 0;
80 
81  delete fragmentNtuple;
82  fragmentNtuple = 0;
83 
84  delete theROOTIonTuple;
85  theROOTIonTuple = 0;
86 
87  delete theROOTNtuple;
88  theROOTNtuple = 0;
89 
90 
91  delete histo16;
92  histo16 = 0;
93 
94  delete histo15;
95  histo15 = 0;
96 
97  delete histo14;
98  histo14 = 0;
99 
100  delete histo13;
101  histo13 = 0;
102 
103  delete histo12;
104  histo12 = 0;
105 
106  delete histo11;
107  histo11 = 0;
108 
109  delete histo10;
110  histo10 = 0;
111 
112  delete histo9;
113  histo9 = 0;
114 
115  delete histo8;
116  histo8 = 0;
117 
118  delete histo7;
119  histo7 = 0;
120 
121  delete histo6;
122  histo6 = 0;
123 
124  delete histo5;
125  histo5 = 0;
126 
127  delete histo4;
128  histo4 = 0;
129 
130  delete histo3;
131  histo3 = 0;
132 
133  delete histo2;
134  histo2 = 0;
135 
136  delete histo1;
137  histo1 = 0;
138  }
139 }
141 
142 void HadrontherapyAnalysisManager::SetAnalysisFileName(G4String aFileName)
143 {
144  this->analysisFileName = aFileName;
145 }
146 
148 G4bool HadrontherapyAnalysisManager::IsTheTFile()
149 {
150  return (theTFile) ? true:false;
151 }
152 void HadrontherapyAnalysisManager::book()
153 {
154  delete theTFile; // this is similar to theTFile->Close() => delete all associated variables created via new, moreover it delete itself.
155 
156  theTFile = new TFile(analysisFileName, "RECREATE");
157 
158  // Create the histograms with the energy deposit along the X axis
159  histo1 = createHistogram1D("braggPeak","slice, energy", 400, 0., 80); //<different waterthicknesses are accoutned for in ROOT-analysis stage
160  histo2 = createHistogram1D("h20","Secondary protons - slice, energy", 400, 0., 400.);
161  histo3 = createHistogram1D("h30","Secondary neutrons - slice, energy", 400, 0., 400.);
162  histo4 = createHistogram1D("h40","Secondary alpha - slice, energy", 400, 0., 400.);
163  histo5 = createHistogram1D("h50","Secondary gamma - slice, energy", 400, 0., 400.);
164  histo6 = createHistogram1D("h60","Secondary electron - slice, energy", 400, 0., 400.);
165  histo7 = createHistogram1D("h70","Secondary triton - slice, energy", 400, 0., 400.);
166  histo8 = createHistogram1D("h80","Secondary deuteron - slice, energy", 400, 0., 400.);
167  histo9 = createHistogram1D("h90","Secondary pion - slice, energy", 400, 0., 400.);
168  histo10 = createHistogram1D("h100","Energy distribution of secondary electrons", 70, 0., 70.);
169  histo11 = createHistogram1D("h110","Energy distribution of secondary photons", 70, 0., 70.);
170  histo12 = createHistogram1D("h120","Energy distribution of secondary deuterons", 70, 0., 70.);
171  histo13 = createHistogram1D("h130","Energy distribution of secondary tritons", 70, 0., 70.);
172  histo14 = createHistogram1D("h140","Energy distribution of secondary alpha particles", 70, 0., 70.);
173  histo15 = createHistogram1D("heliumEnergyAfterPhantom","Energy distribution of secondary helium fragments after the phantom",
174  70, 0., 500.);
175  histo16 = createHistogram1D("hydrogenEnergyAfterPhantom","Energy distribution of secondary helium fragments after the phantom",
176  70, 0., 500.);
177 
178  kinFragNtuple = new TNtuple("kinFragNtuple",
179  "Kinetic energy by voxel & fragment",
180  "i:j:k:A:Z:kineticEnergy");
181  kineticEnergyPrimaryNtuple= new TNtuple("kineticEnergyPrimaryNtuple",
182  "Kinetic energy by voxel of primary",
183  "i:j:k:kineticEnergy");
184  doseFragNtuple = new TNtuple("doseFragNtuple",
185  "Energy deposit by voxel & fragment",
186  "i:j:k:A:Z:energy");
187 
188  fluenceFragNtuple = new TNtuple("fluenceFragNtuple",
189  "Fluence by voxel & fragment",
190  "i:j:k:A:Z:fluence");
191 
192  letFragNtuple = new TNtuple("letFragNtuple",
193  "Let by voxel & fragment",
194  "i:j:k:A:Z:letT:letD");
195 
196  theROOTNtuple = new TNtuple("theROOTNtuple",
197  "Energy deposit by slice",
198  "i:j:k:energy");
199 
200  theROOTIonTuple = new TNtuple("theROOTIonTuple",
201  "Generic ion information",
202  "a:z:occupancy:energy");
203 
204  fragmentNtuple = new TNtuple("fragmentNtuple",
205  "Fragments",
206  "A:Z:energy:posX:posY:posZ");
207 
208  metaData = new TNtuple("metaData",
209  "Metadata",
210  "events:detectorDistance:waterThickness:beamEnergy:energyError:phantomCenterDistance");
211 }
212 
214 void HadrontherapyAnalysisManager::FillEnergyDeposit(G4int i,
215  G4int j,
216  G4int k,
218 {
219  if (theROOTNtuple)
220  {
221  theROOTNtuple->Fill(i, j, k, energy);
222  }
223 }
224 
226 void HadrontherapyAnalysisManager::BraggPeak(G4int slice, G4double energy)
227 {
228  histo1->SetBinContent(slice, energy); //This uses setbincontent instead of fill to get labels correct
229 }
230 
232 void HadrontherapyAnalysisManager::SecondaryProtonEnergyDeposit(G4int slice, G4double energy)
233 {
234  histo2->Fill(slice, energy);
235 }
236 
238 void HadrontherapyAnalysisManager::SecondaryNeutronEnergyDeposit(G4int slice, G4double energy)
239 {
240  histo3->Fill(slice, energy);
241 }
242 
244 void HadrontherapyAnalysisManager::SecondaryAlphaEnergyDeposit(G4int slice, G4double energy)
245 {
246  histo4->Fill(slice, energy);
247 }
248 
250 void HadrontherapyAnalysisManager::SecondaryGammaEnergyDeposit(G4int slice, G4double energy)
251 {
252  histo5->Fill(slice, energy);
253 }
254 
256 void HadrontherapyAnalysisManager::SecondaryElectronEnergyDeposit(G4int slice, G4double energy)
257 {
258  histo6->Fill(slice, energy);
259 }
260 
262 void HadrontherapyAnalysisManager::SecondaryTritonEnergyDeposit(G4int slice, G4double energy)
263 {
264  histo7->Fill(slice, energy);
265 }
266 
268 void HadrontherapyAnalysisManager::SecondaryDeuteronEnergyDeposit(G4int slice, G4double energy)
269 {
270  histo8->Fill(slice, energy);
271 }
272 
274 void HadrontherapyAnalysisManager::SecondaryPionEnergyDeposit(G4int slice, G4double energy)
275 {
276  histo9->Fill(slice, energy);
277 }
278 
280 void HadrontherapyAnalysisManager::electronEnergyDistribution(G4double energy)
281 {
282  histo10->Fill(energy);
283 }
284 
286 void HadrontherapyAnalysisManager::gammaEnergyDistribution(G4double energy)
287 {
288  histo11->Fill(energy);
289 }
290 
292 void HadrontherapyAnalysisManager::deuteronEnergyDistribution(G4double energy)
293 {
294  histo12->Fill(energy);
295 }
296 
298 void HadrontherapyAnalysisManager::tritonEnergyDistribution(G4double energy)
299 {
300  histo13->Fill(energy);
301 }
302 
304 void HadrontherapyAnalysisManager::alphaEnergyDistribution(G4double energy)
305 {
306  histo14->Fill(energy);
307 }
309 void HadrontherapyAnalysisManager::heliumEnergy(G4double secondaryParticleKineticEnergy)
310 {
311  histo15->Fill(secondaryParticleKineticEnergy);
312 }
313 
315 void HadrontherapyAnalysisManager::hydrogenEnergy(G4double secondaryParticleKineticEnergy)
316 {
317  histo16->Fill(secondaryParticleKineticEnergy);
318 }
319 
321  // FillKineticFragmentTuple create an ntuple where the voxel indexs, the atomic number and mass and the kinetic
322  // energy of all the particles interacting with the phantom, are stored
323 void HadrontherapyAnalysisManager::FillKineticFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double kinEnergy)
324 {
325  kinFragNtuple -> Fill(i, j, k, A, Z, kinEnergy);
326 }
327 
329  // FillKineticEnergyPrimaryNTuple creates a ntuple where the voxel indexs and the kinetic
330  // energies of ONLY primary particles interacting with the phantom, are stored
331 void HadrontherapyAnalysisManager::FillKineticEnergyPrimaryNTuple(G4int i, G4int j, G4int k, G4double kinEnergy)
332 {
333  kineticEnergyPrimaryNtuple -> Fill(i, j, k, kinEnergy);
334 }
335 
337  // This function is called only if ROOT is activated.
338  // It is called by the HadrontherapyMatric.cc class file and it is used to create two ntuples containing
339  // the total energy deposited and the fluence values, in each voxel and per any particle (primary
340  // and secondary particles beam)
341 void HadrontherapyAnalysisManager::FillVoxelFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double energy, G4double fluence)
342 {
343  // Fill the ntuple containing the voxel, mass and atomic number and the energy deposited
344  doseFragNtuple -> Fill( i, j, k, A, Z, energy );
345 
346  // Fill the ntuple containing the voxel, mass and atomic number and the fluence
347  if (i==1 && Z==1) {
348  fluenceFragNtuple -> Fill( i, j, k, A, Z, fluence );
349 
350  }
351 }
352 
353 void HadrontherapyAnalysisManager::FillLetFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double letD)
354 {
355  letFragNtuple -> Fill( i, j, k, A, Z, letD/MeV);
356 }
358 void HadrontherapyAnalysisManager::FillFragmentTuple(G4int A, G4double Z, G4double energy, G4double posX, G4double posY, G4double posZ)
359 {
360  fragmentNtuple->Fill(A, Z, energy, posX, posY, posZ);
361 }
362 
364 void HadrontherapyAnalysisManager::genericIonInformation(G4int a,
365  G4double z,
366  G4int electronOccupancy,
367  G4double energy)
368 {
369  if (theROOTIonTuple) {
370  theROOTIonTuple->Fill(a, z, electronOccupancy, energy);
371  }
372 }
373 
375 void HadrontherapyAnalysisManager::startNewEvent()
376 {
377  eventCounter++;
378 }
380 void HadrontherapyAnalysisManager::setGeometryMetaData(G4double endDetectorPosition, G4double waterThickness, G4double phantomCenter)
381 {
382  this->detectorDistance = endDetectorPosition;
383  this->phantomDepth = waterThickness;
384  this->phantomCenterDistance = phantomCenter;
385 }
386 void HadrontherapyAnalysisManager::setBeamMetaData(G4double meanKineticEnergy,G4double sigmaEnergy)
387 {
388  this->beamEnergy = meanKineticEnergy;
389  this->energyError = sigmaEnergy;
390 }
392 // Flush data & close the file
393 void HadrontherapyAnalysisManager::flush()
394 {
395  if (theTFile)
396  {
397  theTFile -> Write();
398  theTFile -> Close();
399  }
400  theTFile = 0;
401  eventCounter = 0;
402 }
403 
404 #endif
static HadrontherapyAnalysisManager * GetInstance()
Get the pointer to the analysis manager.
static const double MeV
Definition: G4SIunits.hh:193
G4double z
Definition: TRTMaterials.hh:39
void Clear(Node *)
static HadrontherapyAnalysisManager * instance
G4double a
Definition: TRTMaterials.hh:39
A messenger object of this class is created by the AnalysisManager.
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
HadrontherapyAnalysisFileMessenger * fMess
static const G4double A[nN]
G4double energy(const ThreeVector &p, const G4double m)
HadrontherapyAnalysisManager()
Analysis manager is a singleton object (there is only one instance).
A class for connecting the simulation to an analysis package.
double G4double
Definition: G4Types.hh:76