LDMX Software
biasing::TargetProcessFilter Class Reference

Biases Geant4 to only process events where PN reaction occurred in the target. More...

#include <TargetProcessFilter.h>

Public Member Functions

 TargetProcessFilter (const std::string &name, framework::config::Parameters &parameters)
 Class constructor.
 
virtual ~TargetProcessFilter ()=default
 Destructor.
 
void stepping (const G4Step *step) override
 Implementmthe stepping action which performs the target volume biasing.
 
void EndOfEventAction (const G4Event *) override
 End of event action.
 
G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *aTrack, const G4ClassificationOfNewTrack &currentTrackClass) override
 Classify a new track which postpones track processing.
 
std::vector< simcore::TYPE > getTypes () override
 Retrieve the type of actions this class defines.
 
- Public Member Functions inherited from simcore::UserAction
 UserAction (const std::string &name, framework::config::Parameters &parameters)
 Constructor.
 
virtual ~UserAction ()=default
 Destructor.
 
virtual void BeginOfEventAction (const G4Event *)
 Method called at the beginning of every event.
 
virtual void BeginOfRunAction (const G4Run *)
 Method called at the beginning of a run.
 
virtual void EndOfRunAction (const G4Run *)
 Method called at the end of a run.
 
virtual void PreUserTrackingAction (const G4Track *)
 Method called before the UserTrackingAction.
 
virtual void PostUserTrackingAction (const G4Track *)
 Method called after the UserTrackingAction.
 
virtual void NewStage ()
 Method called at the beginning of a new stage.
 
virtual void PrepareNewEvent ()
 Method called at the beginning of a new event.
 

Private Attributes

G4Track * currentTrack_ {nullptr}
 Pointer to the current track being processed.
 
std::string process_ {""}
 The process to bias.
 

Additional Inherited Members

- Public Types inherited from simcore::UserAction
using Factory
 factory for user actions
 
- Protected Member Functions inherited from simcore::UserAction
UserEventInformationgetEventInfo () const
 Get a handle to the event information.
 
const std::map< int, ldmx::SimParticle > & getCurrentParticleMap () const
 Get the current particle map.
 
- Protected Attributes inherited from simcore::UserAction
std::string name_ {""}
 Name of the UserAction.
 
framework::config::Parameters parameters_
 The set of parameters used to configure this class.
 
mutable::framework::logging::logger theLog_
 the logging channel user actions can use ldmx_log with
 

Detailed Description

Biases Geant4 to only process events where PN reaction occurred in the target.

Definition at line 27 of file TargetProcessFilter.h.

Constructor & Destructor Documentation

◆ TargetProcessFilter()

biasing::TargetProcessFilter::TargetProcessFilter ( const std::string & name,
framework::config::Parameters & parameters )

Class constructor.

Definition at line 23 of file TargetProcessFilter.cxx.

25 : simcore::UserAction(name, parameters) {
26 process_ = parameters.getParameter<std::string>("process");
27}
std::string process_
The process to bias.
Interface that defines a user action.
Definition UserAction.h:43

References process_.

Member Function Documentation

◆ ClassifyNewTrack()

G4ClassificationOfNewTrack biasing::TargetProcessFilter::ClassifyNewTrack ( const G4Track * aTrack,
const G4ClassificationOfNewTrack & currentTrackClass )
overridevirtual

Classify a new track which postpones track processing.

Track processing resumes normally if a target PN interaction occurred.

Parameters
aTrackThe Geant4 track.
currentTrackClassThe current track classification.

Reimplemented from simcore::UserAction.

Definition at line 29 of file TargetProcessFilter.cxx.

30 {
31 if (track == currentTrack_) {
32 currentTrack_ = nullptr;
33 // std::cout << "[ TargetBremFilter ]: Pushing track to waiting stack." <<
34 // std::endl;
35 return fWaiting;
36 }
37
38 // Use current classification by default so values from other plugins are not
39 // overridden.
40 G4ClassificationOfNewTrack classification = currentTrackClass;
41
42 return classification;
43}
G4Track * currentTrack_
Pointer to the current track being processed.

References currentTrack_.

◆ EndOfEventAction()

void biasing::TargetProcessFilter::EndOfEventAction ( const G4Event * )
overridevirtual

End of event action.

Reimplemented from simcore::UserAction.

Definition at line 167 of file TargetProcessFilter.cxx.

167{}

◆ getTypes()

std::vector< simcore::TYPE > biasing::TargetProcessFilter::getTypes ( )
inlineoverridevirtual

Retrieve the type of actions this class defines.

Implements simcore::UserAction.

Definition at line 60 of file TargetProcessFilter.h.

60 {
61 return {simcore::TYPE::EVENT, simcore::TYPE::STACKING,
62 simcore::TYPE::STEPPING};
63 }

◆ stepping()

void biasing::TargetProcessFilter::stepping ( const G4Step * step)
overridevirtual

Implementmthe stepping action which performs the target volume biasing.

Parameters
stepThe Geant4 step.

Check if the electron will be exiting the target

The 'recoil_PV' volume name is automatically constructed by Geant4's GDML parser and was found by inspecting the geometry using a visualization. This Physical Volume (PV) is associated with the recoil parent volume and so it will break if the recoil parent volume changes its name.

