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
12#include "G4DarkBreM/G4APrime.h" //checking if particles match A'
13#include "G4Electron.hh" //to check if track is electron
14#include "SimCore/UserTrackInformation.h" //make sure A' is saved
15
16namespace biasing {
17
19 const std::string& name, framework::config::Parameters& parameters)
20 : simcore::UserAction(name, parameters) {
21 threshold_ = parameters.getParameter<double>("threshold");
22}
23
25 found_aprime_ = false;
26}
27
28void TargetDarkBremFilter::stepping(const G4Step* step) {
29 // don't process if event is already aborted
30 if (G4EventManager::GetEventManager()->GetConstCurrentEvent()->IsAborted())
31 return;
32
33 // Get the track associated with this step.
34 auto track{step->GetTrack()};
35
36 // Leave if track is not an electron
37 auto particle_def{track->GetParticleDefinition()};
38 if (particle_def != G4Electron::Electron()) {
39 if (particle_def == G4APrime::APrime() and
40 track->GetCurrentStepNumber() == 1) {
53 if (isOutsideTargetRegion(track->GetLogicalVolumeAtVertex()))
54 AbortEvent("A' was not created within target.");
55 // A' was found and originated in correct region
56 found_aprime_ = true;
57 } // first step of A'
58 return;
59 }
60
61 // Leave if track is not primary
62 if (track->GetParentID() != 0) return;
63
64 // Leave if track is not in target region
65 if (isOutsideTargetRegion(track->GetVolume())) return;
66
67 if (isOutsideTargetRegion(track->GetNextVolume()) // leaving target region
68 or
69 track->GetTrackStatus() == fStopAndKill // stopping within target region
70 or track->GetKineticEnergy() == 0. // stopping within target region
71 ) {
72 // Get the electron secondries
73 const std::vector<G4Track*>* secondaries{step->GetSecondary()};
74 if (not secondaries or secondaries->size() == 0) {
75 AbortEvent("Primary electron did not create secondaries.");
76 return;
77 } else {
78 // check secondaries to see if we made a dark brem
79 for (auto& secondary_track : *secondaries) {
80 if (secondary_track->GetParticleDefinition() == G4APrime::APrime()) {
81 // we found an A', woo-hoo!
82
83 if (secondary_track->GetTotalEnergy() < threshold_) {
85 "A' was not created with total energy above input threshold.");
86 return;
87 }
88
89 if (isOutsideTargetRegion(secondary_track->GetVolume())) {
90 AbortEvent("A' was not created within target.");
91 return;
92 }
93
94 /* debug message
95 std::cout << "[ TargetDarkBremFilter ] : "
96 << "(" <<
97 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID()
98 << ") "
99 << "Found an A' weighted " << secondary_track->GetWeight()
100 << " originating in " <<
101 secondary_track->GetVolume()->GetLogicalVolume()->GetName()
102 << " with energy " << secondary_track->GetTotalEnergy() << " MeV."
103 << std::endl;
104 */
105
106 // we found a good A', so we can leave
107 return;
108 } // this secondary is an A'
109 } // loop through secondaries
110 } // are there secondaries to loop through
111
112 // got here without finding A'
113 AbortEvent("Primary electron did not create A'.");
114 return;
115 } // should check secondaries of primary
116
117 // get here if we are the primary in the target region,
118 // but shouldn't check the secondaries yet
119 return;
120}
121
122void TargetDarkBremFilter::EndOfEventAction(const G4Event* event) {
123 if (not event->IsAborted() and not found_aprime_) {
124 AbortEvent("Did not find an A' after entire event was simulated.");
125 }
126}
127
128void TargetDarkBremFilter::AbortEvent(const std::string& reason) const {
129 if (G4RunManager::GetRunManager()->GetVerboseLevel() > 1) {
130 std::cout << "[ TargetDarkBremFilter ]: "
131 << "("
132 << G4EventManager::GetEventManager()
133 ->GetConstCurrentEvent()
134 ->GetEventID()
135 << ") " << reason << " Aborting event." << std::endl;
136 }
137 G4RunManager::GetRunManager()->AbortEvent();
138 return;
139}
140} // namespace biasing
141
142DECLARE_ACTION(biasing, TargetDarkBremFilter)
Class defining a UserActionPlugin that allows a user to filter out events that don't result in a dark...
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 EndOfEventAction(const G4Event *event) override
Check flag signaling finding of A', if false, abort the event so it isn't saved.
void BeginOfEventAction(const G4Event *e) override
Reset flag signaling finding of A' to false.
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:27
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89