Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/hadrontherapy/src/HadrontherapyDetectorSD.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // Hadrontherapy advanced example for Geant4
 27 // See more at: https://twiki.cern.ch/twiki/bin/view/Geant4/AdvancedExamplesHadrontherapy
 28 
 29 #include "HadrontherapyDetectorSD.hh"
 30 #include "HadrontherapyDetectorHit.hh"
 31 #include "G4Step.hh"
 32 #include "G4VTouchable.hh"
 33 #include "G4TouchableHistory.hh"
 34 #include "G4SDManager.hh"
 35 #include "HadrontherapyMatrix.hh"
 36 #include "HadrontherapyLet.hh"
 37 #include "G4Track.hh"
 38 #include "G4SystemOfUnits.hh"
 39 #include "HadrontherapyMatrix.hh"
 40 
 41 
 42 #include "G4SteppingManager.hh"
 43 #include "G4TrackVector.hh"
 44 #include "HadrontherapySteppingAction.hh"
 45 #include "G4ios.hh"
 46 #include "G4SteppingManager.hh"
 47 #include "G4Track.hh"
 48 #include "G4Step.hh"
 49 #include "G4StepPoint.hh"
 50 #include "G4TrackStatus.hh"
 51 #include "G4TrackVector.hh"
 52 #include "G4ParticleDefinition.hh"
 53 #include "G4ParticleTypes.hh"
 54 #include "G4UserEventAction.hh"
 55 #include "G4TransportationManager.hh"
 56 #include "G4VSensitiveDetector.hh"
 57 #include "HadrontherapyRunAction.hh"
 58 #include "G4SystemOfUnits.hh"
 59 #include "HadrontherapyRBE.hh"
 60 #include <G4AccumulableManager.hh>
 61 
 62 
 63 /////////////////////////////////////////////////////////////////////////////
 64 HadrontherapyDetectorSD::HadrontherapyDetectorSD(G4String name):
 65 G4VSensitiveDetector(name)
 66 {
 67     G4String HCname;
 68     collectionName.insert(HCname="HadrontherapyDetectorHitsCollection");
 69     HitsCollection = NULL;
 70     sensitiveDetectorName = name;
 71     
 72 }
 73 
 74 /////////////////////////////////////////////////////////////////////////////
 75 HadrontherapyDetectorSD::~HadrontherapyDetectorSD()
 76 {}
 77 
 78 /////////////////////////////////////////////////////////////////////////////
 79 void HadrontherapyDetectorSD::Initialize(G4HCofThisEvent*)
 80 {
 81     
 82     HitsCollection = new HadrontherapyDetectorHitsCollection(sensitiveDetectorName,
 83                                                              collectionName[0]);
 84 }
 85 
 86 /////////////////////////////////////////////////////////////////////////////
 87 G4bool HadrontherapyDetectorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* )
 88 {
 89     
 90     
 91     if (aStep -> GetPreStepPoint() -> GetPhysicalVolume() -> GetName() != "RODetectorZDivisionPhys") return false;
 92     
 93     
 94     // Get kinetic energy
 95     G4Track * theTrack = aStep  ->  GetTrack();
 96     G4double kineticEnergy = theTrack->GetKineticEnergy();
 97     G4ParticleDefinition *particleDef = theTrack -> GetDefinition();
 98     //Get particle name
 99     G4String particleName =  particleDef -> GetParticleName();
100     
101     // Get particle PDG code
102     G4int pdg = particleDef ->GetPDGEncoding();
103     
104     // Get unique track_id (in an event)
105     G4int trackID = theTrack -> GetTrackID();
106     
107     G4double energyDeposit = aStep -> GetTotalEnergyDeposit();
108     
109     G4double DX = aStep -> GetStepLength();
110     G4int Z = particleDef-> GetAtomicNumber();
111     G4int A = particleDef-> GetAtomicMass();
112     G4StepPoint* PreStep = aStep->GetPreStepPoint();
113     
114     // Read voxel indexes: i is the x index, k is the z index
115     const G4VTouchable* touchable = aStep->GetPreStepPoint()->GetTouchable();
116     G4int k  = touchable->GetReplicaNumber(0);
117     G4int i  = touchable->GetReplicaNumber(2);
118     G4int j  = touchable->GetReplicaNumber(1);
119     
120     G4TouchableHandle touchPreStep = PreStep->GetTouchableHandle();
121     G4VPhysicalVolume* volumePre = touchPreStep->GetVolume();
122     G4String namePre = volumePre->GetName();
123     
124 
125     
126     HadrontherapyMatrix* matrix = HadrontherapyMatrix::GetInstance();
127     HadrontherapyLet* let = HadrontherapyLet::GetInstance();
128     
129     G4int* hitTrack = matrix -> GetHitTrack(i,j,k);
130     
131     
132     //  ******************** let ***************************
133     if (let)
134     {
135         if ( !(Z==0 && A==1) ) // All but not neutrons
136         {
137             if( energyDeposit>0. && DX >0. )// calculate only energy deposit
138             {
139                 if (pdg !=22 && pdg !=11) // not gamma and electrons
140                 {
141                     
142                     // Get the pre-step kinetic energy
143                     G4double eKinPre = aStep -> GetPreStepPoint() -> GetKineticEnergy();
144                     // Get the post-step kinetic energy
145                     G4double eKinPost = aStep -> GetPostStepPoint() -> GetKineticEnergy();
146                     // Get the step average kinetic energy
147                     G4double eKinMean = (eKinPre + eKinPost) * 0.5;
148                     
149                     // get the material
150                     const G4Material * materialStep = aStep -> GetPreStepPoint() -> GetMaterial();
151                     
152                     // get the secondary paticles
153                     G4Step fstep = *theTrack -> GetStep();
154                     // store all the secondary partilce in current step
155                     const std::vector<const G4Track*> * secondary = fstep.GetSecondaryInCurrentStep();
156                     
157                     size_t SecondarySize = (*secondary).size();
158                     G4double EnergySecondary = 0.;
159                     
160                     // get secondary electrons energy deposited
161                     if (SecondarySize) // calculate only secondary particles
162                     {
163                         for (size_t numsec = 0; numsec< SecondarySize ; numsec ++)
164                         {
165                             //Get the PDG code of every secondaty particles in current step
166                             G4int PDGSecondary=(*secondary)[numsec]->GetDefinition()->GetPDGEncoding();
167                             
168                             if(PDGSecondary == 11) // calculate only secondary electrons
169                             {
170                                 // calculate the energy deposit of secondary electrons in current step
171                                 EnergySecondary += (*secondary)[numsec]->GetKineticEnergy();
172                             }
173                         }
174                         
175                     }
176                     
177                     // call the let filldatas function to calculate let
178                     let -> FillEnergySpectrum(trackID, particleDef, eKinMean, materialStep,
179                                               energyDeposit,EnergySecondary,DX, i, j, k);
180                 }
181             }
182         }
183     }
184     
185     
186     if (matrix)
187     {
188         
189         // Increment Fluences & accumulate energy spectra
190         // Hit voxels are marked with track_id throught hitTrack matrix
191         //G4int* hitTrack = matrix -> GetHitTrack(i,j,k); // hitTrack MUST BE cleared at every eventAction!
192         if ( *hitTrack != trackID )
193         {
194             *hitTrack = trackID;
195             /*
196              * Fill FLUENCE data for every single nuclide
197              * Exclude e-, neutrons, gamma, ...
198              */
199             if ( Z >= 1) matrix -> Fill(trackID, particleDef, i, j, k, 0, true);
200             
201         }
202         
203         if(energyDeposit != 0)
204         {
205             /*
206              *  This method will fill a dose matrix for every single nuclide.
207              *  A method of the HadrontherapyMatrix class (StoreDoseFluenceAscii())
208              *  is called automatically at the end of main (or via the macro command /analysis/writeDoseFile.
209              *  It permits to store all dose/fluence data into a single plane ASCII file.
210              */
211             // if (A==1 && Z==1) // primary and sec. protons
212             if ( Z>=1 )    //  exclude e-, neutrons, gamma, ...
213                 matrix -> Fill(trackID, particleDef, i, j, k, energyDeposit);
214             /*
215              * Create a hit with the information of position is in the detector
216              */
217             HadrontherapyDetectorHit* detectorHit = new HadrontherapyDetectorHit();
218             detectorHit -> SetEdepAndPosition(i, j, k, energyDeposit);
219             HitsCollection -> insert(detectorHit);
220         }
221     }
222 
223     auto rbe = HadrontherapyRBE::GetInstance();
224     if (rbe->IsCalculationEnabled())
225     {
226         if (!fRBEAccumulable)
227         {
228             fRBEAccumulable = dynamic_cast<HadrontherapyRBEAccumulable*>(G4AccumulableManager::Instance()->GetAccumulable("RBE"));
229             if (!fRBEAccumulable)
230             {
231                 G4Exception("HadrontherapyDetectorSD::ProcessHits", "NoAccumulable", FatalException, "Accumulable RBE not found.");
232             }
233         }
234   if (A>0) //protect against gammas, e- , etc
235     {
236       fRBEAccumulable->Accumulate(kineticEnergy / A, energyDeposit, DX, Z, i, j, k);
237     }
238     }
239 
240 
241     return true;
242 }
243 
244 /////////////////////////////////////////////////////////////////////////////
245 void HadrontherapyDetectorSD::EndOfEvent(G4HCofThisEvent* HCE)
246 {
247     
248     static G4int HCID = -1;
249     if(HCID < 0)
250     {
251         HCID = GetCollectionID(0);
252     }
253     
254     HCE -> AddHitsCollection(HCID,HitsCollection);
255 }
256 
257