Project
Loading...
Searching...
No Matches
o2::eventgen::AODToHepMC Struct Reference

#include <AODToHepMC.h>

Public Types

The containers we subscribe to
using Headers = o2::aod::McCollisions
 
using Header = o2::aod::McCollision
 
using Tracks = o2::aod::McParticles
 
using Track = typename Tracks::iterator
 
using XSections = o2::aod::HepMCXSections
 
using PdfInfos = o2::aod::HepMCPdfInfos
 
using HeavyIons = o2::aod::HepMCHeavyIons
 
using XSection = typename XSections::iterator
 
using PdfInfo = typename PdfInfos::iterator
 
using HeavyIon = typename HeavyIons::iterator
 
Types from HepMC3
using FourVector = HepMC3::FourVector
 
using ParticlePtr = HepMC3::GenParticlePtr
 
using VertexPtr = HepMC3::GenVertexPtr
 
using Event = HepMC3::GenEvent
 
using HeavyIonPtr = HepMC3::GenHeavyIonPtr
 
using CrossSectionPtr = HepMC3::GenCrossSectionPtr
 
using PdfInfoPtr = HepMC3::GenPdfInfoPtr
 
using ParticleMap = std::unordered_map< long, ParticlePtr >
 
using ParticleVector = std::vector< ParticlePtr >
 
using ParticleList = std::list< ParticlePtr >
 
using VertexVector = std::vector< VertexPtr >
 
using VertexList = std::list< VertexPtr >
 
using WriterAscii = HepMC3::WriterAscii
 
using WriterAsciiPtr = std::shared_ptr< WriterAscii >
 

Public Member Functions

Interface member functions
virtual void init ()
 
virtual void startEvent ()
 
virtual void process (Header const &collision, Tracks const &tracks)
 
virtual void process (Header const &collision, XSections const &xsections, PdfInfos const &pdfs, HeavyIons const &heavyions)
 
virtual void endEvent ()
 
virtual bool postRun ()
 

Public Attributes

struct { 
 
   std::string   dump {""} 
 
   bool   onlyGen {false} 
 
   bool   useTree {false} 
 
   int   precision {8} 
 
   bool   recenter {false} 
 
configs 
 
HepMC3 objects
Event mEvent
 
CrossSectionPtr mCrossSec = nullptr
 
HeavyIonPtr mIon = nullptr
 
PdfInfoPtr mPdf = nullptr
 
Containers etc.
ParticleMap mParticles
 
VertexList mVertices
 Cache of particles.
 
ParticleVector mBeams
 Cache of vertices.
 
ParticleList mOrphans
 Cache of beam particles.
 
Options and such

Cache of particles w/o mothers

WriterAsciiPtr mWriter = nullptr
 
int mEventNo = 0
 
int mLastBC = -1
 
bool mOnlyGen = false
 
bool mUseTree = true
 
std::ofstream * mOutput = nullptr
 
int mPrecision = 16
 
bool mRecenter = false
 

Protected Member Functions

void enableDump (const std::string &dump)
 
Actual worker member functions
virtual void makeEvent (Header const &collision, Tracks const &tracks)
 
virtual void makeHeader (Header const &header)
 
virtual void makeXSection (XSections const &xsection)
 
virtual void makePdfInfo (PdfInfos const &pdf)
 
virtual void makeHeavyIon (HeavyIons const &heavyion, Header const &header)
 
virtual void makeBeams (Header const &header, const VertexPtr ip)
 
virtual void makeParticles (Tracks const &tracks)
 
virtual ParticlePtr getParticle (Track const &ref) const
 
virtual bool isIgnored (Track const &track) const
 
virtual ParticlePtr makeParticleRecursive (Track const &track, Tracks const &tracks, bool force=false)
 
virtual ParticlePtr makeParticle (const Track &track, Int_t &motherStatus, bool force) const
 
virtual VertexPtr makeVertex (const Track &track)
 
virtual void fleshOut (Tracks const &tracks)
 
virtual void fleshOutParticle (Track const &track, Tracks const &tracks)
 
virtual void recenter (Header const &collision)
 

Detailed Description

Convert AOD tables of MC information into a HepMC event structure.

The event structure is kept in memory.

The conventions used here are