We also check for 'World_PV' because in later geometries, there is an air gap between the target region and the recoil tracker.

Reimplemented from simcore::UserAction.

Definition at line 45 of file TargetProcessFilter.cxx.

45 {
46 // Get the track associated with this step.
47 auto track{step->GetTrack()};
48
49 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted())
50 return;
51
52 // Get the track info and check if this track is a brem candidate
53 auto trackInfo{simcore::UserTrackInformation::get(track)};
54 if ((trackInfo != nullptr) && !trackInfo->isBremCandidate()) return;
55
56 // Get the particles daughters.
57 auto secondaries{step->GetSecondary()};
58
59 // Get the region the particle is currently in. Continue processing
60 // the particle only if it's in the target region.
61 static auto target_region =
62 simcore::g4user::ptrretrieval::getRegion("target");
63 if (!target_region) {
64 ldmx_log(warn) << "Region 'target' not found in Geant4 region store";
65 }
66 auto phy_vol{track->GetVolume()};
67 auto log_vol{phy_vol ? phy_vol->GetLogicalVolume() : nullptr};
68 auto current_region{log_vol ? log_vol->GetRegion() : nullptr};
69 if (current_region != target_region) {
70 // If secondaries were produced outside of the volume of interest,
71 // and there aren't additional brems to process, abort the event.
72 // Otherwise, suspend the track and move on to the next brem.
73 if (secondaries->size() != 0) {
74 if (getEventInfo()->bremCandidateCount() == 1) {
75 track->SetTrackStatus(fKillTrackAndSecondaries);
76 G4RunManager::GetRunManager()->AbortEvent();
77 currentTrack_ = nullptr;
78 } else {
79 currentTrack_ = track;
80 track->SetTrackStatus(fSuspend);
82 trackInfo->tagBremCandidate(false);
83 }
84 }
85 return;
86 }
87
88 // If the brem photon doesn't undergo any reaction in the target, stop
89 // processing the rest of the event if the particle is exiting the
90 // target region.
91 if (secondaries->size() == 0) {
104 static auto recoil_physical_volume =
105 simcore::g4user::ptrretrieval::getPhysicalVolume("recoil_PV");
106 static auto world_physical_volume =
107 simcore::g4user::ptrretrieval::getPhysicalVolume("World_PV");
108
109 if (!recoil_physical_volume) {
110 ldmx_log(warn) << "Volume 'recoil_PV' not found in Geant4 volume store";
111 }
112 if (!world_physical_volume) {
113 ldmx_log(warn) << "Volume 'World_PV' not found in Geant4 volume store";
114 }
115
116 auto current_volume = track->GetNextVolume();
117 if (current_volume == recoil_physical_volume or
118 current_volume == world_physical_volume) {
119 if (getEventInfo()->bremCandidateCount() == 1) {
120 track->SetTrackStatus(fKillTrackAndSecondaries);
121 G4RunManager::GetRunManager()->AbortEvent();
122 currentTrack_ = nullptr;
123 } else {
124 currentTrack_ = track;
125 track->SetTrackStatus(fSuspend);
127 trackInfo->tagBremCandidate(false);
128 }
129 }
130 return;
131 } else {
132 // If the brem gamma interacts and produced secondaries, get the
133 // process used to create them.
134 G4String processName =
135 secondaries->at(0)->GetCreatorProcess()->GetProcessName();
136
137 // Only record the process that is being biased
138 if (!processName.contains(process_)) {
139 if (getEventInfo()->bremCandidateCount() == 1) {
140 track->SetTrackStatus(fKillTrackAndSecondaries);
141 G4RunManager::GetRunManager()->AbortEvent();
142 currentTrack_ = nullptr;
143 } else {
144 currentTrack_ = track;
145 track->SetTrackStatus(fSuspend);
147 trackInfo->tagBremCandidate(false);
148 }
149 return;
150 }
151
152 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
153 std::cout << "[ TargetProcessFilter ]: "
154 << G4EventManager::GetEventManager()
155 ->GetConstCurrentEvent()
156 ->GetEventID()
157 << " Brem photon produced " << secondaries->size()
158 << " particle via " << processName << " process." << std::endl;
159 }
160 trackInfo->tagBremCandidate(false);
161 trackInfo->setSaveFlag(true);
162 trackInfo->tagPNGamma();
164 }
165}
UserEventInformation * getEventInfo() const
Get a handle to the event information.
void decBremCandidateCount()
Decrease the number of brem candidates in an event.
static UserTrackInformation * get(const G4Track *track)
get

References currentTrack_, simcore::UserEventInformation::decBremCandidateCount(), simcore::UserTrackInformation::get(), simcore::UserAction::getEventInfo(), and process_.

Member Data Documentation

◆ currentTrack_

G4Track* biasing::TargetProcessFilter::currentTrack_ {nullptr}
private

Pointer to the current track being processed.

Definition at line 67 of file TargetProcessFilter.h.

67{nullptr};

Referenced by ClassifyNewTrack(), and stepping().

◆ process_

std::string biasing::TargetProcessFilter::process_ {""}
private

The process to bias.

Definition at line 70 of file TargetProcessFilter.h.

70{""};

Referenced by stepping(), and TargetProcessFilter().


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