Geant4  10.01
HadrontherapyDetectorSD.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 
31 
34 #include "G4Step.hh"
35 #include "G4VTouchable.hh"
36 #include "G4TouchableHistory.hh"
37 #include "G4SDManager.hh"
38 #include "HadrontherapyMatrix.hh"
39 #include "HadrontherapyLet.hh"
40 #include "G4Track.hh"
42 #include "G4SystemOfUnits.hh"
43 
47 {
48  G4String HCname;
49  collectionName.insert(HCname="HadrontherapyDetectorHitsCollection");
50  HitsCollection = NULL;
52 
53 }
54 
57 {
58 }
59 
62 {
63 
65  collectionName[0]);
66 }
67 
70 {
71 
72 
73  if (aStep -> GetPreStepPoint() -> GetPhysicalVolume() -> GetName() != "RODetectorZDivisionPhys") return false;
74 
75 
76  // Get kinetic energy
77  G4Track * theTrack = aStep -> GetTrack();
78  G4double kineticEnergy = theTrack -> GetKineticEnergy();
79 
80  G4ParticleDefinition *particleDef = theTrack -> GetDefinition();
81  //Get particle name
82  G4String particleName = particleDef -> GetParticleName();
83 
84  G4int pdg = particleDef ->GetPDGEncoding();
85 
86  // Get unique track_id (in an event)
87  G4int trackID = theTrack -> GetTrackID();
88 
89  G4double energyDeposit = aStep -> GetTotalEnergyDeposit();
90 
91  G4double DX = aStep -> GetStepLength();
92  G4int Z = particleDef-> GetAtomicNumber();
93  G4int A = particleDef-> GetAtomicMass();
94 
95 
96  // Read voxel indexes: i is the x index, k is the z index
97  const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
98  G4int k = touchable->GetReplicaNumber(0);
99  G4int i = touchable->GetReplicaNumber(2);
100  G4int j = touchable->GetReplicaNumber(1);
101 
102 
103 
104 #ifdef G4ANALYSIS_USE_ROOT
106 #endif
107 
110 
111 
112  // ******************** let ***************************
113  if (let)
114  {
115  if ( !(Z==0 && A==1) ) // All but not neutrons
116  {
117  if( energyDeposit>0. && DX >0. )
118  {
119  if (pdg !=22) // not gamma
120  {
121  let -> FillEnergySpectrum(trackID, particleDef,energyDeposit, DX, i, j, k);
122  }
123  else if (kineticEnergy > 50.*keV) // gamma cut
124  {
125  let -> FillEnergySpectrum(trackID, particleDef,energyDeposit, DX, i , j, k);
126  }
127 
128  }
129  }
130  }
131 
132 
133 
134  // ******************** let ***************************
135 
136 
137  if (matrix)
138  {
139 
140  // Increment Fluences & accumulate energy spectra
141  // Hit voxels are marked with track_id throught hitTrack matrix
142  G4int* hitTrack = matrix -> GetHitTrack(i,j,k); // hitTrack MUST BE cleared at every eventAction!
143  if ( *hitTrack != trackID )
144  {
145  *hitTrack = trackID;
146  /*
147  * Fill FLUENCE data for every single nuclide
148  * Exclude e-, neutrons, gamma, ...
149  */
150  if ( Z >= 1) matrix -> Fill(trackID, particleDef, i, j, k, 0, true);
151 
152 
153 #ifdef G4ANALYSIS_USE_ROOT
154 /*
155  // Fragments kinetic energy (ntuple)
156  if (trackID !=1 && Z>=1)
157  {
158  // First step kinetic energy for every fragment
159  analysis -> FillKineticFragmentTuple(i, j, k, A, Z, kineticEnergy/MeV);
160  }
161  // Kinetic energy spectra for primary particles
162 
163  if ( trackID == 1 && i == 0)
164  {
165  // First step kinetic energy for primaries only
166  analysis -> FillKineticEnergyPrimaryNTuple(i, j, k, kineticEnergy/MeV);
167  }
168 */
169 #endif
170  }
171 
172  if(energyDeposit != 0)
173  {
174 /*
175  * This method will fill a dose matrix for every single nuclide.
176  * A method of the HadrontherapyMatrix class (StoreDoseFluenceAscii())
177  * is called automatically at the end of main (or via the macro command /analysis/writeDoseFile.
178  * It permits to store all dose/fluence data into a single plane ASCII file.
179 */
180  // if (A==1 && Z==1) // primary and sec. protons
181  if ( Z>=1 ) // exclude e-, neutrons, gamma, ...
182  matrix -> Fill(trackID, particleDef, i, j, k, energyDeposit);
183  /*
184  * Create a hit with the information of position is in the detector
185  */
187  detectorHit -> SetEdepAndPosition(i, j, k, energyDeposit);
188  HitsCollection -> insert(detectorHit);
189  }
190  }
191 
192 #ifdef G4ANALYSIS_USE_ROOT
193  if(energyDeposit != 0)
194  {
195  if(trackID != 1)
196  {
197  if (particleName == "proton")
198  analysis -> SecondaryProtonEnergyDeposit(i, energyDeposit/MeV);
199 
200  else if (particleName == "neutron")
201  analysis -> SecondaryNeutronEnergyDeposit(i, energyDeposit/MeV);
202 
203  else if (particleName == "alpha")
204  analysis -> SecondaryAlphaEnergyDeposit(i, energyDeposit/MeV);
205 
206  else if (particleName == "gamma")
207  analysis -> SecondaryGammaEnergyDeposit(i, energyDeposit/MeV);
208 
209  else if (particleName == "e-")
210  analysis -> SecondaryElectronEnergyDeposit(i, energyDeposit/MeV);
211 
212  else if (particleName == "triton")
213  analysis -> SecondaryTritonEnergyDeposit(i, energyDeposit/MeV);
214 
215  else if (particleName == "deuteron")
216  analysis -> SecondaryDeuteronEnergyDeposit(i, energyDeposit/MeV);
217 
218  else if (particleName == "pi+" || particleName == "pi-" || particleName == "pi0")
219  analysis -> SecondaryPionEnergyDeposit(i, energyDeposit/MeV);
220  }
221  }
222 #endif
223 
224  return true;
225 }
226 
229 {
230 
231  static G4int HCID = -1;
232  if(HCID < 0)
233  {
234  HCID = GetCollectionID(0);
235  }
236 
237  HCE -> AddHitsCollection(HCID,HitsCollection);
238 }
239 
static HadrontherapyAnalysisManager * GetInstance()
Get the pointer to the analysis manager.
static const double MeV
Definition: G4SIunits.hh:193
static HadrontherapyLet * GetInstance()
G4String name
Definition: TRTMaterials.hh:40
const G4VTouchable * GetTouchable() const
int G4int
Definition: G4Types.hh:78
virtual G4int GetCollectionID(G4int i)
G4StepPoint * GetPreStepPoint() const
static HadrontherapyMatrix * GetInstance()
bool G4bool
Definition: G4Types.hh:79
G4THitsCollection< HadrontherapyDetectorHit > HadrontherapyDetectorHitsCollection
Definition: G4Step.hh:76
static const G4double A[nN]
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
void Initialize(G4HCofThisEvent *)
A class for connecting the simulation to an analysis package.
HadrontherapyDetectorHitsCollection * HitsCollection
void EndOfEvent(G4HCofThisEvent *HCE)
G4CollectionNameVector collectionName
static const double keV
Definition: G4SIunits.hh:195
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist)
double G4double
Definition: G4Types.hh:76