  • A track is a kinematics particle stored in the o2::aod::McParticles table (aliased as AODToHepMC::Tracks), and correspond to a track used during simulation transport. These are of type o2::aod::McParticle aliased to AODToHepMC::Track.
    • A generated track is a particle originally made by the event generator. These can be recognised by Track::producedByGenerator()
  • A particle is a particle stored in the HepMC event structure of the class HepMC3::GenParticle.
  • A vertex is where a particle is either produced or disappears. A vertex with incoming particles will have a number of outgoing particles. Thus a vertex is the production vertex for out-going particles of that vertex, and the end vertex of incoming particles of that vertex. Vertexes are of the type HepMC3::GenVertex.
    • The relationship between mother and daughter tracks in AODToHepMC::Tracks is encoded by indexes.
      • At most two mothers can be stored per track. If a track has no mothers, then both indexes are -1. If a track has a single mother track, then both indexes are the same.
      • A track can have any number of daughters. The first daughter track is identified by the first index and the last daughter by the second stored index. All tracks indexes with in the range from the first to the second index (both inclusive) are daughter tracks. If a track has no daughters, then both indexes are -1. If a track has one daughter, then both indexes are equal. The number of daughters can be obtained by
        last - first + 1
        
        if both last and first are zero or greater.
  • An event header is the information stored in a row of o2::aod::McCollisions table (aliased to AODToHepMC::Headers). Each event has one such header (aliased as AODToHepMC::Header).

    The header stores information on the event such as event number, weight, primary interaction point (IP), and so on.

    In addition, auxiliary information (impact parameter, $N_{\mathrm{part}}$, and so on) may be stored separate tables.

    • The table o2::aod::HepMCXSections (aliased to AODToHepMC::XSections) stores the event cross section and weight.
    • The table o2::aod::HepMCPdfInfos (aliased to AODToHepMC::PdfInfos) stores the parton distribution function parameters used in the event.
    • The table o2::aod::HepMCHeavyIons (aliased to AODToHepMC::PdfHeavyIons) stores the heavy-ion auxiliary information, such as the event impact parameter $b$, $N_{\mathrm{part}}$, $N_{\mathrm{coll}}$, and so on.
  • A HepMC event has a simple header which contains the event number, number of particles, and number of vertexes in the event.

    Other event information is stored in specific structures.

The conversion is done as follows:

