Back to home page

EIC code displayed by LXR

 
 

    


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

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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0051 
0052 ScoreLET::ScoreLET(G4String name) : G4VPrimitiveScorer(name), G4UImessenger(), fEvtMap(0)
0053 {
0054   fpLETDir = new G4UIdirectory("/scorer/LET/");
0055   fpLETDir->SetGuidance("LET scorer commands");
0056 
0057   fpCutoff = new G4UIcmdWithADoubleAndUnit("/scorer/LET/cutoff", this);
0058 }
0059 
0060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0061 
0062 ScoreLET::~ScoreLET()
0063 {
0064   delete fpLETDir;
0065   delete fpCutoff;
0066 }
0067 
0068 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0069 
0070 void ScoreLET::SetNewValue(G4UIcommand* command, G4String newValue)
0071 {
0072   if (command == fpCutoff) fCutoff = atof(newValue);
0073 }
0074 
0075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0076 
0077 G4bool ScoreLET::ProcessHits(G4Step* aStep, G4TouchableHistory* /*TH*/)
0078 {
0079   // In order to follow the primary track
0080   // regardless charge increasing or decreasing
0081   if (aStep->GetTrack()->GetTrackID() != 1
0082       && aStep->GetTrack()->GetParticleDefinition()->GetPDGEncoding() != 11)
0083   {
0084     G4int subType = aStep->GetTrack()->GetCreatorProcess()->GetProcessSubType();
0085     if (subType == 56 || subType == 57) {
0086       fTrackID = aStep->GetTrack()->GetTrackID();
0087     }
0088   }
0089 
0090   // Ignore the step if it is not primary.
0091   if (aStep->GetTrack()->GetTrackID() != fTrackID)
0092     return false;
0093   else {
0094     fStepL += aStep->GetStepLength() / um;
0095     fEdep += aStep->GetTotalEnergyDeposit() / keV;
0096 
0097     G4int subType = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessSubType();
0098 
0099     // Don't add the kinetic energy of primary particle
0100     if (subType == 56 || subType == 57) return false;
0101 
0102     const std::vector<const G4Track*>* secondary = aStep->GetSecondaryInCurrentStep();
0103 
0104     size_t nbtrk = (*secondary).size();
0105 
0106     if (nbtrk) {
0107       for (size_t lp = 0; lp < nbtrk; lp++) {
0108         // Store the kinetic energy of secondaries
0109         // which less than cutoff energy.
0110         if ((*secondary)[lp]->GetKineticEnergy() / eV < fCutoff) {
0111           fEdep += (*secondary)[lp]->GetKineticEnergy() / keV;
0112         }
0113       }
0114     }
0115   }
0116   return true;
0117 }
0118 
0119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0120 
0121 void ScoreLET::Initialize(G4HCofThisEvent* HCE)
0122 {
0123   fEvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(), GetName());
0124 
0125   fEvtMap->clear();
0126   fNEvent = 0;
0127   static G4int HCID = -1;
0128   if (HCID < 0) {
0129     HCID = GetCollectionID(0);
0130   }
0131   HCE->AddHitsCollection(HCID, (G4VHitsCollection*)fEvtMap);
0132 }
0133 
0134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0135 
0136 void ScoreLET::EndOfEvent(G4HCofThisEvent*)
0137 {
0138   if (fStepL > 0) {
0139     fLET = fEdep / fStepL;
0140     if (!G4RunManager::GetRunManager()->GetCurrentEvent()->IsAborted()) {
0141       fEvtMap->add(fNEvent, fLET);
0142       fNEvent++;
0143     }
0144   }
0145   fTrackID = 1;
0146   fLET = 0;
0147   fEdep = 0;
0148   fStepL = 0;
0149 }
0150 
0151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0152 
0153 G4int ScoreLET::GetIndex(G4Step* /*aStep*/)
0154 {
0155   return 0;
0156 }
0157 
0158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....