Geant4 Cross Reference |
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 // Any report or published results obtained using the Geant4-DNA software 28 // shall cite the following Geant4-DNA collaboration publication: 29 // Med. Phys. 37 (2010) 4692-4708 30 // J. Comput. Phys. 274 (2014) 841-882 31 // Phys. Med. Biol. 63(10) (2018) 105014-12pp 32 // The Geant4-DNA web site is available at http://geant4-dna.org 33 // 34 // 35 #include "PrimaryKiller.hh" 36 37 #include <G4Event.hh> 38 #include <G4EventManager.hh> 39 #include <G4RunManager.hh> 40 #include <G4SystemOfUnits.hh> 41 #include <G4UIcmdWithADoubleAndUnit.hh> 42 #include <G4UnitsTable.hh> 43 #include <globals.hh> 44 45 /** \file PrimaryKiller.cc 46 \class PrimaryKiller 47 48 Kill the primary particle: 49 - either after a given energy loss 50 - or after the primary particle has reached a given energy 51 */ 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 PrimaryKiller::PrimaryKiller(G4String name, G4int depth) 56 : G4VPrimitiveScorer(name, depth), G4UImessenger() 57 { 58 fELoss = 0.; // cumulated energy for current event 59 60 fELossRange_Min = DBL_MAX; // fELoss from which the primary is killed 61 fELossRange_Max = DBL_MAX; // fELoss from which the event is aborted 62 fKineticE_Min = 0; // kinetic energy below which the primary is killed 63 64 fpELossUI = new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin", this); 65 fpAbortEventIfELossUpperThan = new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMax", this); 66 fpMinKineticE = new G4UIcmdWithADoubleAndUnit("/primaryKiller/minKineticE", this); 67 } 68 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 70 71 PrimaryKiller::~PrimaryKiller() 72 { 73 ; 74 } 75 76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 77 78 void PrimaryKiller::SetNewValue(G4UIcommand* command, G4String newValue) 79 { 80 if (command == fpELossUI) { 81 fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue); 82 } 83 else if (command == fpAbortEventIfELossUpperThan) { 84 fELossRange_Max = fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue); 85 } 86 } 87 88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 89 90 G4bool PrimaryKiller::ProcessHits(G4Step* aStep, G4TouchableHistory*) 91 { 92 const G4Track* track = aStep->GetTrack(); 93 94 if (track->GetTrackID() != 1) return FALSE; 95 96 //------------------- 97 98 double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy(); 99 100 G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy() - kineticE; 101 102 if (eLoss == 0.) return FALSE; 103 104 //------------------- 105 106 fELoss += eLoss; 107 108 if (fELoss > fELossRange_Max) { 109 G4RunManager::GetRunManager()->AbortEvent(); 110 int eventID = G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID(); 111 G4cout << "Aborts event " << eventID 112 << ", energy loss " 113 "is too large. \n" 114 << " * Energy loss by primary is: " << G4BestUnit(fELoss, "Energy") 115 << ". Event is aborted if the Eloss is > " << G4BestUnit(fELossRange_Max, "Energy") 116 << G4endl; 117 } 118 119 if (fELoss >= fELossRange_Min) { //|| kineticE <= fKineticE_Min){ 120 ((G4Track*)track)->SetTrackStatus(fStopAndKill); 121 G4cout << "Kill track at : " << G4BestUnit(kineticE, "Energy") 122 << ", Energy loss by primary is: " << G4BestUnit(fELoss, "Energy") 123 << ", primary is terminated as Eloss is >: " << G4BestUnit(fELossRange_Min, "Energy") 124 << G4endl; //", EThreshold is: " 125 // << G4BestUnit(fEThreshold, "Energy") 126 // << G4endl; 127 } 128 129 return true; 130 } 131 132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 133 134 void PrimaryKiller::Initialize(G4HCofThisEvent* /*HCE*/) 135 { 136 Clear(); 137 } 138 139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 140 141 void PrimaryKiller::EndOfEvent(G4HCofThisEvent*) {} 142 143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 144 145 void PrimaryKiller::Clear() 146 { 147 fELoss = 0.; 148 } 149 150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 151 152 void PrimaryKiller::DrawAll() 153 { 154 ; 155 } 156 157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 158 159 void PrimaryKiller::PrintAll() {} 160 161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 162