Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:50

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 // This example is provided by the Geant4-DNA collaboration
0027 // chem6 example is derived from chem4 and chem5 examples
0028 //
0029 // Any report or published results obtained using the Geant4-DNA software
0030 // shall cite the following Geant4-DNA collaboration publication:
0031 // J. Appl. Phys. 125 (2019) 104301
0032 // Med. Phys. 45 (2018) e722-e739
0033 // J. Comput. Phys. 274 (2014) 841-882
0034 // Med. Phys. 37 (2010) 4692-4708
0035 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
0036 // The Geant4-DNA web site is available at http://geant4-dna.org
0037 //
0038 // Authors: W. G. Shin and S. Incerti (CENBG, France)
0039 //
0040 // $Id$
0041 //
0042 /// \file ScoreLET.cc
0043 /// \brief Implementation of the ScoreLET class
0044 
0045 #include "ScoreLET.hh"
0046 
0047 #include "G4RunManager.hh"
0048 #include "G4SystemOfUnits.hh"
0049 
0050 ScoreLET::ScoreLET(G4String name) : G4VPrimitiveScorer(name), G4UImessenger(), fEvtMap(0)
0051 {
0052   fLET = 0;
0053   fEdep = 0;
0054   fStepL = 0;
0055   fNEvent = 0;
0056   fTrackID = 1;
0057 
0058   fpLETDir = new G4UIdirectory("/scorer/LET/");
0059   fpLETDir->SetGuidance("LET scorer commands");
0060 
0061   fpCutoff = new G4UIcmdWithADoubleAndUnit("/scorer/LET/cutoff", this);
0062 
0063   fCutoff = DBL_MAX;
0064 }
0065 
0066 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0067 
0068 ScoreLET::~ScoreLET()
0069 {
0070   delete fpLETDir;
0071   delete fpCutoff;
0072 }
0073 
0074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0075 
0076 void ScoreLET::SetNewValue(G4UIcommand* command, G4String newValue)
0077 {
0078   if (command == fpCutoff) fCutoff = atof(newValue);
0079 }
0080 
0081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0082 
0083 G4bool ScoreLET::ProcessHits(G4Step* aStep, G4TouchableHistory* /*TH*/)
0084 {
0085   // In order to follow the primary track
0086   // regardless charge increasing or decreasing
0087   if (aStep->GetTrack()->GetTrackID() != 1
0088       && aStep->GetTrack()->GetParticleDefinition()->GetPDGEncoding() != 11)
0089   {
0090     G4int subType = aStep->GetTrack()->GetCreatorProcess()->GetProcessSubType();
0091     if (subType == 56 || subType == 57) {
0092       fTrackID = aStep->GetTrack()->GetTrackID();
0093     }
0094   }
0095 
0096   // Ignore the step if it is not primary.
0097   if (aStep->GetTrack()->GetTrackID() != fTrackID)
0098     return false;
0099   else {
0100     fStepL += aStep->GetStepLength() / um;
0101     fEdep += aStep->GetTotalEnergyDeposit() / keV;
0102 
0103     G4int subType = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessSubType();
0104 
0105     // Don't add the kinetic energy of primary particle
0106     if (subType == 56 || subType == 57) return false;
0107 
0108     const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep();
0109 
0110     size_t nbtrk = (*secondary).size();
0111 
0112     if (nbtrk) {
0113       for (size_t lp = 0; lp < nbtrk; lp++) {
0114         // Store the kinetic energy of secondaries
0115         // which less than cutoff energy.
0116         if ((*secondary)[lp]->GetKineticEnergy() / eV < fCutoff) {
0117           fEdep += (*secondary)[lp]->GetKineticEnergy() / keV;
0118         }
0119       }
0120     }
0121   }
0122   return true;
0123 }
0124 
0125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0126 
0127 void ScoreLET::Initialize(G4HCofThisEvent* HCE)
0128 {
0129   fEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
0130 
0131   fEvtMap->clear();
0132   fNEvent = 0;
0133   static G4int HCID = -1;
0134   if (HCID < 0) {
0135     HCID = GetCollectionID(0);
0136   }
0137   HCE->AddHitsCollection(HCID, (G4VHitsCollection*)fEvtMap);
0138 }
0139 
0140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0141 
0142 void ScoreLET::EndOfEvent(G4HCofThisEvent*)
0143 {
0144   fLET = fEdep / fStepL;
0145   if (!G4RunManager::GetRunManager()->GetCurrentEvent()->IsAborted()) {
0146     fEvtMap->add(fNEvent, fLET);
0147     fNEvent++;
0148   }
0149   fTrackID = 1;
0150   fLET = 0;
0151   fEdep = 0;
0152   fStepL = 0;
0153 }
0154 
0155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0156 
0157 G4int ScoreLET::GetIndex(G4Step* /*aStep*/)
0158 {
0159   return 0;
0160 }
0161 
0162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....