EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ActsExamples::HepMC3Event Struct Reference

Helper struct to convert HepMC3 event to the internal format. More...

#include <acts/blob/sPHENIX/Examples/Io/HepMC3/include/ActsExamples/Plugins/HepMC3/HepMC3Event.hpp>

Public Member Functions

void momentumUnit (std::shared_ptr< HepMC3::GenEvent > event, const double momentumUnit)
 Sets new units for momentums.
 
void lengthUnit (std::shared_ptr< HepMC3::GenEvent > event, const double lengthUnit)
 Sets new units for lengths.
 
void shiftPositionBy (std::shared_ptr< HepMC3::GenEvent > event, const Acts::Vector3D &deltaPos, const double deltaTime)
 Shifts the positioning of an event in space and time.
 
void shiftPositionTo (std::shared_ptr< HepMC3::GenEvent > event, const Acts::Vector3D &pos, const double time)
 Shifts the positioning of an event to a paint in space and time.
 
void shiftPositionTo (std::shared_ptr< HepMC3::GenEvent > event, const Acts::Vector3D &pos)
 Shifts the positioning of an event to a paint in space.
 
void shiftPositionTo (std::shared_ptr< HepMC3::GenEvent > event, const double time)
 Shifts the positioning of an event to a paint in time.
 
void addParticle (std::shared_ptr< HepMC3::GenEvent > event, std::shared_ptr< SimParticle > particle)
 Adds a new particle.
 
void addVertex (std::shared_ptr< HepMC3::GenEvent > event, const std::shared_ptr< SimVertex > vertex)
 Adds a new vertex.
 
void removeParticle (std::shared_ptr< HepMC3::GenEvent > event, const std::shared_ptr< SimParticle > &particle)
 Removes a particle from the record.
 
void removeVertex (std::shared_ptr< HepMC3::GenEvent > event, const std::shared_ptr< SimVertex > &vertex)
 Removes a vertex from the record.
 
double momentumUnit (const std::shared_ptr< HepMC3::GenEvent > event)
 Getter of the unit of momentum used.
 
double lengthUnit (const std::shared_ptr< HepMC3::GenEvent > event)
 Getter of the unit of length used.
 
Acts::Vector3D eventPos (const std::shared_ptr< HepMC3::GenEvent > event)
 Getter of the position of the event.
 
double eventTime (const std::shared_ptr< HepMC3::GenEvent > event)
 Getter of the time of the event.
 
std::vector< std::unique_ptr
< SimParticle > > 
particles (const std::shared_ptr< HepMC3::GenEvent > event)
 Get list of particles.
 
std::vector< std::unique_ptr
< SimVertex > > 
vertices (const std::shared_ptr< HepMC3::GenEvent > event)
 Get list of vertices.
 
std::vector< std::unique_ptr
< SimParticle > > 
beams (const std::shared_ptr< HepMC3::GenEvent > event)
 Get beam particles.
 
std::vector< std::unique_ptr
< SimParticle > > 
finalState (const std::shared_ptr< HepMC3::GenEvent > event)
 Get final state particles.
 

Private Member Functions

HepMC3::GenParticlePtr actsParticleToGen (std::shared_ptr< SimParticle > actsParticle)
 Converts an SimParticle into HepMC3::GenParticle.
 
HepMC3::GenVertexPtr createGenVertex (const std::shared_ptr< SimVertex > &actsVertex)
 Converts an Acts vertex to a HepMC3::GenVertexPtr.
 
bool compareVertices (const std::shared_ptr< SimVertex > &actsVertex, const HepMC3::GenVertexPtr &genVertex)
 Compares an Acts vertex with a HepMC3::GenVertex.
 

Detailed Description

Helper struct to convert HepMC3 event to the internal format.

Definition at line 24 of file HepMC3Event.hpp.

View newest version in sPHENIX GitHub at line 24 of file HepMC3Event.hpp

Member Function Documentation

HepMC3::GenParticlePtr ActsExamples::HepMC3Event::actsParticleToGen ( std::shared_ptr< SimParticle actsParticle)
private

Converts an SimParticle into HepMC3::GenParticle.

Note
The conversion ignores HepMC status codes
Parameters
actsParticleActs particle that will be converted
Returns
converted particle

Adder

Definition at line 91 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 91 of file HepMC3Event.cpp

