File indexing completed on 2025-04-10 08:06:46
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "PrimaryKiller.hh"
0030
0031 #include "G4Event.hh"
0032 #include "G4RunManager.hh"
0033 #include "G4UIcmdWith3VectorAndUnit.hh"
0034 #include "G4UIcmdWithADoubleAndUnit.hh"
0035 #include "G4UnitsTable.hh"
0036
0037 namespace scavenger
0038 {
0039
0040
0041
0042 PrimaryKiller::PrimaryKiller(const G4String& name, const G4int& depth)
0043 : G4VPrimitiveScorer(name, depth),
0044 G4UImessenger(),
0045 fpELossUI(new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin", this)),
0046 fpAbortEventIfELossUpperThan(new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMax", this)),
0047 fpMinKineticE(new G4UIcmdWithADoubleAndUnit("/primaryKiller/minKineticE", this)),
0048 fpSizeUI(new G4UIcmdWith3VectorAndUnit("/primaryKiller/setSize", this))
0049 {
0050 fpSizeUI->SetDefaultUnit("um");
0051 }
0052
0053
0054
0055 void PrimaryKiller::SetNewValue(G4UIcommand* command, G4String newValue)
0056 {
0057 if (command == fpELossUI.get()) {
0058 fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue);
0059 }
0060 else if (command == fpAbortEventIfELossUpperThan.get()) {
0061 fELossRange_Max = fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue);
0062 }
0063 else if (command == fpSizeUI.get()) {
0064 fPhantomSize = fpSizeUI->GetNew3VectorValue(newValue);
0065 }
0066 }
0067
0068
0069
0070 G4bool PrimaryKiller::ProcessHits(G4Step* aStep, G4TouchableHistory*)
0071 {
0072 const G4Track* track = aStep->GetTrack();
0073 G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
0074 if (std::abs(pos.x()) > fPhantomSize.getX() / 2 || std::abs(pos.y()) > fPhantomSize.getY() / 2
0075 || std::abs(pos.z()) > fPhantomSize.getZ() / 2)
0076 {
0077 ((G4Track*)track)->SetTrackStatus(fStopAndKill);
0078 return false;
0079 }
0080 if (track->GetTrackID() != 1) {
0081 return FALSE;
0082 }
0083 G4double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy();
0084 G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy() - kineticE;
0085 if (eLoss == 0.) {
0086 return FALSE;
0087 }
0088 fELoss += eLoss;
0089 if (fELoss > fELossRange_Max) {
0090 G4RunManager::GetRunManager()->AbortEvent();
0091 }
0092 if (fELoss >= fELossRange_Min || kineticE <= fKineticE_Min) {
0093 ((G4Track*)track)->SetTrackStatus(fStopAndKill);
0094 }
0095 return TRUE;
0096 }
0097
0098
0099
0100 void PrimaryKiller::Initialize(G4HCofThisEvent* )
0101 {
0102 fELoss = 0.;
0103 }
0104
0105
0106
0107 }