Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/medical/dna/chem6/src/ScoreLET.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 // This example is provided by the Geant4-DNA collaboration
 27 // chem6 example is derived from chem4 and chem5 examples
 28 //
 29 // Any report or published results obtained using the Geant4-DNA software
 30 // shall cite the following Geant4-DNA collaboration publication:
 31 // J. Appl. Phys. 125 (2019) 104301
 32 // Med. Phys. 45 (2018) e722-e739
 33 // J. Comput. Phys. 274 (2014) 841-882
 34 // Med. Phys. 37 (2010) 4692-4708
 35 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
 36 // The Geant4-DNA web site is available at http://geant4-dna.org
 37 //
 38 // Authors: W. G. Shin and S. Incerti (CENBG, France)
 39 //
 40 // $Id$
 41 //
 42 /// \file ScoreLET.cc
 43 /// \brief Implementation of the ScoreLET class
 44 
 45 #include "ScoreLET.hh"
 46 
 47 #include "G4RunManager.hh"
 48 #include "G4SystemOfUnits.hh"
 49 
 50 ScoreLET::ScoreLET(G4String name) : G4VPrimitiveScorer(name), G4UImessenger(), fEvtMap(0)
 51 {
 52   fLET = 0;
 53   fEdep = 0;
 54   fStepL = 0;
 55   fNEvent = 0;
 56   fTrackID = 1;
 57 
 58   fpLETDir = new G4UIdirectory("/scorer/LET/");
 59   fpLETDir->SetGuidance("LET scorer commands");
 60 
 61   fpCutoff = new G4UIcmdWithADoubleAndUnit("/scorer/LET/cutoff", this);
 62 
 63   fCutoff = DBL_MAX;
 64 }
 65 
 66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 67 
 68 ScoreLET::~ScoreLET()
 69 {
 70   delete fpLETDir;
 71   delete fpCutoff;
 72 }
 73 
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 75 
 76 void ScoreLET::SetNewValue(G4UIcommand* command, G4String newValue)
 77 {
 78   if (command == fpCutoff) fCutoff = atof(newValue);
 79 }
 80 
 81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 82 
 83 G4bool ScoreLET::ProcessHits(G4Step* aStep, G4TouchableHistory* /*TH*/)
 84 {
 85   // In order to follow the primary track
 86   // regardless charge increasing or decreasing
 87   if (aStep->GetTrack()->GetTrackID() != 1
 88       && aStep->GetTrack()->GetParticleDefinition()->GetPDGEncoding() != 11)
 89   {
 90     G4int subType = aStep->GetTrack()->GetCreatorProcess()->GetProcessSubType();
 91     if (subType == 56 || subType == 57) {
 92       fTrackID = aStep->GetTrack()->GetTrackID();
 93     }
 94   }
 95 
 96   // Ignore the step if it is not primary.
 97   if (aStep->GetTrack()->GetTrackID() != fTrackID)
 98     return false;
 99   else {
100     fStepL += aStep->GetStepLength() / um;
101     fEdep += aStep->GetTotalEnergyDeposit() / keV;
102 
103     G4int subType = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessSubType();
104 
105     // Don't add the kinetic energy of primary particle
106     if (subType == 56 || subType == 57) return false;
107 
108     const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep();
109 
110     size_t nbtrk = (*secondary).size();
111 
112     if (nbtrk) {
113       for (size_t lp = 0; lp < nbtrk; lp++) {
114         // Store the kinetic energy of secondaries
115         // which less than cutoff energy.
116         if ((*secondary)[lp]->GetKineticEnergy() / eV < fCutoff) {
117           fEdep += (*secondary)[lp]->GetKineticEnergy() / keV;
118         }
119       }
120     }
121   }
122   return true;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
126 
127 void ScoreLET::Initialize(G4HCofThisEvent* HCE)
128 {
129   fEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
130 
131   fEvtMap->clear();
132   fNEvent = 0;
133   static G4int HCID = -1;
134   if (HCID < 0) {
135     HCID = GetCollectionID(0);
136   }
137   HCE->AddHitsCollection(HCID, (G4VHitsCollection*)fEvtMap);
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
141 
142 void ScoreLET::EndOfEvent(G4HCofThisEvent*)
143 {
144   fLET = fEdep / fStepL;
145   if (!G4RunManager::GetRunManager()->GetCurrentEvent()->IsAborted()) {
146     fEvtMap->add(fNEvent, fLET);
147     fNEvent++;
148   }
149   fTrackID = 1;
150   fLET = 0;
151   fEdep = 0;
152   fStepL = 0;
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
156 
157 G4int ScoreLET::GetIndex(G4Step* /*aStep*/)
158 {
159   return 0;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
163