LDMX Software
TargetProcessFilter.cxx
1
2#include "Biasing/TargetProcessFilter.h"
3
4/*~~~~~~~~~~~~*/
5/* Geant4 */
6/*~~~~~~~~~~~~*/
7#include "G4Event.hh"
8#include "G4EventManager.hh"
9#include "G4RunManager.hh"
10#include "G4Step.hh"
11#include "G4Track.hh"
12
13/*~~~~~~~~~~~~~*/
14/* Biasing */
15/*~~~~~~~~~~~~~*/
16#include "Biasing/TargetBremFilter.h"
17
18/*~~~~~~~~~~~~~*/
19/* SimCore */
20/*~~~~~~~~~~~~~*/
21#include "SimCore/UserTrackInformation.h"
22
23namespace biasing {
24
26 const std::string& name, framework::config::Parameters& parameters)
27 : simcore::UserAction(name, parameters) {
28 process_ = parameters.getParameter<std::string>("process");
29}
30
31G4ClassificationOfNewTrack TargetProcessFilter::ClassifyNewTrack(
32 const G4Track* track, const G4ClassificationOfNewTrack& currentTrackClass) {
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}
46
47void TargetProcessFilter::stepping(const G4Step* step) {
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}
150
152} // namespace biasing
153
154DECLARE_ACTION(biasing, TargetProcessFilter)
G4Track * currentTrack_
Pointer to the current track being processed.
G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack, const G4ClassificationOfNewTrack &currentTrackClass) override
Classify a new track which postpones track processing.
std::string process_
The process to bias.
TargetProcessFilter(const std::string &name, framework::config::Parameters &parameters)
Class constructor.
void EndOfEventAction(const G4Event *) override
End of event action.
void stepping(const G4Step *step) override
Implementmthe stepping action which performs the target volume biasing.
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
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