void ActsExamples::HepMC3Event::addParticle ( std::shared_ptr< HepMC3::GenEvent >  event,
std::shared_ptr< SimParticle particle 
)

Adds a new particle.

Adder

Parameters
eventevent in HepMC data type
particlenew particle that will be added

Definition at line 103 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 103 of file HepMC3Event.cpp

void ActsExamples::HepMC3Event::addVertex ( std::shared_ptr< HepMC3::GenEvent >  event,
const std::shared_ptr< SimVertex vertex 
)

Adds a new vertex.

Parameters
eventevent in HepMC data type
vertexnew vertex that will be added
Note
The statuses are not represented in Acts and therefore set to 0

Definition at line 134 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 134 of file HepMC3Event.cpp

std::vector< std::unique_ptr< ActsExamples::SimParticle > > ActsExamples::HepMC3Event::beams ( const std::shared_ptr< HepMC3::GenEvent >  event)

Get beam particles.

Parameters
eventevent in HepMC data type
Returns
List of beam particles

Definition at line 263 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 263 of file HepMC3Event.cpp

References ActsExamples::HepMC3Particle::particle().

Referenced by main().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

bool ActsExamples::HepMC3Event::compareVertices ( const std::shared_ptr< SimVertex > &  actsVertex,
const HepMC3::GenVertexPtr &  genVertex 
)
private

Compares an Acts vertex with a HepMC3::GenVertex.

Note
An Acts vertex does not store a barcode. Therefore the content of both vertices is compared. The position, time and number of incoming and outgoing particles will be compared. Since a second vertex could exist in the record with identical informations (although unlikely), this comparison could lead to false positive results. On the other hand, a numerical deviation of the parameters could lead to a false negative.
Parameters
actsVertexActs vertex
genVertexHepMC3::GenVertex
Returns
boolean result if both vertices are identical

Definition at line 160 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 160 of file HepMC3Event.cpp

HepMC3::GenVertexPtr ActsExamples::HepMC3Event::createGenVertex ( const std::shared_ptr< SimVertex > &  actsVertex)
private

Converts an Acts vertex to a HepMC3::GenVertexPtr.

Note
The conversion ignores HepMC status codes
Parameters
actsVertexActs vertex that will be converted
Returns
Converted Acts vertex to HepMC3::GenVertexPtr

Definition at line 110 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 110 of file HepMC3Event.cpp

References particle.

Acts::Vector3D ActsExamples::HepMC3Event::eventPos ( const std::shared_ptr< HepMC3::GenEvent >  event)

Getter of the position of the event.

Parameters
eventevent in HepMC data type
Returns
vector to the location of the event

Definition at line 214 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 214 of file HepMC3Event.cpp

Referenced by main().

+ Here is the caller graph for this function:

double ActsExamples::HepMC3Event::eventTime ( const std::shared_ptr< HepMC3::GenEvent >  event)

Getter of the time of the event.

Parameters
eventevent in HepMC data type
Returns
time of the event

Definition at line 224 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 224 of file HepMC3Event.cpp

Referenced by main().

+ Here is the caller graph for this function:

std::vector< std::unique_ptr< ActsExamples::SimParticle > > ActsExamples::HepMC3Event::finalState ( const std::shared_ptr< HepMC3::GenEvent >  event)

Get final state particles.

Parameters
eventevent in HepMC data type
Returns
List of final state particles

Definition at line 278 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 278 of file HepMC3Event.cpp

References ActsExamples::HepMC3Particle::particle(), and particle.

Referenced by main().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void ActsExamples::HepMC3Event::lengthUnit ( std::shared_ptr< HepMC3::GenEvent >  event,
const double  lengthUnit 
)

Sets new units for lengths.

Note
The allowed units are mm and cm
Parameters
eventevent in HepMC data type
lengthUnitnew unit of length

Definition at line 36 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 36 of file HepMC3Event.cpp

References Acts::UnitConstants::cm, and Acts::UnitConstants::mm.

Referenced by main().

+ Here is the caller graph for this function:

double ActsExamples::HepMC3Event::lengthUnit ( const std::shared_ptr< HepMC3::GenEvent >  event)

Getter of the unit of length used.

Parameters
eventevent in HepMC data type
Returns
unit of length

Definition at line 206 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 206 of file HepMC3Event.cpp

References Acts::UnitConstants::cm, and Acts::UnitConstants::mm.

