LDMX Software
PhotonuclearTracker.cxx
2
3#include "G4Element.hh"
4#include "G4EventManager.hh"
5#include "G4Material.hh"
6#include "G4RunManager.hh"
7#include "G4Step.hh"
8#include "G4StepPoint.hh"
9#include "G4Track.hh"
10#include "G4VProcess.hh"
12#include "SimCore/G4User/UserEventInformation.h"
13
14namespace simcore {
15
16PhotonuclearTracker* PhotonuclearTracker::instance = nullptr;
17
19 const std::string& name, framework::config::Parameters& parameters)
20 : UserAction(name, parameters) {
21 instance = this; // Set the static instance pointer
22}
23
24void PhotonuclearTracker::BeginOfEventAction(const G4Event* event) {
25 // Clear all data structures for the new event
26 pn_interactions_.clear();
29}
30
31void PhotonuclearTracker::EndOfEventAction(const G4Event* event) {
32 ldmx_log(debug) << "Event complete, recorded " << pn_interactions_.size()
33 << " PN interactions";
34}
35
36void PhotonuclearTracker::stepping(const G4Step* step) {
37 // Check if any secondaries were created in this step
38 const std::vector<const G4Track*>* secondaries =
39 step->GetSecondaryInCurrentStep();
40 if (!secondaries || secondaries->empty()) return;
41
42 // Check if any secondary was created by a photonNuclear process
43 bool is_pn_interaction = false;
44 for (const G4Track* secondary : *secondaries) {
45 const G4VProcess* creator = secondary->GetCreatorProcess();
46 if (creator) {
47 const G4String& creator_name = creator->GetProcessName();
48 if (creator_name.contains("photonNuclear")) {
49 is_pn_interaction = true;
50 break;
51 }
52 }
53 }
54
55 if (!is_pn_interaction) return;
56
57 // We have a PN interaction! Create a new PhotonuclearInteraction object
58 ldmx::PhotonuclearInteraction pn_interaction;
59
60 // Get the incident photon information from the track that took the step
61 const G4Track* incident_track = step->GetTrack();
62 const G4StepPoint* pre_step = step->GetPreStepPoint();
63
64 auto incident_pos = pre_step->GetPosition();
65 auto incident_momentum = pre_step->GetMomentum();
66
67 pn_interaction.setIncidentPhoton(
68 incident_track->GetTrackID(),
69 incident_track->GetDefinition()->GetPDGEncoding(),
70 pre_step->GetKineticEnergy() +
71 incident_track->GetDefinition()->GetPDGMass(),
72 incident_momentum.x(), incident_momentum.y(), incident_momentum.z(),
73 incident_pos.x(), incident_pos.y(), incident_pos.z(),
74 pre_step->GetGlobalTime());
75
76 // Get target nucleus information
77 int target_z, target_a;
78 std::string target_material;
79 extractTargetInfo(step, target_z, target_a, target_material);
80 pn_interaction.setTarget(target_z, target_a, target_material);
81
82 // Set interaction metadata
83 pn_interaction.setInteractionVolume(incident_track->GetVolume()->GetName());
84 pn_interaction.setProcessName("photonNuclear");
85
86 // Record all immediate secondaries
87 int pn_index = pn_interactions_.size();
88
89 // Set target Z/A in UserEventInformation for the first PN interaction
90 if (pn_index == 0) {
91 auto event_info = static_cast<UserEventInformation*>(
92 G4EventManager::GetEventManager()->GetUserInformation());
93 if (event_info) {
94 event_info->setPNTargetZ(target_z);
95 event_info->setPNTargetA(target_a);
96 }
97 }
98 for (const G4Track* secondary : *secondaries) {
99 auto particle_info = createParticleInfo(secondary);
100 pn_interaction.addImmediateSecondary(particle_info);
101
102 // Map this secondary to the current PN interaction
103 int secondary_id = secondary->GetTrackID();
104 secondary_to_pn_index_[secondary_id] = pn_index;
105 // Initialize ancestry: immediate secondaries map to themselves
106 descendant_ancestry_[secondary_id] = secondary_id;
107
108 ldmx_log(debug) << " Immediate secondary: track_id=" << secondary_id
109 << " pdg=" << particle_info.pdg_id_
110 << " E=" << particle_info.energy_ << " MeV";
111 }
112
113 // Add the completed PN interaction to the collection
114 pn_interactions_.push_back(pn_interaction);
115
116 ldmx_log(debug) << "Detected PN interaction #" << pn_index << " with "
117 << secondaries->size() << " secondaries";
118}
119
121 // Propagate ancestry information to new tracks
122 int track_id = track->GetTrackID();
123 int parent_id = track->GetParentID();
124
125 // If the parent is a descendant of a PN secondary, this track is too
126 auto parent_it = descendant_ancestry_.find(parent_id);
127 if (parent_it != descendant_ancestry_.end()) {
128 // Propagate the immediate secondary ID to this descendant
129 descendant_ancestry_[track_id] = parent_it->second;
130
131 ldmx_log(debug) << " Track " << track_id
132 << " is descendant of PN secondary " << parent_it->second
133 << " (parent=" << parent_id << ")";
134 }
135}
136
138 int track_id = track->GetTrackID();
139
140 // Check if this track is a descendant of a PN secondary
141 auto ancestry_it = descendant_ancestry_.find(track_id);
142 if (ancestry_it == descendant_ancestry_.end()) {
143 return; // Not related to any PN interaction
144 }
145
146 // This track is a descendant! Find which immediate secondary it came from
147 int immediate_secondary_id = ancestry_it->second;
148
149 // Find which PN interaction this secondary belongs to
150 auto pn_index_it = secondary_to_pn_index_.find(immediate_secondary_id);
151 if (pn_index_it == secondary_to_pn_index_.end()) {
152 ldmx_log(warn) << "Inconsistent ancestry mapping for track " << track_id;
153 return;
154 }
155
156 int pn_index = pn_index_it->second;
157
158 // Record ALL descendants (intermediate + final state) for full cascade
159 // genealogy This provides information beyond SimParticles which only has
160 // final state
161 pn_interactions_[pn_index].addDescendant(immediate_secondary_id, track_id);
162
163 ldmx_log(debug) << " Recording descendant: track " << track_id
164 << " -> secondary " << immediate_secondary_id
165 << " in PN interaction " << pn_index;
166}
167
168void PhotonuclearTracker::extractTargetInfo(const G4Step* step, int& Z, int& A,
169 std::string& material) {
170 const G4Material* mat = step->GetPreStepPoint()->GetMaterial();
171 material = mat->GetName();
172
173 // Get the element - for compound materials, we take the first element
174 // (This is a simplification; the actual nucleus struck may vary)
175 const G4Element* element = mat->GetElement(0);
176 Z = static_cast<int>(element->GetZ());
177 A = static_cast<int>(element->GetA() /
178 (CLHEP::g / CLHEP::mole)); // Convert to mass number
179}
180
184
185 info.track_id_ = track->GetTrackID();
186 info.pdg_id_ = track->GetDefinition()->GetPDGEncoding();
187
188 auto position = track->GetPosition();
189 info.x_ = position.x();
190 info.y_ = position.y();
191 info.z_ = position.z();
192
193 auto momentum = track->GetMomentum();
194 info.px_ = momentum.x();
195 info.py_ = momentum.y();
196 info.pz_ = momentum.z();
197
198 info.energy_ =
199 track->GetKineticEnergy() + track->GetDefinition()->GetPDGMass();
200 info.time_ = track->GetGlobalTime();
201
202 return info;
203}
204
205} // namespace simcore
206
UserAction for tracking detailed photonuclear interaction information.
Class which implements the user tracking action.
#define DECLARE_ACTION(CLASS)
register a new UserAction with its factory
Definition UserAction.h:206
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
Stores detailed information about a photonuclear interaction.
void setProcessName(const std::string &process)
Set the process name.
void setInteractionVolume(const std::string &volume)
Set the volume where the interaction occurred.
void addImmediateSecondary(const ParticleInfo &particle)
Add an immediate secondary particle from the PN interaction.
void setIncidentPhoton(int track_id, int pdg_id, double energy, double px, double py, double pz, double x, double y, double z, double time)
Set the incident photon information.
void setTarget(int Z, int A, const std::string &material)
Set the target nucleus information.
Tracks detailed information about photonuclear interactions.
static PhotonuclearTracker * instance
Static pointer to the current instance.
void BeginOfEventAction(const G4Event *event) override
Called at the beginning of each event.
ldmx::PhotonuclearInteraction::ParticleInfo createParticleInfo(const G4Track *track)
Create a ParticleInfo struct from a G4Track.
std::map< int, int > descendant_ancestry_
Maps descendant track ID -> immediate secondary track ID Propagates through the genealogy to trace ba...
void PreUserTrackingAction(const G4Track *track) override
Called before tracking a new track.
PhotonuclearTracker(const std::string &name, framework::config::Parameters &parameters)
Constructor.
std::map< int, int > secondary_to_pn_index_
Maps secondary track ID -> index in pn_interactions_ Used to quickly find which PN interaction a seco...
void extractTargetInfo(const G4Step *step, int &Z, int &A, std::string &material)
Extract target nucleus information from the current step.
void PostUserTrackingAction(const G4Track *track) override
Called after tracking a track.
void EndOfEventAction(const G4Event *event) override
Called at the end of each event.
void stepping(const G4Step *step) override
Called after each simulation step.
std::vector< ldmx::PhotonuclearInteraction > pn_interactions_
Collection of PN interactions in this event.
Interface that defines a user action.
Definition UserAction.h:47
Encapsulates user defined information associated with a Geant4 event.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...
Stores kinematic and identity information for a particle.