LDMX Software
Public Member Functions | Private Attributes | List of all members
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 = ::simcore::Factory< UserAction, std::shared_ptr< UserAction >, const std::string &, framework::config::Parameters & >
 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.
 

Detailed Description

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

Definition at line 26 of file TargetProcessFilter.h.

Constructor & Destructor Documentation

◆ TargetProcessFilter()

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

Class constructor.

Definition at line 25 of file TargetProcessFilter.cxx.

27 : simcore::UserAction(name, parameters) {
28 process_ = parameters.getParameter<std::string>("process");
29}
std::string process_
The process to bias.
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
Interface that defines a user action.
Definition UserAction.h:42

References framework::config::Parameters::getParameter(), and 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 31 of file TargetProcessFilter.cxx.

32 {
33 if (track == currentTrack_) {
34 currentTrack_ = nullptr;
35 // std::cout << "[ TargetBremFilter ]: Pushing track to waiting stack." <<
36 // std::endl;
37 return fWaiting;
38 }
39
40 // Use current classification by default so values from other plugins are not
41 // overridden.
42 G4ClassificationOfNewTrack classification = currentTrackClass;
43
44 return classification;
45}
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 151 of file TargetProcessFilter.cxx.

151{}

◆ getTypes()

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

Retrieve the type of actions this class defines.

Implements simcore::UserAction.

Definition at line 59 of file TargetProcessFilter.h.

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

◆ 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 47 of file TargetProcessFilter.cxx.

47 {
48 // Get the track associated with this step.
49 auto track{step->GetTrack()};
50
51 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted())
52 return;
53
54 // Get the track info and check if this track is a brem candidate
55 auto trackInfo{simcore::UserTrackInformation::get(track)};
56 if ((trackInfo != nullptr) && !trackInfo->isBremCandidate()) return;
57
58 // Get the particles daughters.
59 auto secondaries{step->GetSecondary()};
60
61 // Get the region the particle is currently in. Continue processing
62 // the particle only if it's in the target region.
63 if (auto region{
64 track->GetVolume()->GetLogicalVolume()->GetRegion()->GetName()};
65 region.compareTo("target") != 0) {
66 // If secondaries were produced outside of the volume of interest,
67 // and there aren't additional brems to process, abort the event.
68 // Otherwise, suspend the track and move on to the next brem.
69 if (secondaries->size() != 0) {
70 if (getEventInfo()->bremCandidateCount() == 1) {
71 track->SetTrackStatus(fKillTrackAndSecondaries);
72 G4RunManager::GetRunManager()->AbortEvent();
73 currentTrack_ = nullptr;
74 } else {
75 currentTrack_ = track;
76 track->SetTrackStatus(fSuspend);
78 trackInfo->tagBremCandidate(false);
79 }
80 }
81 return;
82 }
83
84 // If the brem photon doesn't undergo any reaction in the target, stop
85 // processing the rest of the event if the particle is exiting the
86 // target region.
87 if (secondaries->size() == 0) {
100 if (auto volume{track->GetNextVolume()->GetName()};
101 volume.compareTo("recoil_PV") == 0 or
102 volume.compareTo("World_PV") == 0) {
103 if (getEventInfo()->bremCandidateCount() == 1) {
104 track->SetTrackStatus(fKillTrackAndSecondaries);
105 G4RunManager::GetRunManager()->AbortEvent();
106 currentTrack_ = nullptr;
107 } else {
108 currentTrack_ = track;
109 track->SetTrackStatus(fSuspend);
111 trackInfo->tagBremCandidate(false);
112 }
113 }
114 return;
115 } else {
116 // If the brem gamma interacts and produced secondaries, get the
117 // process used to create them.
118 G4String processName =
119 secondaries->at(0)->GetCreatorProcess()->GetProcessName();
120
121 // Only record the process that is being biased
122 if (!processName.contains(process_)) {
123 if (getEventInfo()->bremCandidateCount() == 1) {
124 track->SetTrackStatus(fKillTrackAndSecondaries);
125 G4RunManager::GetRunManager()->AbortEvent();
126 currentTrack_ = nullptr;
127 } else {
128 currentTrack_ = track;
129 track->SetTrackStatus(fSuspend);
131 trackInfo->tagBremCandidate(false);
132 }
133 return;
134 }
135
136 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
137 std::cout << "[ TargetProcessFilter ]: "
138 << G4EventManager::GetEventManager()
139 ->GetConstCurrentEvent()
140 ->GetEventID()
141 << " Brem photon produced " << secondaries->size()
142 << " particle via " << processName << " process." << std::endl;
143 }
144 trackInfo->tagBremCandidate(false);
145 trackInfo->setSaveFlag(true);
146 trackInfo->tagPNGamma();
148 }
149}
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 66 of file TargetProcessFilter.h.

66{nullptr};

Referenced by ClassifyNewTrack(), and stepping().

◆ process_

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

The process to bias.

Definition at line 69 of file TargetProcessFilter.h.

69{""};

Referenced by stepping(), and TargetProcessFilter().


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