LDMX Software
simcore::hepmc::HepMCEvent Class Reference

Wrapper for HepMC3::GenEvent with convenience methods. More...

#include <HepMCEvent.h>

Public Member Functions

 HepMCEvent (std::shared_ptr< HepMC3::GenEvent > event)
 Class constructor.
 
virtual ~HepMCEvent ()=default
 Class destructor.
 
int getNumParticles () const
 Get the number of particles in the event.
 
double getEventWeight () const
 Get the event weight.
 
const double * getVertex () const
 Get the vertex location (in mm, as expected by Geant4).
 
double getVertexTime () const
 Get the vertex time (in ns, as expected by Geant4).
 
const std::vector< std::unique_ptr< HepMCParticle > > & getParticles () const
 Get the list of final state particles in the event.
 
std::shared_ptr< HepMC3::GenEvent > getGenEvent () const
 Get the underlying HepMC3 GenEvent.
 

Private Member Functions

void extractParticles () const
 Extract final state particles from the HepMC event.
 

Private Attributes

std::shared_ptr< HepMC3::GenEvent > event_
 The underlying HepMC3 GenEvent.
 
double vtx_ [3]
 Vertex location (in mm)
 
double vtxt_ {0.}
 Vertex time (in ns)
 
std::vector< std::unique_ptr< HepMCParticle > > particles_
 The list of final state particles to be tracked.
 
bool particles_extracted_ {false}
 Flag to indicate if particles have been extracted.
 

Detailed Description

Wrapper for HepMC3::GenEvent with convenience methods.

Note
Detailed information on the HepMC3 format is provided here: HepMC3: A modernized Monte Carlo event structure library.

Definition at line 37 of file HepMCEvent.h.

Constructor & Destructor Documentation

◆ HepMCEvent()

simcore::hepmc::HepMCEvent::HepMCEvent ( std::shared_ptr< HepMC3::GenEvent > event)

Class constructor.

Parameters
eventShared pointer to HepMC3 GenEvent

Definition at line 6 of file HepMCEvent.cxx.

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}
double vtxt_
Vertex time (in ns)
Definition HepMCEvent.h:101
double vtx_[3]
Vertex location (in mm)
Definition HepMCEvent.h:96
std::shared_ptr< HepMC3::GenEvent > event_
The underlying HepMC3 GenEvent.
Definition HepMCEvent.h:91

References event_, vtx_, and vtxt_.

Member Function Documentation

◆ extractParticles()

void simcore::hepmc::HepMCEvent::extractParticles ( ) const
private

Extract final state particles from the HepMC event.

Definition at line 70 of file HepMCEvent.cxx.

70 {
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}
bool particles_extracted_
Flag to indicate if particles have been extracted.
Definition HepMCEvent.h:111
std::vector< std::unique_ptr< HepMCParticle > > particles_
The list of final state particles to be tracked.
Definition HepMCEvent.h:106

References event_, particles_, and particles_extracted_.

Referenced by getNumParticles(), and getParticles().

◆ getEventWeight()

double simcore::hepmc::HepMCEvent::getEventWeight ( ) const

Get the event weight.

Returns
The event weight.

Definition at line 48 of file HepMCEvent.cxx.

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

References event_.

◆ getGenEvent()

std::shared_ptr< HepMC3::GenEvent > simcore::hepmc::HepMCEvent::getGenEvent ( ) const

Get the underlying HepMC3 GenEvent.

Returns
Shared pointer to the HepMC3 GenEvent

Definition at line 66 of file HepMCEvent.cxx.

66 {
67 return event_;
68}

References event_.

◆ getNumParticles()

int simcore::hepmc::HepMCEvent::getNumParticles ( ) const

Get the number of particles in the event.

Returns
The number of particles in event.

Definition at line 41 of file HepMCEvent.cxx.

41 {
44 }
45 return particles_.size();
46}
void extractParticles() const
Extract final state particles from the HepMC event.

References extractParticles(), particles_, and particles_extracted_.

◆ getParticles()

const std::vector< std::unique_ptr< HepMCParticle > > & simcore::hepmc::HepMCEvent::getParticles ( ) const

Get the list of final state particles in the event.

These are particles that should be tracked by Geant4.

Returns
The list of final state particles in the event.

Definition at line 58 of file HepMCEvent.cxx.

59 {
62 }
63 return particles_;
64}

References extractParticles(), particles_, and particles_extracted_.

◆ getVertex()

const double * simcore::hepmc::HepMCEvent::getVertex ( ) const

Get the vertex location (in mm, as expected by Geant4).

Returns
Array double[3] with x, y, z ordering

Definition at line 54 of file HepMCEvent.cxx.

54{ return vtx_; }

References vtx_.

◆ getVertexTime()

double simcore::hepmc::HepMCEvent::getVertexTime ( ) const

Get the vertex time (in ns, as expected by Geant4).

Returns
time of primary vertex.

Definition at line 56 of file HepMCEvent.cxx.

56{ return vtxt_; }

References vtxt_.

Member Data Documentation

◆ event_

std::shared_ptr<HepMC3::GenEvent> simcore::hepmc::HepMCEvent::event_
private

The underlying HepMC3 GenEvent.

Definition at line 91 of file HepMCEvent.h.

Referenced by extractParticles(), getEventWeight(), getGenEvent(), and HepMCEvent().

◆ particles_

std::vector<std::unique_ptr<HepMCParticle> > simcore::hepmc::HepMCEvent::particles_
mutableprivate

The list of final state particles to be tracked.

Definition at line 106 of file HepMCEvent.h.

Referenced by extractParticles(), getNumParticles(), and getParticles().

◆ particles_extracted_

bool simcore::hepmc::HepMCEvent::particles_extracted_ {false}
mutableprivate

Flag to indicate if particles have been extracted.

Definition at line 111 of file HepMCEvent.h.

111{false};

Referenced by extractParticles(), getNumParticles(), and getParticles().

◆ vtx_

double simcore::hepmc::HepMCEvent::vtx_[3]
mutableprivate

Vertex location (in mm)

Definition at line 96 of file HepMCEvent.h.

Referenced by getVertex(), and HepMCEvent().

◆ vtxt_

double simcore::hepmc::HepMCEvent::vtxt_ {0.}
mutableprivate

Vertex time (in ns)

Definition at line 101 of file HepMCEvent.h.

101{0.};

Referenced by getVertexTime(), and HepMCEvent().


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