Loading...
Searching...
No Matches
PrimaryKiller.cc
Go to the documentation of this file.
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/// \file scavenger/src/PrimaryKiller.cc
27/// \brief Implementation of the scavenger::PrimaryKiller class
28
29#include "PrimaryKiller.hh"
30#include "G4UnitsTable.hh"
31#include "G4Event.hh"
32#include "G4UIcmdWithADoubleAndUnit.hh"
33#include "G4UIcmdWith3VectorAndUnit.hh"
34#include "G4RunManager.hh"
35
36namespace scavenger
37{
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
40
41PrimaryKiller::PrimaryKiller(const G4String &name, const G4int &depth)
42 : G4VPrimitiveScorer(name, depth),
44 fpELossUI(new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin", this)),
45 fpAbortEventIfELossUpperThan(
46 new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMax", this)),
47 fpMinKineticE(new G4UIcmdWithADoubleAndUnit("/primaryKiller/minKineticE", this)),
48 fpSizeUI(new G4UIcmdWith3VectorAndUnit("/primaryKiller/setSize", this)) {
49 fpSizeUI->SetDefaultUnit("um");
50}
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
53
55 G4String newValue) {
56 if (command == fpELossUI.get()) {
57 fELossRange_Min = fpELossUI->GetNewDoubleValue(newValue);
58 } else if (command == fpAbortEventIfELossUpperThan.get()) {
60 fpAbortEventIfELossUpperThan->GetNewDoubleValue(newValue);
61 } else if (command == fpSizeUI.get()) {
62 fPhantomSize = fpSizeUI->GetNew3VectorValue(newValue);
63 }
64}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
67
69 const G4Track *track = aStep->GetTrack();
70 G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
71 if (std::abs(pos.x()) > fPhantomSize.getX() / 2 ||
72 std::abs(pos.y()) > fPhantomSize.getY() / 2 ||
73 std::abs(pos.z()) > fPhantomSize.getZ() / 2) {
74 ((G4Track *) track)->SetTrackStatus(fStopAndKill);
75 return false;
76 }
77 if (track->GetTrackID() != 1) {
78 return FALSE;
79 }
80 G4double kineticE = aStep->GetPostStepPoint()->GetKineticEnergy();
81 G4double eLoss = aStep->GetPreStepPoint()->GetKineticEnergy()
82 - kineticE;
83 if (eLoss == 0.) {
84 return FALSE;
85 }
86 fELoss += eLoss;
87 if (fELoss > fELossRange_Max) {
88 G4RunManager::GetRunManager()->AbortEvent();
89 }
90 if (fELoss >= fELossRange_Min || kineticE <= fKineticE_Min) {
91 ((G4Track *) track)->SetTrackStatus(fStopAndKill);
92 }
93 return TRUE;
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
97
99 fELoss = 0.;
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
103
104}
Definition of the scavenger::PrimaryKiller class.
void SetNewValue(G4UIcommand *command, G4String newValue) override
Method related to G4UImessenger used to control energy cuts through macro file.
PrimaryKiller(const G4String &name, const G4int &depth=0)
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
std::unique_ptr< G4UIcmdWithADoubleAndUnit > fpELossUI
void Initialize(G4HCofThisEvent *) override
std::unique_ptr< G4UIcmdWithADoubleAndUnit > fpAbortEventIfELossUpperThan
std::unique_ptr< G4UIcmdWith3VectorAndUnit > fpSizeUI

Applications | User Support | Publications | Collaboration