Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
scavenger::PrimaryKiller Class Reference

Kill the primary particle: More...

#include <Doxymodules_medical.h>

Inheritance diagram for scavenger::PrimaryKiller:
G4VPrimitiveScorer G4UImessenger

Public Member Functions

 PrimaryKiller (const G4String &name, const G4int &depth=0)
 
 ~PrimaryKiller () override=default
 
void SetMinLossEnergyLimit (const G4double &energy)
 Set the energy loss from which the primary is killed.
 
void SetMaxLossEnergyLimit (const G4double &energy)
 Set the energy loss from which the event is aborted.
 
void SetNewValue (G4UIcommand *command, G4String newValue) override
 Method related to G4UImessenger used to control energy cuts through macro file.
 
void Initialize (G4HCofThisEvent *) override
 
void EndOfEvent (G4HCofThisEvent *) override
 
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override
 

Private Attributes

G4double fELoss = 0
 
G4double fELossRange_Min = DBL_MAX
 
G4double fELossRange_Max = DBL_MAX
 
G4double fKineticE_Min = 0
 
G4ThreeVector fPhantomSize = G4ThreeVector(1 * km, 1 * km, 1 * km)
 
std::unique_ptr< G4UIcmdWithADoubleAndUnitfpELossUI
 
std::unique_ptr< G4UIcmdWithADoubleAndUnitfpAbortEventIfELossUpperThan
 
std::unique_ptr< G4UIcmdWithADoubleAndUnitfpMinKineticE
 
std::unique_ptr< G4UIcmdWith3VectorAndUnitfpSizeUI
 

Detailed Description

Kill the primary particle:

Definition at line 350 of file Doxymodules_medical.h.

Constructor & Destructor Documentation

◆ PrimaryKiller()

scavenger::PrimaryKiller::PrimaryKiller ( const G4String name,
const G4int &  depth = 0 
)
explicit

Definition at line 41 of file PrimaryKiller.cc.

42 : G4VPrimitiveScorer(name, depth),
44 fpELossUI(new G4UIcmdWithADoubleAndUnit("/primaryKiller/eLossMin", this)),
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}
std::unique_ptr< G4UIcmdWithADoubleAndUnit > fpMinKineticE
std::unique_ptr< G4UIcmdWithADoubleAndUnit > fpELossUI
std::unique_ptr< G4UIcmdWithADoubleAndUnit > fpAbortEventIfELossUpperThan
std::unique_ptr< G4UIcmdWith3VectorAndUnit > fpSizeUI

◆ ~PrimaryKiller()

scavenger::PrimaryKiller::~PrimaryKiller ( )
overridedefault

Member Function Documentation

◆ SetMinLossEnergyLimit()

void scavenger::PrimaryKiller::SetMinLossEnergyLimit ( const G4double &  energy)
inline

Set the energy loss from which the primary is killed.

Definition at line 59 of file PrimaryKiller.hh.

59 {
60 fELossRange_Min = energy;
61 }

◆ SetMaxLossEnergyLimit()

void scavenger::PrimaryKiller::SetMaxLossEnergyLimit ( const G4double &  energy)
inline

Set the energy loss from which the event is aborted.

Definition at line 65 of file PrimaryKiller.hh.

65 {
66 fELossRange_Max = energy;
67 }

◆ SetNewValue()

void scavenger::PrimaryKiller::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
override

Method related to G4UImessenger used to control energy cuts through macro file.

Definition at line 54 of file PrimaryKiller.cc.

55 {
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}

◆ Initialize()

void scavenger::PrimaryKiller::Initialize ( G4HCofThisEvent )
override

Definition at line 98 of file PrimaryKiller.cc.

98 {
99 fELoss = 0.;
100}

◆ EndOfEvent()

void scavenger::PrimaryKiller::EndOfEvent ( G4HCofThisEvent )
inlineoverride

Definition at line 75 of file PrimaryKiller.hh.

75{};

◆ ProcessHits()

G4bool scavenger::PrimaryKiller::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
override

Definition at line 68 of file PrimaryKiller.cc.

68 {
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}

Member Data Documentation

◆ fELoss

G4double scavenger::PrimaryKiller::fELoss = 0
private

Definition at line 80 of file PrimaryKiller.hh.

◆ fELossRange_Min

G4double scavenger::PrimaryKiller::fELossRange_Min = DBL_MAX
private

Definition at line 81 of file PrimaryKiller.hh.

◆ fELossRange_Max

G4double scavenger::PrimaryKiller::fELossRange_Max = DBL_MAX
private

Definition at line 82 of file PrimaryKiller.hh.

◆ fKineticE_Min

G4double scavenger::PrimaryKiller::fKineticE_Min = 0
private

Definition at line 83 of file PrimaryKiller.hh.

◆ fPhantomSize

G4ThreeVector scavenger::PrimaryKiller::fPhantomSize = G4ThreeVector(1 * km, 1 * km, 1 * km)
private

Definition at line 84 of file PrimaryKiller.hh.

◆ fpELossUI

std::unique_ptr<G4UIcmdWithADoubleAndUnit> scavenger::PrimaryKiller::fpELossUI
private

Definition at line 85 of file PrimaryKiller.hh.

◆ fpAbortEventIfELossUpperThan

std::unique_ptr<G4UIcmdWithADoubleAndUnit> scavenger::PrimaryKiller::fpAbortEventIfELossUpperThan
private

Definition at line 86 of file PrimaryKiller.hh.

◆ fpMinKineticE

std::unique_ptr<G4UIcmdWithADoubleAndUnit> scavenger::PrimaryKiller::fpMinKineticE
private

Definition at line 87 of file PrimaryKiller.hh.

◆ fpSizeUI

std::unique_ptr<G4UIcmdWith3VectorAndUnit> scavenger::PrimaryKiller::fpSizeUI
private

Definition at line 88 of file PrimaryKiller.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration