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