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