void ActsExamples::HepMC3Event::momentumUnit ( std::shared_ptr< HepMC3::GenEvent >  event,
const double  momentumUnit 
)

Sets new units for momentums.

Setter

Note
The allowed units are MeV and Gev
Parameters
eventevent in HepMC data type
momentumUnitnew unit of momentum

Setter

Definition at line 18 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 18 of file HepMC3Event.cpp

References Acts::UnitConstants::GeV, Acts::UnitConstants::MeV, and charm_jet_coverage::mom.

Referenced by main().

+ Here is the caller graph for this function:

double ActsExamples::HepMC3Event::momentumUnit ( const std::shared_ptr< HepMC3::GenEvent >  event)

Getter of the unit of momentum used.

Getter

Parameters
eventevent in HepMC data type
Returns
unit of momentum

Getter

Definition at line 198 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 198 of file HepMC3Event.cpp

References Acts::UnitConstants::GeV, and Acts::UnitConstants::MeV.

std::vector< std::unique_ptr< ActsExamples::SimParticle > > ActsExamples::HepMC3Event::particles ( const std::shared_ptr< HepMC3::GenEvent >  event)

Get list of particles.

Parameters
eventevent in HepMC data type
Returns
List of particles

Definition at line 231 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 231 of file HepMC3Event.cpp

References ActsExamples::HepMC3Particle::particle().

Referenced by main().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void ActsExamples::HepMC3Event::removeParticle ( std::shared_ptr< HepMC3::GenEvent >  event,
const std::shared_ptr< SimParticle > &  particle 
)

Removes a particle from the record.

Remover

Parameters
eventevent in HepMC data type
particleparticle that will be removed

Remover

Definition at line 145 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 145 of file HepMC3Event.cpp

void ActsExamples::HepMC3Event::removeVertex ( std::shared_ptr< HepMC3::GenEvent >  event,
const std::shared_ptr< SimVertex > &  vertex 
)

Removes a vertex from the record.

Note
The identification of the vertex is potentially unstable (c.f. HepMC3Event::compareVertices())
Parameters
eventevent in HepMC data type
vertexvertex that will be removed

Definition at line 181 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 181 of file HepMC3Event.cpp

void ActsExamples::HepMC3Event::shiftPositionBy ( std::shared_ptr< HepMC3::GenEvent >  event,
const Acts::Vector3D deltaPos,
const double  deltaTime 
)

Shifts the positioning of an event in space and time.

Parameters
eventevent in HepMC data type
deltaPosrelative spatial shift that will be applied
deltaTimerelative time shift that will be applied

Definition at line 55 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 55 of file HepMC3Event.cpp

void ActsExamples::HepMC3Event::shiftPositionTo ( std::shared_ptr< HepMC3::GenEvent >  event,
const Acts::Vector3D pos,
const double  time 
)

Shifts the positioning of an event to a paint in space and time.

Parameters
eventevent in HepMC data type
posnew position of the event
timenew time of the event

Definition at line 64 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 64 of file HepMC3Event.cpp

References pos().

+ Here is the call graph for this function:

void ActsExamples::HepMC3Event::shiftPositionTo ( std::shared_ptr< HepMC3::GenEvent >  event,
const Acts::Vector3D pos 
)

Shifts the positioning of an event to a paint in space.

Parameters
eventevent in HepMC data type
posnew position of the event

Definition at line 72 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 72 of file HepMC3Event.cpp

References pos().

+ Here is the call graph for this function:

void ActsExamples::HepMC3Event::shiftPositionTo ( std::shared_ptr< HepMC3::GenEvent >  event,
const double  time 
)

Shifts the positioning of an event to a paint in time.

Parameters
eventevent in HepMC data type
timenew time of the event

Definition at line 79 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 79 of file HepMC3Event.cpp

References Acts::Test::time.

std::vector< std::unique_ptr< ActsExamples::SimVertex > > ActsExamples::HepMC3Event::vertices ( const std::shared_ptr< HepMC3::GenEvent >  event)

Get list of vertices.

Parameters
eventevent in HepMC data type
Returns
List of vertices

Definition at line 247 of file HepMC3Event.cpp.

View newest version in sPHENIX GitHub at line 247 of file HepMC3Event.cpp

References ActsExamples::HepMC3Vertex::processVertex().

Referenced by main().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:


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