LDMX Software
simcore::PhotonuclearTracker Class Reference

Tracks detailed information about photonuclear interactions. More...

#include <PhotonuclearTracker.h>

Public Member Functions

 PhotonuclearTracker (const std::string &name, framework::config::Parameters &parameters)
 Constructor.
 
virtual ~PhotonuclearTracker ()=default
 Destructor.
 
std::vector< TYPEgetTypes () override
 Get the types of actions this class handles.
 
void BeginOfEventAction (const G4Event *event) override
 Called at the beginning of each event.
 
void EndOfEventAction (const G4Event *event) override
 Called at the end of each event.
 
void stepping (const G4Step *step) override
 Called after each simulation step.
 
void PreUserTrackingAction (const G4Track *track) override
 Called before tracking a new track.
 
void PostUserTrackingAction (const G4Track *track) override
 Called after tracking a track.
 
std::vector< ldmx::PhotonuclearInteractiongetInteractions () const
 Get the collection of PN interactions for this event.
 
- Public Member Functions inherited from simcore::UserAction
 UserAction (const std::string &name, framework::config::Parameters &parameters)
 Constructor.
 
 DECLARE_FACTORY (UserAction, std::shared_ptr< UserAction >, const std::string &, framework::config::Parameters &)
 factory for user actions
 
virtual ~UserAction ()=default
 Destructor.
 
virtual void BeginOfRunAction (const G4Run *)
 Method called at the beginning of a run.
 
virtual void EndOfRunAction (const G4Run *)
 Method called at the end of a run.
 
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *, const G4ClassificationOfNewTrack &cl)
 Method called when a track is updated.
 
virtual void NewStage ()
 Method called at the beginning of a new stage.
 
virtual void PrepareNewEvent ()
 Method called at the beginning of a new event.
 

Static Public Member Functions

static PhotonuclearTrackerget ()
 Get the current PhotonuclearTracker instance.
 

Private Member Functions

void extractTargetInfo (const G4Step *step, int &Z, int &A, std::string &material)
 Extract target nucleus information from the current step.
 
ldmx::PhotonuclearInteraction::ParticleInfo createParticleInfo (const G4Track *track)
 Create a ParticleInfo struct from a G4Track.
 

Private Attributes

std::vector< ldmx::PhotonuclearInteractionpn_interactions_
 Collection of PN interactions in this event.
 
std::map< int, int > secondary_to_pn_index_
 Maps secondary track ID -> index in pn_interactions_ Used to quickly find which PN interaction a secondary belongs to.
 
std::map< int, int > descendant_ancestry_
 Maps descendant track ID -> immediate secondary track ID Propagates through the genealogy to trace back to the original PN secondary.
 

Static Private Attributes

static PhotonuclearTrackerinstance = nullptr
 Static pointer to the current instance.
 

Additional Inherited Members

- Protected Member Functions inherited from simcore::UserAction
UserEventInformationgetEventInfo () const
 Get a handle to the event information.
 
const std::map< int, ldmx::SimParticle > & getCurrentParticleMap () const
 Get the current particle map.
 
- Protected Attributes inherited from simcore::UserAction
std::string name_ {""}
 Name of the UserAction.
 
framework::config::Parameters parameters_
 The set of parameters used to configure this class.
 
mutable::framework::logging::logger the_log_
 the logging channel user actions can use ldmx_log with
 

Detailed Description

Tracks detailed information about photonuclear interactions.

This UserAction captures comprehensive details about photonuclear (PN) interactions during simulation, including:

  • Incident photon kinematics
  • Target nucleus information
  • All immediate secondary particles from the cascade
  • Mapping from cascade products to their final state descendants

The tracker operates by:

  1. Detecting PN interactions in the stepping action (when secondaries are created)
  2. Capturing all immediate secondaries and their G4 track IDs
  3. Tracking descendants through the event via ancestry propagation
  4. Building final state maps when tracks complete

The collected interactions are saved to the event bus as a vector: "PhotonuclearInteractions"

Definition at line 52 of file PhotonuclearTracker.h.

Constructor & Destructor Documentation

◆ PhotonuclearTracker()

simcore::PhotonuclearTracker::PhotonuclearTracker ( const std::string & name,
framework::config::Parameters & parameters )

Constructor.

Parameters
nameName of this user action
parametersConfiguration parameters

Definition at line 18 of file PhotonuclearTracker.cxx.

20 : UserAction(name, parameters) {
21 instance = this; // Set the static instance pointer
22}
static PhotonuclearTracker * instance
Static pointer to the current instance.
UserAction(const std::string &name, framework::config::Parameters &parameters)
Constructor.

References instance.

Member Function Documentation

◆ BeginOfEventAction()

void simcore::PhotonuclearTracker::BeginOfEventAction ( const G4Event * event)
overridevirtual

Called at the beginning of each event.

Clears the interaction collection and ancestry maps for the new event.

Parameters
eventGeant4 event object

Reimplemented from simcore::UserAction.

Definition at line 24 of file PhotonuclearTracker.cxx.

24 {
25 // Clear all data structures for the new event
26 pn_interactions_.clear();
29}
std::map< int, int > descendant_ancestry_
Maps descendant track ID -> immediate secondary track ID Propagates through the genealogy to trace ba...
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...
std::vector< ldmx::PhotonuclearInteraction > pn_interactions_
Collection of PN interactions in this event.

References descendant_ancestry_, pn_interactions_, and secondary_to_pn_index_.

◆ createParticleInfo()

ldmx::PhotonuclearInteraction::ParticleInfo simcore::PhotonuclearTracker::createParticleInfo ( const G4Track * track)
private

Create a ParticleInfo struct from a G4Track.

Parameters
trackGeant4 track
Returns
ParticleInfo with kinematic data

