LDMX Software
TargetDarkBremFilter.cxx
1/*
2 * @file TargetDarkBremFilter.cxx
3 * @class TargetDarkBremFilter
4 * @brief Class defining a UserActionPlugin that allows a user to filter out
5 * events that don't result in a dark brem within a given volume
6 * @author Michael Revering, University of Minnesota
7 * @author Tom Eichlersmith, University of Minnesota
8 */
9
11
12namespace biasing {
13
15 const std::string& name, framework::config::Parameters& parameters)
16 : simcore::UserAction(name, parameters) {
17 threshold_ = parameters.get<double>("threshold");
18}
19
21 found_aprime_ = false;
22}
23
24void TargetDarkBremFilter::stepping(const G4Step* step) {
25 // don't process if event is already aborted
26 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted())
27 return;
28
29 // Get the track associated with this step.
30 auto track{step->GetTrack()};
31
32 // Leave if track is not an electron
33 auto particle_def{track->GetParticleDefinition()};
34 if (particle_def != G4Electron::Electron()) {
35 if (particle_def == G4APrime::APrime() and
36 track->GetCurrentStepNumber() == 1) {
49 if (isOutsideTargetRegion(track->GetLogicalVolumeAtVertex()))
50 AbortEvent("A' was not created within target.");
51 // A' was found and originated in correct region
52 found_aprime_ = true;
53 } // first step of A'
54 return;
55 }
56
57 // Leave if track is not primary
58 if (track->GetParentID() != 0) return;
59
60 // Leave if track is not in target region
61 if (isOutsideTargetRegion(track->GetVolume())) return;
62
63 if (isOutsideTargetRegion(track->GetNextVolume()) // leaving target region
64 or
65 track->GetTrackStatus() == fStopAndKill // stopping within target region
66 or track->GetKineticEnergy() == 0. // stopping within target region
67 ) {
68 // Get the electron secondries
69 const std::vector<G4Track*>* secondaries{step->GetSecondary()};
70 if (not secondaries or secondaries->size() == 0) {
71 AbortEvent("Primary electron did not create secondaries.");
72 return;
73 } else {
74 // check secondaries to see if we made a dark brem
75 for (auto& secondary_track : *secondaries) {
76 if (secondary_track->GetParticleDefinition() == G4APrime::APrime()) {
77 // we found an A', woo-hoo!
78
79 if (secondary_track->GetTotalEnergy() < threshold_) {
81 "A' was not created with total energy above input threshold.");
82 return;
83 }
84
85 if (isOutsideTargetRegion(secondary_track->GetVolume())) {
86 AbortEvent("A' was not created within target.");
87 return;
88 }
89
90 /* debug message
91 std::cout << "[ TargetDarkBremFilter ] : "
92 << "(" <<
93 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID()
94 << ") "
95 << "Found an A' weighted " << secondary_track->GetWeight()
96 << " originating in " <<
97 secondary_track->GetVolume()->GetLogicalVolume()->GetName()
98 << " with energy " << secondary_track->GetTotalEnergy() << " MeV."
99 << std::endl;
100 */
101
102 // we found a good A', so we can leave
103 return;
104 } // this secondary is an A'
105 } // loop through secondaries
106 } // are there secondaries to loop through
107
108 // got here without finding A'
109 AbortEvent("Primary electron did not create A'.");
110 return;
111 } // should check secondaries of primary
112
113 // get here if we are the primary in the target region,
114 // but shouldn't check the secondaries yet
115 return;
116}
117
118void TargetDarkBremFilter::EndOfEventAction(const G4Event* event) {
119 if (not event->IsAborted() and not found_aprime_) {
120 AbortEvent("Did not find an A' after entire event was simulated.");
121 }
122}
123
124void TargetDarkBremFilter::AbortEvent(const std::string& reason) const {
125 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
126 std::cout << "[ TargetDarkBremFilter ]: " << "("
127 << G4EventManager::GetEventManager()
128 ->GetConstCurrentEvent()
129 ->GetEventID()
130 << ") " << reason << " Aborting event." << std::endl;
131 }
132 G4RunManager::GetRunManager()->AbortEvent();
133 return;
134}
135} // namespace biasing
136
#define DECLARE_ACTION(CLASS)
register a new UserAction with its factory
Definition UserAction.h:206
This class is meant to filter for events that produce a dark brem matching originating in the target ...
void stepping(const G4Step *step) override
Looking for A' while primary is stepping.
TargetDarkBremFilter(const std::string &name, framework::config::Parameters &parameters)
Class constructor.
double threshold_
Minimum energy [MeV] that the A' should have to keep the event.
bool isOutsideTargetRegion(const G4VPhysicalVolume *vol) const
Check if the volume is outside the target region.
void BeginOfEventAction(const G4Event *e) override
Reset flag signaling finding of A' to false.
void EndOfEventAction(const G4Event *event) override
Check flag signaling finding of A', if false, abort the event so it isn't saved.
void AbortEvent(const std::string &reason) const
Helper to abort an event with a message.
bool found_aprime_
flag to signal that we saw an A'
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...