Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-02 07:37:06

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file ScoreLET.cc
0027 /// \brief Implementation of the ScoreLET class
0028 
0029 // This example is provided by the Geant4-DNA collaboration
0030 // chem6 example is derived from chem4 and chem5 examples
0031 //
0032 // Any report or published results obtained using the Geant4-DNA software
0033 // shall cite the following Geant4-DNA collaboration publication:
0034 // J. Appl. Phys. 125 (2019) 104301
0035 // Med. Phys. 45 (2018) e722-e739
0036 // J. Comput. Phys. 274 (2014) 841-882
0037 // Med. Phys. 37 (2010) 4692-4708
0038 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
0039 // The Geant4-DNA web site is available at http://geant4-dna.org
0040 //
0041 // Authors: W. G. Shin and S. Incerti (CENBG, France)
0042 
0043 #include "ScoreLET.hh"
0044 
0045 #include "G4RunManager.hh"
0046 #include "G4SystemOfUnits.hh"
0047 
0048 ScoreLET::ScoreLET(G4String name) : G4VPrimitiveScorer(name), G4UImessenger() {
0049   fpLETDir = std::make_unique<G4UIdirectory>("/scorer/LET/");
0050   fpLETDir->SetGuidance("LET scorer commands");
0051   fpCutoff = std::make_unique<G4UIcmdWithADoubleAndUnit>("/scorer/LET/cutoff", this);
0052 }
0053 
0054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0055 
0056 void ScoreLET::SetNewValue(G4UIcommand *command, const G4String newValue) {
0057   if (command == fpCutoff.get()) {
0058     fCutoff = G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(newValue);
0059   }
0060 }
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0063 
0064 G4bool ScoreLET::ProcessHits(G4Step *aStep, G4TouchableHistory * /*TH*/) {
0065   // In order to follow the primary track
0066   // regardless charge increasing or decreasing
0067   if (aStep->GetTrack()->GetTrackID() != 1
0068       && aStep->GetTrack()->GetParticleDefinition()->GetPDGEncoding() != 11) {
0069     G4int subType = aStep->GetTrack()->GetCreatorProcess()->GetProcessSubType();
0070     if (subType == 56 || subType == 57) {
0071       fTrackID = aStep->GetTrack()->GetTrackID();
0072     }
0073   }
0074 
0075   // Ignore the step if it is not primary.
0076   if (aStep->GetTrack()->GetTrackID() != fTrackID) { return false; }
0077   fStepL += aStep->GetStepLength() / um;
0078   fEdep += aStep->GetTotalEnergyDeposit() / keV;
0079 
0080   G4int subType = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessSubType();
0081 
0082   // Don't add the kinetic energy of primary particle
0083   if (subType == 56 || subType == 57) { return false; }
0084 
0085   const std::vector<const G4Track *> *secondary = aStep->GetSecondaryInCurrentStep();
0086 
0087   size_t nbtrk = secondary->size();
0088 
0089   if (nbtrk) {
0090     for (size_t lp = 0; lp < nbtrk; lp++) {
0091       // Store the kinetic energy of secondaries
0092       // which less than cutoff energy.
0093       if ((*secondary)[lp]->GetKineticEnergy() / eV < fCutoff / eV) {
0094         fEdep += (*secondary)[lp]->GetKineticEnergy() / keV;
0095       }
0096     }
0097   }
0098   return true;
0099 }
0100 
0101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0102 
0103 void ScoreLET::Initialize(G4HCofThisEvent *HCE) {
0104   fEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
0105 
0106   fEvtMap->clear();
0107   fNEvent = 0;
0108   static G4int HCID = -1;
0109   if (HCID < 0) {
0110     HCID = GetCollectionID(0);
0111   }
0112   HCE->AddHitsCollection(HCID, fEvtMap);
0113 }
0114 
0115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0116 
0117 void ScoreLET::EndOfEvent(G4HCofThisEvent *) {
0118   fLET = fEdep / fStepL;
0119   if (!G4RunManager::GetRunManager()->GetCurrentEvent()->IsAborted()) {
0120     fEvtMap->add(fNEvent, fLET);
0121     fNEvent++;
0122   }
0123   fTrackID = 1;
0124   fLET = 0;
0125   fEdep = 0;
0126   fStepL = 0;
0127 }
0128 
0129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....