LDMX Software
TargetProcessFilter.cxx
1/*~~~~~~~~~~~~~*/
2/* Biasing */
3/*~~~~~~~~~~~~~*/
4#include "Biasing/TargetProcessFilter.h"
5
6/*~~~~~~~~~~~~*/
7/* Geant4 */
8/*~~~~~~~~~~~~*/
9#include "G4Event.hh"
10#include "G4EventManager.hh"
11#include "G4RunManager.hh"
12#include "G4Step.hh"
13#include "G4Track.hh"
14
15/*~~~~~~~~~~~~~*/
16/* SimCore */
17/*~~~~~~~~~~~~~*/
18#include "SimCore/G4User/PtrRetrieval.h"
19#include "SimCore/UserTrackInformation.h"
20
21namespace biasing {
22
24 const std::string& name, framework::config::Parameters& parameters)
25 : simcore::UserAction(name, parameters) {
26 process_ = parameters.getParameter<std::string>("process");
27}
28
29G4ClassificationOfNewTrack TargetProcessFilter::ClassifyNewTrack(
30 const G4Track* track, const G4ClassificationOfNewTrack& currentTrackClass) {
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}
44
45void TargetProcessFilter::stepping(const G4Step* step) {
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}
166
168} // namespace biasing
169
170DECLARE_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:29
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