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

Last change on this file since 1312 was 1230, checked in by garnier, 16 years ago

update to geant4.9.3

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