LDMX Software
PhotonuclearInteraction.h
1#ifndef SIMCORE_EVENT_PHOTONUCLEARINTERACTION_H
2#define SIMCORE_EVENT_PHOTONUCLEARINTERACTION_H
3
4/*~~~~~~~~~~*/
5/* ROOT */
6/*~~~~~~~~~~*/
7#include "TObject.h"
8
9/*~~~~~~~~~~~~~~~~*/
10/* C++ StdLib */
11/*~~~~~~~~~~~~~~~~*/
12#include <map>
13#include <string>
14#include <vector>
15
16namespace ldmx {
17
33 public:
38 struct ParticleInfo {
40 int track_id_{0};
41
43 int pdg_id_{0};
44
46 double energy_{0.0};
47
49 double px_{0.0};
50 double py_{0.0};
51 double pz_{0.0};
52
54 double x_{0.0};
55 double y_{0.0};
56 double z_{0.0};
57
59 double time_{0.0};
60
62 ParticleInfo() = default;
63
65 virtual ~ParticleInfo() = default;
66
69 };
70
73
75 virtual ~PhotonuclearInteraction() = default;
76
78 void clear();
79
93 void setIncidentPhoton(int track_id, int pdg_id, double energy, double px,
94 double py, double pz, double x, double y, double z,
95 double time);
96
103
111 void setTarget(int Z, int A, const std::string& material);
112
118 int getTargetZ() const { return target_z_; }
119
125 int getTargetA() const { return target_a_; }
126
132 std::string getTargetMaterial() const { return target_material_; }
133
139 void addImmediateSecondary(const ParticleInfo& particle);
140
146 std::vector<ParticleInfo> getImmediateSecondaries() const {
148 }
149
156 return immediate_secondaries_.size();
157 }
158
168 void addDescendant(int immediate_secondary_id, int final_state_track_id);
169
177 std::map<int, std::vector<int>> getDescendantMap() const {
178 return descendant_map_;
179 }
180
187 std::vector<int> getDescendants(int immediate_secondary_id) const;
188
194 void setInteractionVolume(const std::string& volume) {
195 interaction_volume_ = volume;
196 }
197
203 std::string getInteractionVolume() const { return interaction_volume_; }
204
210 void setProcessName(const std::string& process) { process_name_ = process; }
211
217 std::string getProcessName() const { return process_name_; }
218
222 void print() const;
223
227 friend std::ostream& operator<<(std::ostream& o,
228 const PhotonuclearInteraction& pn);
229
230 private:
233
235 int target_z_{0};
236
238 int target_a_{0};
239
241 std::string target_material_{""};
242
244 std::vector<ParticleInfo> immediate_secondaries_;
245
247 std::map<int, std::vector<int>> descendant_map_;
248
250 std::string interaction_volume_{""};
251
253 std::string process_name_{""};
254
257
258}; // PhotonuclearInteraction
259
260} // namespace ldmx
261
262#endif // SIMCORE_EVENT_PHOTONUCLEARINTERACTION_H
Stores detailed information about a photonuclear interaction.
int getTargetZ() const
Get the target atomic number.
std::string target_material_
Target material name.
virtual ~PhotonuclearInteraction()=default
Destructor.
std::vector< ParticleInfo > getImmediateSecondaries() const
Get all immediate secondary particles.
friend std::ostream & operator<<(std::ostream &o, const PhotonuclearInteraction &pn)
Stream operator for printing.
void setProcessName(const std::string &process)
Set the process name.
void setInteractionVolume(const std::string &volume)
Set the volume where the interaction occurred.
void addDescendant(int immediate_secondary_id, int final_state_track_id)
Add a final state descendant track for an immediate secondary.
void print() const
Print information about this interaction.
ParticleInfo incident_photon_
Incident photon information.
void clear()
Clear the interaction data.
std::vector< int > getDescendants(int immediate_secondary_id) const
Get descendants for a specific immediate secondary.
void addImmediateSecondary(const ParticleInfo &particle)
Add an immediate secondary particle from the PN interaction.
std::string process_name_
Process name (e.g., "photonNuclear")
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.
std::map< int, std::vector< int > > getDescendantMap() const
Get the descendant map.
std::string interaction_volume_
Volume where interaction occurred.
std::string getTargetMaterial() const
Get the target material name.
std::string getProcessName() const
Get the process name.
PhotonuclearInteraction()=default
Constructor.
void setTarget(int Z, int A, const std::string &material)
Set the target nucleus information.
int getTargetA() const
Get the target mass number.
std::vector< ParticleInfo > immediate_secondaries_
Immediate secondary particles from the cascade.
std::map< int, std::vector< int > > descendant_map_
Maps immediate secondary track ID -> final state track IDs.
std::string getInteractionVolume() const
Get the interaction volume name.
int getNumImmediateSecondaries() const
Get the number of immediate secondary particles.
ClassDef(PhotonuclearInteraction, 1)
ROOT dictionary generation.
ParticleInfo getIncidentPhoton() const
Get the incident photon information.
Stores kinematic and identity information for a particle.
ParticleInfo()=default
Default constructor.
ClassDef(ParticleInfo, 1)
ROOT dictionary generation.
virtual ~ParticleInfo()=default
Virtual destructor (required by ClassDef)