Definition at line 182 of file PhotonuclearTracker.cxx.

182 {
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}
Stores kinematic and identity information for a particle.

Referenced by stepping().

◆ EndOfEventAction()

void simcore::PhotonuclearTracker::EndOfEventAction ( const G4Event * event)
overridevirtual

Called at the end of each event.

Final cleanup and validation of collected data.

Parameters
eventGeant4 event object

Reimplemented from simcore::UserAction.

Definition at line 31 of file PhotonuclearTracker.cxx.

31 {
32 ldmx_log(debug) << "Event complete, recorded " << pn_interactions_.size()
33 << " PN interactions";
34}

References pn_interactions_.

◆ extractTargetInfo()

void simcore::PhotonuclearTracker::extractTargetInfo ( const G4Step * step,
int & Z,
int & A,
std::string & material )
private

Extract target nucleus information from the current step.

Parameters
stepCurrent Geant4 step
ZOutput: atomic number
AOutput: mass number
materialOutput: material name

Definition at line 168 of file PhotonuclearTracker.cxx.

169 {
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}

Referenced by stepping().

◆ get()

static PhotonuclearTracker * simcore::PhotonuclearTracker::get ( )
inlinestatic

Get the current PhotonuclearTracker instance.

Returns
Pointer to the current instance, or nullptr if not created

Definition at line 139 of file PhotonuclearTracker.h.

139{ return instance; }

References instance.

◆ getInteractions()

std::vector< ldmx::PhotonuclearInteraction > simcore::PhotonuclearTracker::getInteractions ( ) const
inline

Get the collection of PN interactions for this event.

Returns
Vector of PhotonuclearInteraction objects

Definition at line 130 of file PhotonuclearTracker.h.

130 {
131 return pn_interactions_;
132 }

References pn_interactions_.

◆ getTypes()

std::vector< TYPE > simcore::PhotonuclearTracker::getTypes ( )
inlineoverridevirtual

Get the types of actions this class handles.

Returns
Vector of action types (EVENT, STEPPING, TRACKING)

Implements simcore::UserAction.

Definition at line 73 of file PhotonuclearTracker.h.

73 {
74 return {TYPE::EVENT, TYPE::STEPPING, TYPE::TRACKING};
75 }

◆ PostUserTrackingAction()

void simcore::PhotonuclearTracker::PostUserTrackingAction ( const G4Track * track)
overridevirtual

Called after tracking a track.

Records final state particles when tracks complete, building the descendant map for PN interactions.

Parameters
trackCurrent Geant4 track

Reimplemented from simcore::UserAction.

Definition at line 137 of file PhotonuclearTracker.cxx.

137 {
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}

References descendant_ancestry_, pn_interactions_, and secondary_to_pn_index_.

◆ PreUserTrackingAction()

void simcore::PhotonuclearTracker::PreUserTrackingAction ( const G4Track * track)
overridevirtual

Called before tracking a new track.

Propagates ancestry information for descendants of PN secondaries.

Parameters
trackCurrent Geant4 track

Reimplemented from simcore::UserAction.

Definition at line 120 of file PhotonuclearTracker.cxx.

120 {
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}

References descendant_ancestry_.

◆ stepping()

void simcore::PhotonuclearTracker::stepping ( const G4Step * step)
overridevirtual

Called after each simulation step.

Detects PN interactions by examining secondaries created in the step. When a PN interaction is found, captures all details and initializes genealogy tracking.

Parameters
stepCurrent Geant4 step

Reimplemented from simcore::UserAction.

Definition at line 36 of file PhotonuclearTracker.cxx.

36 {
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}
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.
ldmx::PhotonuclearInteraction::ParticleInfo createParticleInfo(const G4Track *track)
Create a ParticleInfo struct from a G4Track.
void extractTargetInfo(const G4Step *step, int &Z, int &A, std::string &material)
Extract target nucleus information from the current step.

References ldmx::PhotonuclearInteraction::addImmediateSecondary(), createParticleInfo(), descendant_ancestry_, extractTargetInfo(), pn_interactions_, secondary_to_pn_index_, ldmx::PhotonuclearInteraction::setIncidentPhoton(), ldmx::PhotonuclearInteraction::setInteractionVolume(), ldmx::PhotonuclearInteraction::setProcessName(), and ldmx::PhotonuclearInteraction::setTarget().

Member Data Documentation

◆ descendant_ancestry_

std::map<int, int> simcore::PhotonuclearTracker::descendant_ancestry_
private

Maps descendant track ID -> immediate secondary track ID Propagates through the genealogy to trace back to the original PN secondary.

Definition at line 173 of file PhotonuclearTracker.h.

Referenced by BeginOfEventAction(), PostUserTrackingAction(), PreUserTrackingAction(), and stepping().

◆ instance

PhotonuclearTracker * simcore::PhotonuclearTracker::instance = nullptr
staticprivate

Static pointer to the current instance.

Definition at line 176 of file PhotonuclearTracker.h.

Referenced by get(), and PhotonuclearTracker().

◆ pn_interactions_

std::vector<ldmx::PhotonuclearInteraction> simcore::PhotonuclearTracker::pn_interactions_
private

Collection of PN interactions in this event.

Definition at line 164 of file PhotonuclearTracker.h.

Referenced by BeginOfEventAction(), EndOfEventAction(), getInteractions(), PostUserTrackingAction(), and stepping().

◆ secondary_to_pn_index_

std::map<int, int> simcore::PhotonuclearTracker::secondary_to_pn_index_
private

Maps secondary track ID -> index in pn_interactions_ Used to quickly find which PN interaction a secondary belongs to.

Definition at line 168 of file PhotonuclearTracker.h.

Referenced by BeginOfEventAction(), PostUserTrackingAction(), and stepping().


The documentation for this class was generated from the following files: