source: trunk/examples/advanced/hadrontherapy/include/HadrontherapyAnalysisManager.hh@ 1330

Last change on this file since 1330 was 1313, checked in by garnier, 15 years ago

geant4.9.4 beta rc0

File size: 7.8 KB
Line 
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 // HadrontherapyAnalysisManager.hh; May 2005
27 // See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28
29#ifndef HADRONTHERAPYANALYSISMANAGER_HH
30#define HADRONTHERAPYANALYSISMANAGER_HH 1
31
32#include "globals.hh"
33
34
35#ifdef G4ANALYSIS_USE_ROOT ///< If analysis is done directly with ROOT
36#include "TROOT.h"
37#include "TFile.h"
38#include "TNtuple.h"
39#include "TH1F.h"
40#endif
41/**
42 * Messenger class for analysis-settings for HadronTherapyAnalysisManager
43 */
44class HadrontherapyAnalysisFileMessenger;
45
46/**
47 * A class for connecting the simulation to an analysis package.
48 */
49class HadrontherapyAnalysisManager
50{
51private:
52 /**
53 * Analysis manager is a singleton object (there is only one instance).
54 * The pointer to this object is available through the use of the method GetInstance();
55 *
56 * @see GetInstance
57 */
58 HadrontherapyAnalysisManager();
59
60public:
61 ~HadrontherapyAnalysisManager();
62
63 /**
64 * Get the pointer to the analysis manager.
65 */
66 static HadrontherapyAnalysisManager* GetInstance();
67
68#ifdef G4ANALYSIS_USE_ROOT
69 /**
70 * Clear analysis manager heap.
71 */
72 void Clear();
73
74 /**
75 * Book the histograms and ntuples in an AIDA or ROOT file.
76 */
77 void book();
78 /**
79 * Set name for the analysis file .root (used by macro)
80 */
81 void SetAnalysisFileName(G4String);
82
83 /**
84 * Fill the ntuple with the energy deposit in the phantom
85 */
86 void FillEnergyDeposit(G4int voxelXId, G4int voxelYId, G4int voxelZId,
87 G4double energyDeposit);
88
89 void BraggPeak(G4int, G4double); ///< Fill 1D histogram with the Bragg peak in the phantom
90
91 void SecondaryProtonEnergyDeposit(G4int slice, G4double energy);
92 ///< Fill 1D histogram with the energy deposit of secondary protons
93
94 void SecondaryNeutronEnergyDeposit(G4int slice, G4double energy);
95 ///< Fill 1D histogram with the energy deposit of secondary neutrons
96
97 void SecondaryAlphaEnergyDeposit(G4int slice, G4double energy);
98 ///< Fill 1D histogram with the energy deposit of secondary alpha particles
99
100 void SecondaryGammaEnergyDeposit(G4int slice, G4double energy);
101 ///< Fill 1D histogram with the energy deposit of secondary gamma
102
103 void SecondaryElectronEnergyDeposit(G4int slice, G4double energy);
104 ///< Fill 1D histogram with the energy deposit of secondary electrons
105
106 void SecondaryTritonEnergyDeposit(G4int slice, G4double energy);
107 ///< Fill 1D histogram with the energy deposit of secondary tritons
108
109 void SecondaryDeuteronEnergyDeposit(G4int slice, G4double energy);
110 ///< Fill 1D histogram with the energy deposit of secondary deuterons
111
112 void SecondaryPionEnergyDeposit(G4int slice, G4double energy);
113 ///< Fill 1D histogram with the energy deposit of secondary pions
114
115 void electronEnergyDistribution(G4double secondaryParticleKineticEnergy);
116 ///< Energy distribution of secondary electrons originated in the phantom
117
118 void gammaEnergyDistribution(G4double secondaryParticleKineticEnergy);
119 ///< Energy distribution of secondary gamma originated in the phantom
120
121 void deuteronEnergyDistribution(G4double secondaryParticleKineticEnergy);
122 ///< Energy distribution of secondary deuterons originated in the phantom
123
124 void tritonEnergyDistribution(G4double secondaryParticleKineticEnergy);
125 ///< Energy distribution of secondary tritons originated in the phantom
126
127 void alphaEnergyDistribution(G4double secondaryParticleKineticEnergy);
128 ///< Energy distribution of secondary alpha originated in the phantom
129
130 void heliumEnergy(G4double secondaryParticleKineticEnergy);
131 ///< Energy distribution of the helium (He3 and alpha) particles after the phantom
132
133 void hydrogenEnergy(G4double secondaryParticleKineticEnergy);
134 ///< Energy distribution of the hydrogen (proton, d, t) particles after the phantom
135
136 //Kinetic energy by voxel, mass number A and atomic number Z.
137 void FillKineticFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double kinEnergy);
138
139 //Kinetic energy by voxel, mass number A and atomic number Z of only primary particles
140 void FillKineticEnergyPrimaryNTuple(G4int i, G4int j, G4int k, G4double kinEnergy);
141
142 ///< Energy by voxel, mass number A and atomic number Z.
143 void FillVoxelFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double energy, G4double fluence);
144
145 void FillFragmentTuple(G4int A, G4double Z, G4double energy, G4double posX, G4double posY, G4double posZ);
146 ///< Energy ntuple
147
148 void FillLetFragmentTuple(G4int i, G4int j, G4int k, G4int A, G4double Z, G4double letT, G4double letD);
149 ///< let ntuple
150
151 void genericIonInformation(G4int, G4double, G4int, G4double);
152
153 void ThintargetBeamDisp(G4double,G4double);
154
155 void startNewEvent();
156 ///< Tell the analysis manager that a new event is starting
157
158 void setGeometryMetaData(G4double, G4double, G4double);
159 ///< from the detector construction information about the geometry can be written as metadata
160
161 void setBeamMetaData(G4double, G4double);
162 ///< metadata about the beam can be written this way
163
164 void flush();
165 ///< Close the .hbk file with the histograms and the ntuples
166private:
167 TH1F *createHistogram1D(const TString name, const TString title, int bins, double xmin, double xmax) {
168 TH1F *histo = new TH1F(name, title, bins, xmin, xmax);
169 histo->SetLineWidth(2);
170 return histo;
171 }
172
173private:
174#endif
175 static HadrontherapyAnalysisManager* instance;
176 HadrontherapyAnalysisFileMessenger* fMess;
177#ifdef G4ANALYSIS_USE_ROOT
178 G4String analysisFileName;
179 TFile *theTFile;
180 TH1F *histo1;
181 TH1F *histo2;
182 TH1F *histo3;
183 TH1F *histo4;
184 TH1F *histo5;
185 TH1F *histo6;
186 TH1F *histo7;
187 TH1F *histo8;
188 TH1F *histo9;
189 TH1F *histo10;
190 TH1F *histo11;
191 TH1F *histo12;
192 TH1F *histo13;
193 TH1F *histo14;
194 TH1F *histo15;
195 TH1F *histo16;
196
197 TNtuple *kinFragNtuple;
198 TNtuple *kineticEnergyPrimaryNtuple;
199
200 // ntuple containing the fluence of all the particle in any voxel
201 TNtuple *doseFragNtuple;
202
203 // ntuple containing the fluence of all the particle in any voxel
204 TNtuple *fluenceFragNtuple;
205
206 // ntuple containing the fluence of all the particle in any voxel
207 TNtuple *letFragNtuple;
208
209 TNtuple *theROOTNtuple;
210 TNtuple *theROOTIonTuple;
211 TNtuple *fragmentNtuple; // fragments
212 TNtuple *metaData;
213 G4long eventCounter; // Simulation metadata
214 G4double detectorDistance;
215 G4double phantomDepth;
216 G4double beamEnergy;
217 G4double energyError;
218 G4double phantomCenterDistance;
219#endif
220};
221#endif
222
223
224
Note: See TracBrowser for help on using the repository browser.