  • For all MC tracks, create the correspond particle (AODToHepMC::makeParticleRecursive)
    • Check if we have already created the corresponding particle (AODToHepMC::getParticle). If so, then we go on to the next track.
    • If we are asked to only keep generated tracks (AODToHepMC::mOnlyGen set by option --hepmc-only-generated), i.e., tracks that were posted by the event generator, and this track is not such a track (AODToHepMC::isIgnored), nor forced (argument force) then return nothing. Note that mothers of particles that we keep are always made (force is set to true).
    • If we do not have a particle yet for the track (look-up in AODToHepMC::mParticles) and it is not excluded, then create a particle that corresponds to the track, and store the mapping from track to particle (in AODToHepMC::mParticles).
    • For all mother tracks of the current track, create the corresponding particle, following this algorithm (recursion, AODToHepMC::makeParticleRecursive).
    • If a mother particle has an end vertex, set that vertex as the production vertex of the current particle.
    • If no mother particle has an end vertex, and this particle is not a beam particle and it does have mothers, then create a vertex (AODToHepMC::makeVertex) corresponding to the track production position, and add this particle as an outgoing particle of that vertex.
      • In this case, if some mother does not have an outgoing vertex, add that mother as an incoming particle to the created vertex.
    • If the particle is a beam particle (status=4) then store this particle as a beam particle.
    • If not a beam particle, and the particle has no mothers, mark this particle as an orphan.
  • Once we have made particle for all tracks, we flesh-out all particles. For all tracks (AODToHepMC::fleshOutParticle)
    • If this track is ignored (AODToHepMC::isIgnored), or we have no particle (AODToHepMC::getParticle) corresponding to this track, go on to the next track.
    • Get the end vertex of the particle. If any. Set the candidate end vertex to this vertex, whether it exists or not.
    • Then for each daughter track of the current track, check if it is consistent with the end vertex of the particle.
      • Check if the dauther is ignored (AODToHepMC::isIgnored). If so, move on to the next daughter track.
      • Get the particle corresponding to the daughter track (AODToHepMC::getParticle.). If no particle is found, move on to the next daughter track.
      • Check that the daughter particle has an end vertex. If it doesn't, mark it as a head-less particle.
      • If the production vertex of the daughter doesn't match the end vertex of the current particle, or the current candidate end vertex, then issue a warning, and move on to the next daughter.
      • Update the candidate end vertex to the daughter end vertex.
    • After processing all daughter particles, and if we have no end vertex for the current particle, then
      • if we found no candiate end vertex, and the particle is either a beam (status=4) or decayed (status=2) particle, issue a warning.
        • if we do have a candidate end vertex, set that as the end vertex of the current particle.
    • If, after this, the current particle does have an end vertex, loop over all daughters previsouly as head-less and set their production vertex to this end vertex.
  • At this point, we should have the particles in a proper HepMC event structure.
  • During simulation transport, the interaction point vertex (IP) may not be at (0,0,0,0). Since some consumers of the the HepMC event structure may expect the IP to be at (0,0,0,0), we can recenter (AODMCToHepMC::recenter) all vertexes of the event. This is governed by the member AODMCToHepMC::mRecenter set by the option --hepmc-recenter
  • We can then fill in remaing information into the HepMC event header.
  • Once all the above is done, we have a complete HepMC event (AODToHepMC::mEvent). This event structure is kept in memory.
    • Optionally (option --hepmc-dumb filename) we may write the events to a file (AODToHepMC::mOutput). The event is still kept in memory.
    • Clients of this class (e.g., o2::pwgmm::RivetWrapper) may access the event structure for further processing.

The utility o2-aod-mc-to-hepmc will read in AODs and write out a HepMC event file (plain ASCII).

Definition at line 253 of file AODToHepMC.h.

Member Typedef Documentation

◆ CrossSectionPtr

using o2::eventgen::AODToHepMC::CrossSectionPtr = HepMC3::GenCrossSectionPtr

Alias (smart-)pointer to HepMC cross section object

Definition at line 317 of file AODToHepMC.h.

◆ Event

using o2::eventgen::AODToHepMC::Event = HepMC3::GenEvent

Alias HepMC eventt structure

Definition at line 313 of file AODToHepMC.h.

◆ FourVector

using o2::eventgen::AODToHepMC::FourVector = HepMC3::FourVector

Alias HepMC four-vector

Definition at line 307 of file AODToHepMC.h.

◆ Header

Alias of MC collisions table row

Definition at line 282 of file AODToHepMC.h.

◆ Headers

Alias of MC collisions table

Definition at line 280 of file AODToHepMC.h.

◆ HeavyIon

using o2::eventgen::AODToHepMC::HeavyIon = typename HeavyIons::iterator

Alias row of auxiliary MC table of heavy-ion parameters

Definition at line 300 of file AODToHepMC.h.

◆ HeavyIonPtr

using o2::eventgen::AODToHepMC::HeavyIonPtr = HepMC3::GenHeavyIonPtr

Alias (smart-)pointer to HepMC heavy-ion object

Definition at line 315 of file AODToHepMC.h.

◆ HeavyIons

using o2::eventgen::AODToHepMC::HeavyIons = o2::aod::HepMCHeavyIons

Alias auxiliary MC table of heavy-ion parameters

Definition at line 293 of file AODToHepMC.h.

◆ ParticleList

A container of pointers to particles

Definition at line 326 of file AODToHepMC.h.

◆ ParticleMap

Type used to map tracks to HepMC particles

Definition at line 322 of file AODToHepMC.h.

◆ ParticlePtr

using o2::eventgen::AODToHepMC::ParticlePtr = HepMC3::GenParticlePtr

Alias (smart-)pointer to HepMC particle

Definition at line 309 of file AODToHepMC.h.

◆ ParticleVector

A container of pointers to particles

Definition at line 324 of file AODToHepMC.h.

◆ PdfInfo

using o2::eventgen::AODToHepMC::PdfInfo = typename PdfInfos::iterator

Alias row of auxiliary MC table of parton distribution function parameters

Definition at line 298 of file AODToHepMC.h.

◆ PdfInfoPtr

using o2::eventgen::AODToHepMC::PdfInfoPtr = HepMC3::GenPdfInfoPtr

Alias (smart-)pointer to HepMC parton distribution function object

Definition at line 320 of file AODToHepMC.h.

◆ PdfInfos

using o2::eventgen::AODToHepMC::PdfInfos = o2::aod::HepMCPdfInfos

Alias auxiliary MC table of parton distribution function parameters

Definition at line 291 of file AODToHepMC.h.

◆ Track

using o2::eventgen::AODToHepMC::Track = typename Tracks::iterator

Alias MC particles (tracks) table row

Definition at line 286 of file AODToHepMC.h.

◆ Tracks

Alias MC particles (tracks) table

Definition at line 284 of file AODToHepMC.h.

◆ VertexList

A container of pointers to vertexes

Definition at line 330 of file AODToHepMC.h.

◆ VertexPtr

using o2::eventgen::AODToHepMC::VertexPtr = HepMC3::GenVertexPtr

Alias (smart-)pointer to HepMC vertex

Definition at line 311 of file AODToHepMC.h.

◆ VertexVector

A container of pointers to vertexes

Definition at line 328 of file AODToHepMC.h.

◆ WriterAscii

using o2::eventgen::AODToHepMC::WriterAscii = HepMC3::WriterAscii

Alias of HepMC writer class

Definition at line 332 of file AODToHepMC.h.

◆ WriterAsciiPtr

The of pointer to HepMC writer class

Definition at line 334 of file AODToHepMC.h.

◆ XSection

using o2::eventgen::AODToHepMC::XSection = typename XSections::iterator

Alias row of auxiliary MC table of cross-sections

Definition at line 295 of file AODToHepMC.h.

◆ XSections

using o2::eventgen::AODToHepMC::XSections = o2::aod::HepMCXSections

Alias auxiliary MC table of cross-sections

Definition at line 288 of file AODToHepMC.h.

Member Function Documentation

◆ enableDump()

void o2::eventgen::AODToHepMC::enableDump ( const std::string &  dump)
protected

Open dump output, or close if an empty string was given.

Parameters
dump

Definition at line 561 of file AODToHepMC.cxx.

◆ endEvent()

void o2::eventgen::AODToHepMC::endEvent ( )
virtual

Call after process an. Thisf finalises the event and optionally outputs to dump.

Definition at line 71 of file AODToHepMC.cxx.

◆ fleshOut()

void o2::eventgen::AODToHepMC::fleshOut ( Tracks const &  tracks)
protectedvirtual

Visit all tracks, but this time look for daughters and attach mothers as incoming particles if not already done.

Definition at line 436 of file AODToHepMC.cxx.

◆ fleshOutParticle()

void o2::eventgen::AODToHepMC::fleshOutParticle ( Track const &  track,
Tracks const &  tracks 
)
protectedvirtual

Flesh out a single particle

Parameters
trackThe track for which we should flesh out the corresponding particle.

Definition at line 443 of file AODToHepMC.cxx.

◆ getParticle()

AODToHepMC::ParticlePtr o2::eventgen::AODToHepMC::getParticle ( Track const &  ref) const
protectedvirtual

Get particle corresponding to track no from particle cache

Parameters
refTrack reference
Returns
Pointer to HepMC3::GenParticle or null

Definition at line 292 of file AODToHepMC.cxx.

◆ init()

void o2::eventgen::AODToHepMC::init ( )
virtual

Initialize the converter. Sets internal parameters based on the configurables.

Definition at line 21 of file AODToHepMC.cxx.

◆ isIgnored()

bool o2::eventgen::AODToHepMC::isIgnored ( Track const &  track) const
protectedvirtual

Check if we are ignoring this track

Parameters
trackTrack to check

Definition at line 302 of file AODToHepMC.cxx.

◆ makeBeams()

virtual void o2::eventgen::AODToHepMC::makeBeams ( Header const &  header,
const VertexPtr  ip 
)
inlineprotectedvirtual

This is supposed to make the beam particles from the information available. However, it seems like we really don't have enough information, so for now we will do nothing here. Perhaps the user will be forced to provide that information - either via the analyser configurables or somehow from somewhere else.

Definition at line 475 of file AODToHepMC.h.

◆ makeEvent()

void o2::eventgen::AODToHepMC::makeEvent ( Header const &  collision,
Tracks const &  tracks 
)
protectedvirtual

Generate the final event, including fleshing out the vertexes, and so on

Parameters
collisionHeader information
tracksParticle tracks

Definition at line 82 of file AODToHepMC.cxx.

◆ makeHeader()

void o2::eventgen::AODToHepMC::makeHeader ( Header const &  header)
protectedvirtual

Set the various fields in the header of the HepMC3::GenEvent object

Parameters
headerHeader object

Definition at line 130 of file AODToHepMC.cxx.

◆ makeHeavyIon()

void o2::eventgen::AODToHepMC::makeHeavyIon ( HeavyIons const &  heavyion,
Header const &  header 
)
protectedvirtual

Make heavy-ion collision information. If no header given, then fill in other reasonable values

Definition at line 204 of file AODToHepMC.cxx.

◆ makeParticle()

AODToHepMC::ParticlePtr o2::eventgen::AODToHepMC::makeParticle ( const Track track,
Int_t &  motherStatus,
bool  force 
) const
protectedvirtual

Generate a HepMC particle from a track. Note that the job here is simply to make the object. The more complicated job of adding the track to the tree is done above in makeParticleRecursive

Parameters
trackMC track
mstMother status code - updated on return
forceForce generation of particle even if onlyGen
Returns
Shared pointer to new HepMC particle

Definition at line 398 of file AODToHepMC.cxx.

◆ makeParticleRecursive()

AODToHepMC::ParticlePtr o2::eventgen::AODToHepMC::makeParticleRecursive ( Track const &  track,
Tracks const &  tracks,
bool  force = false 
)
protectedvirtual

Truely make a particle, and its mother particles if any. We add vertexes as needed here too.

Note that this traverses the tree from the bottom up. That is, we check if a particle has any mothers, and if it does, we create those.

However, this can be a bit problematic, since the Kinematic tree (and thus the McParticles table) only allows for at most 2 mothers. In the case a vertex has 3 or more incoming particles, then some of the intermediate particles will be lost.

We remedy that problem by traversing the tree again, but this time from the bottom up - that is, we look for daughters of all particles, and if they are not registered with the out-going vertex, then we reattch the parent to the incoming vertex of the daughters. For this to work, we need to map from track index to HepMC3::GenParticle.

Parameters
trackCurrent track
tracksMC tracks
forceIf true, do make the particle even if it would otherwise be skipped.
Returns
Pointer to created particle.

Definition at line 313 of file AODToHepMC.cxx.

◆ makeParticles()

void o2::eventgen::AODToHepMC::makeParticles ( Tracks const &  tracks)
protectedvirtual

Make all particles. We loop through the MC particles, and for each of them create that particle and any possible mother particles (recursively). This allows us to traverse the data rather straight-forwardly.

Parameters
tracksThe MC tracks

Definition at line 282 of file AODToHepMC.cxx.

◆ makePdfInfo()

void o2::eventgen::AODToHepMC::makePdfInfo ( PdfInfos const &  pdf)
protectedvirtual

Make parton-distribition function information. If no entry in the table, then make dummy information

Definition at line 175 of file AODToHepMC.cxx.

◆ makeVertex()

AODToHepMC::VertexPtr o2::eventgen::AODToHepMC::makeVertex ( const Track track)
protectedvirtual

Generate vertex from production vertex of track

Parameters
trackMC track
Returns
Shared pointer to new HepMC vertex

Definition at line 591 of file AODToHepMC.cxx.

◆ makeXSection()

void o2::eventgen::AODToHepMC::makeXSection ( XSections const &  xsection)
protectedvirtual

Make cross-section information. If no entry in the table, then make dummy information

Definition at line 151 of file AODToHepMC.cxx.

◆ postRun()

virtual bool o2::eventgen::AODToHepMC::postRun ( )
inlinevirtual

End of run - closes output file if enabled. This is called via specialisation of o2::framework::OutputManager<AODToHepMC>.

Definition at line 426 of file AODToHepMC.h.

◆ process() [1/2]

void o2::eventgen::AODToHepMC::process ( Header const &  collision,
Tracks const &  tracks 
)
virtual

Process the collision header and tracks

Parameters
collisionHeader information
tracksParticle tracks

Definition at line 47 of file AODToHepMC.cxx.

◆ process() [2/2]

void o2::eventgen::AODToHepMC::process ( Header const &  collision,
XSections const &  xsections,
PdfInfos const &  pdfs,
HeavyIons const &  heavyions 
)
virtual

Process collision header and HepMC auxiliary information

Parameters
collisionHeader information
xsectionsCross-section table (possible no rows)
pdfsParton-distribution function table (possible no rows)
heavyionsHeavy ion collision table (possible no rows)

Definition at line 60 of file AODToHepMC.cxx.

◆ recenter()

virtual void o2::eventgen::AODToHepMC::recenter ( Header const &  collision)
protectedvirtual

Recenter event to (0,0,0,0). This will use the vertex information from the event header to translate all vertexes in the event.

Parameters
headerEvent header

◆ startEvent()

void o2::eventgen::AODToHepMC::startEvent ( )
virtual

Call before starting to process an event. This clears the current event and internal data structures.

Definition at line 37 of file AODToHepMC.cxx.

Member Data Documentation

◆ [struct]

struct { ... } o2::eventgen::AODToHepMC::configs

Group of configurables which will be added to an program that uses this class. Note that it is really the specialisation of framework::OptionManager<AODToHepMC> that propagates the options to the program.

◆ dump

std::string o2::eventgen::AODToHepMC::dump {""}

Option for dumping HepMC event structures to disk. Takes one argument - the name of the file to write to.

Definition at line 263 of file AODToHepMC.h.

◆ mBeams

ParticleVector o2::eventgen::AODToHepMC::mBeams

Cache of vertices.

List of beam particles

Definition at line 358 of file AODToHepMC.h.

◆ mCrossSec

CrossSectionPtr o2::eventgen::AODToHepMC::mCrossSec = nullptr

Pointer to cross section-ion information

Definition at line 343 of file AODToHepMC.h.

◆ mEvent

Event o2::eventgen::AODToHepMC::mEvent

The result of processing

Definition at line 341 of file AODToHepMC.h.

◆ mEventNo

int o2::eventgen::AODToHepMC::mEventNo = 0

Current sequential event number

Definition at line 369 of file AODToHepMC.h.

◆ mIon

HeavyIonPtr o2::eventgen::AODToHepMC::mIon = nullptr

Pointer to heavy-ion information

Definition at line 345 of file AODToHepMC.h.

◆ mLastBC

int o2::eventgen::AODToHepMC::mLastBC = -1

The last bunch crossing identifier

Definition at line 371 of file AODToHepMC.h.

◆ mOnlyGen

bool o2::eventgen::AODToHepMC::mOnlyGen = false

If true, only store particles from the generator

Definition at line 373 of file AODToHepMC.h.

◆ mOrphans

ParticleList o2::eventgen::AODToHepMC::mOrphans

Cache of beam particles.

Particles without a mother

Definition at line 360 of file AODToHepMC.h.

◆ mOutput

std::ofstream* o2::eventgen::AODToHepMC::mOutput = nullptr

Output stream if enabled

Definition at line 377 of file AODToHepMC.h.

◆ mParticles

ParticleMap o2::eventgen::AODToHepMC::mParticles

Maps tracks to particles

Definition at line 354 of file AODToHepMC.h.

◆ mPdf

PdfInfoPtr o2::eventgen::AODToHepMC::mPdf = nullptr

Pointer to parton distribution function information

Definition at line 347 of file AODToHepMC.h.

◆ mPrecision

int o2::eventgen::AODToHepMC::mPrecision = 16

Precision used on the output stream

Definition at line 379 of file AODToHepMC.h.

◆ mRecenter

bool o2::eventgen::AODToHepMC::mRecenter = false

If true, recenter IP to (0,0,0,0)

Definition at line 381 of file AODToHepMC.h.

◆ mUseTree

bool o2::eventgen::AODToHepMC::mUseTree = true

If true, use HepMC tree parser

Definition at line 375 of file AODToHepMC.h.

◆ mVertices

VertexList o2::eventgen::AODToHepMC::mVertices

Cache of particles.

List of vertexes made

Definition at line 356 of file AODToHepMC.h.

◆ mWriter

WriterAsciiPtr o2::eventgen::AODToHepMC::mWriter = nullptr

Output writer, if enabled

Definition at line 367 of file AODToHepMC.h.

◆ onlyGen

bool o2::eventgen::AODToHepMC::onlyGen {false}

Option for only storing particles from the event generator. Note, if a particle is stored down, then its mothers will also be stored.

Definition at line 267 of file AODToHepMC.h.

◆ precision

int o2::eventgen::AODToHepMC::precision {8}

Floating point precision used when writing to disk

Definition at line 271 of file AODToHepMC.h.

◆ recenter

void o2::eventgen::AODToHepMC::recenter {false}

Recenter event at IP=(0,0,0,0).

Definition at line 273 of file AODToHepMC.h.

◆ useTree

bool o2::eventgen::AODToHepMC::useTree {false}

Use HepMC's tree parsing for building event structure

Definition at line 269 of file AODToHepMC.h.


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