LDMX Software
HepMCEvent.cxx
2
3namespace simcore {
4namespace hepmc {
5
6HepMCEvent::HepMCEvent(std::shared_ptr<HepMC3::GenEvent> event)
7 : event_(event) {
8 if (!event_) {
9 EXCEPTION_RAISE("NullEvent",
10 "Attempted to create HepMCEvent with null GenEvent");
11 }
12
13 // Extract vertex information from the signal process vertex
14 // or the first vertex if no signal process is defined
15 std::shared_ptr<HepMC3::GenVertex> vertex;
16
17 // Try to get the signal process vertex first
18 if (event_->vertices().size() > 0) {
19 vertex = event_->vertices()[0];
20 }
21
22 if (vertex) {
23 // HepMC3 uses mm for length and mm/c for time by default
24 const HepMC3::FourVector& pos = vertex->position();
25 vtx_[0] = pos.x(); // mm
26 vtx_[1] = pos.y(); // mm
27 vtx_[2] = pos.z(); // mm
28 vtxt_ = pos.t(); // mm/c, needs conversion to ns for Geant4
29
30 // Convert time from mm/c to ns: t[ns] = t[mm/c] / c[mm/ns]
31 // c = 299.792458 mm/ns
32 vtxt_ = vtxt_ / 299.792458;
33 } else {
34 vtx_[0] = 0.0;
35 vtx_[1] = 0.0;
36 vtx_[2] = 0.0;
37 vtxt_ = 0.0;
38 }
39}
40
44 }
45 return particles_.size();
46}
47
49 // HepMC3 events can have multiple weights, we take the first one
50 const std::vector<double>& weights = event_->weights();
51 return weights.empty() ? 1.0 : weights[0];
52}
53
54const double* HepMCEvent::getVertex() const { return vtx_; }
55
56double HepMCEvent::getVertexTime() const { return vtxt_; }
57
58const std::vector<std::unique_ptr<HepMCParticle>>& HepMCEvent::getParticles()
59 const {
62 }
63 return particles_;
64}
65
66std::shared_ptr<HepMC3::GenEvent> HepMCEvent::getGenEvent() const {
67 return event_;
68}
69
71 particles_.clear();
72
73 // Get all particles from the event
74 for (const auto& particle : event_->particles()) {
75 // In HepMC3, status code 1 typically means final state particle
76 // Status codes vary by generator, but generally:
77 // status = 1: final state particle
78 // status = 2: intermediate/decayed particle
79 // status > 2: generator-specific codes
80
81 // We include particles with status code 1 (final state)
82 if (particle->status() == 1) {
83 particles_.push_back(std::make_unique<HepMCParticle>(particle));
84 }
85 }
86
88}
89
90} // namespace hepmc
91} // namespace simcore
Class defining a HepMC event with a list of particles.
double vtxt_
Vertex time (in ns)
Definition HepMCEvent.h:101
const std::vector< std::unique_ptr< HepMCParticle > > & getParticles() const
Get the list of final state particles in the event.
bool particles_extracted_
Flag to indicate if particles have been extracted.
Definition HepMCEvent.h:111
const double * getVertex() const
Get the vertex location (in mm, as expected by Geant4).
HepMCEvent(std::shared_ptr< HepMC3::GenEvent > event)
Class constructor.
Definition HepMCEvent.cxx:6
void extractParticles() const
Extract final state particles from the HepMC event.
double getEventWeight() const
Get the event weight.
double vtx_[3]
Vertex location (in mm)
Definition HepMCEvent.h:96
double getVertexTime() const
Get the vertex time (in ns, as expected by Geant4).
std::vector< std::unique_ptr< HepMCParticle > > particles_
The list of final state particles to be tracked.
Definition HepMCEvent.h:106
std::shared_ptr< HepMC3::GenEvent > getGenEvent() const
Get the underlying HepMC3 GenEvent.
int getNumParticles() const
Get the number of particles in the event.
std::shared_ptr< HepMC3::GenEvent > event_
The underlying HepMC3 GenEvent.
Definition HepMCEvent.h:91
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...