Output formats#

In CepGen, the event (particles, their parentage and kinematics) is handled through the cepgen::Event object described in the dedicated event format page.

To ease the user interaction with this object, a few output writers (defined here as “handlers”) are given as examples. All handlers are defined as modules derivating from the following abstract base class:

class EventExporter : public EventHandler#

Subclassed by DelphesHandler, EventHarvester, HepMC2Handler< T >, HepMC3Handler< T >, LHEFHepMCHandler, ProMCHandler, ROOTHistsHandler, TextEventHandler, TextVariablesHandler, YODAHistsHandler< T >, LHEFHandler, EventExporter

class EventExporter : public EventHandler

Output format handler for events export.

Author

Laurent Forthomme laurent.forthomme@cern.ch

Date

Sep 2016

Subclassed by DelphesHandler, EventHarvester, HepMC2Handler< T >, HepMC3Handler< T >, LHEFHepMCHandler, ProMCHandler, ROOTHistsHandler, TextEventHandler, TextVariablesHandler, YODAHistsHandler< T >, LHEFHandler, EventExporter

Public Functions

inline virtual void setCrossSection(const Value&)

Specify the cross section value, in pb.

inline void setEventNumber(unsigned long long ev_id)

Set event number.

virtual bool operator<<(const Event&) = 0

Writer operator.

void initialise(const RunParameters&)

Initialise the handler and its inner parameterisation.

const RunParameters &runParameters() const

List of run parameters.

template<typename T>
inline T *engine()

Retrieve the engine object.

inline const std::string &name() const

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const

Equality operator.

inline bool operator!=(const SteeredObject &oth) const

Inequality operator.

inline virtual const ParametersList &parameters() const override

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)

Set (documented) module parameters.

A full list of the output modules currently supported in CepGen addons, along with their user-steerable parameters, can be found here.


Event harvester#

Added in version 1.0.0.

class EventHarvester : public EventExporter#

Subclassed by GnuplotHarvester, MatplotlibHarvester, TextHarvester, TopdrawerHarvester

This simplest case of an output module allows to generate integrated histograms of kinematic variables, fully configurable by the user. Using the Python steering cards definition, a dictionary histVariables of variable-indexed cepgen.Parameters objects is fed to the output module.

A valid implementation of such objects requires a set of attributes depending on the type of distribution requested by the user. If the variable string contains one :, a 2D distribution is automatically booked for the two variables around it. Otherwise, a 1D distribution is assumed.

These attributes are, namely for 1-dimensional histograms:

  • a number of bins nbins, or nbinsX, and

  • a range (low and high) of interest for the variable, or a set of bins in a xbins Python list

and for 2-dimensional distributions:

  • the two nbinsX and nbinsY number of bins, and

  • the two ranges (lowX and highX, and lowY and highY) of interest for the variables, or equivalently one or two sets of bins in xbins/ybins lists.

As an example, equivalently to vars output defined above, the following output block may be used for a text output histogram with kinematics equivalent to the lpair process:

cepgen.Module('text',
    histVariables = {
       # 1D histogram (pt of central system)
       'pt(4)': cepgen.Parameters(nbins=10, xrange=(0., 20.)),
       # 1D histogram (outgoing proton mass)
       'm(5)': cepgen.Parameters(xbins=[float(bin) for bin in range(0, 20, 1)]),
       # 2D histogram (central system rapidity vs. mass)
       'y(cs):m(cs)': cepgen.Parameters(nbinsX=10, xrange=(-10., 10.),
                                        nbinxY=10, yrange=(0., 400.)),
    },
    save = False,
    show = True
)

To quote a few introduced since v1 of CepGen, the following plotters are handled for harvesters:

  • gnuplot, on systems which handle a working version of the gnuplot executable documented here, and for which the command line piper works,

  • matplotlib, in case the C++ wrapper of the matplotlib,

  • text, for the CepGen-specific terminal-based text histogram drawer, introduced as an overlay of the GSL histogramming capability.

Other types of output modules#

In this page you will find a list of all currently supported output formats, covering a broad spectrum of usages, both in the phenomenological and experimental communities. Please note that this list is under constant evolution, you may contact us with requests for additional interfacing capabilities.

dump#

A simple text-based event dumper, useful for debugging the process and its kinematics, is steered using the cepgen::TextEventHandler module.

class TextEventHandler : public EventExporter#

lhef#

This output format handles the conversion into the Les Houches standard definition. Currently, two implementations of this export module exist:

  • a Pythia 8 LHEF output module (described here) as the default handler, cepgen::LHEFPythiaHandler,

    Warning

    doxygenclass: Cannot find class “cepgen::LHEFPythiaHandler” in doxygen xml output for project “CepGen” from directory: /Documentation/build/xml

  • a HepMC (v≥3) implementation, cepgen::LHEFHepMCHandler.

    class LHEFHepMCHandler : public EventExporter#

hepmc2, hepmc2_ascii, …#

template<typename T>
class HepMC2Handler : public EventExporter#

This handler allows to translate the CepGen event record into one (or multiple) implementation(s) of the version 2 of the HepMC [DH01] ASCII output format. By default, this version is used in older releases. It allows a hepmc2 output format to be supported.

hepmc, hepmc_root, hepevt, …#

template<typename T>
class HepMC3Handler : public EventExporter#

This handler allows to translate the CepGen event record into one (or multiple) implementation(s) of the version 3 of the HepMC [DH01] ASCII output format.

By default, the version 3 of the file format is chosen for versions of HepMC starting from v3.1.0. It may be updated with future derivatives of the HepMC writer base class.

Alternatively, as from this version 3.1.0 of HepMC, the following output formats are also handled:

  • a hepevt ASCII format using the HepMC3::WriterHEPEVT handler,

  • a hepmc_root format using the HepMC3::WriterRoot export module,

  • a hepmc_root_tree using the HepMC3::WriterRootTree module.

promc#

Added in version 0.9.8.

class ProMCHandler : public EventExporter#

The support has been added for the ProMC highly compressed output format.

vars#

Added in version 1.0.0.

class TextVariablesHandler : public EventExporter#

This simplest case of an output module allows to generate a generic (ASCII) output format along with raw text histograms of kinematic variables, fully configurable by the user. Using the Python steering cards definition, a list of variables to be stored is defined through the variables list/array of string-typed definition.

For this text output format, the default behaviour is storing one event per line with variables separated with an user-parameterisable separator (separator string parameter, default is the standard tabulation \t).

The variable (here, var is used as an example) may be defined using the cepgen::utils::EventBrowser string-to-method conventions described in the Event content definition page. Two extra boolean parameters may also be fed to the module configuration:

  • saveBanner, to enable/disable the CepGen banner printout (containing useful information about the process and cuts definition), and

  • saveVariables, to show/hide the list of variables used in this file.

As an example, the following output block may be used for a 2-to-4 process such as lpair:

cepgen.Module('vars',
    filename = 'test.txt',
    variables = [
        'm(4)',    # central two-photon/central system mass
        'pt(cs)',  # central two-photon/central system transverse momentum
        'pt(6)'    # first outgoing central system particle transverse momentum
    ],
    saveBanner = False,
    saveVariables = True,
    separator = ' ',  # single space
)

root_hist, root_tree#

Added in version 0.9.7.

Note

Previously used in dedicated test executables, resp. test_distributions and cepgen-root.

These two modules module allow to produce a ROOT [BR97] file containing either:

  • a list of histograms (stored as ROOT TH1D objects) provided as an input for the earlier:

    class ROOTHistsHandler : public EventExporter#
  • or a set of events and run information (stored as ROOT TTree objects) for the latter:

    Warning

    doxygenclass: Cannot find class “cepgen::ROOTTreeHandler” in doxygen xml output for project “CepGen” from directory: /Documentation/build/xml

The histogramming utilitary follows the same procedure as introduced for the cepgen::TextHandler module above to define the histograms list.

As an example, the following output block may be used:

cepgen.Module('root_hist',
    filename = 'output.hists.root',
    variables = {
       'pt(4)': cepgen.Parameters(nbins=10, xrange=(0., 20.)),
       'm(5)': cepgen.Parameters(nbins=10, xrange=(0., 100.)),
       'y(cs)': cepgen.Parameters(nbins=10, xrange=(-10., 10.)),
    },
)

The tree handler may be used in parallel to the two ROOT::CepGenRun and ROOT::CepGenEvent helper reader objects for a compact analysis workflow:

class CepGenRun#

All useful information about a generation run.

Public Functions

void clear()#

Reinitialise the run tree.

void create()#

Populate the run tree.

inline TTree *tree()#

Retrieve the ROOT tree.

void fill()#

Fill the run tree.

void attach(TFile *file, const std::string &run_tree = TREE_NAME)#

Attach the run tree reader to a given tree Attach the run tree reader to a given file.

Public Members

double sqrt_s = {-1.}#

Centre of mass energy for beam particles.

double xsect = {-1.}#

Process cross section, in pb.

double errxsect = {-1.}#

Uncertainty on process cross section, in pb.

unsigned int num_events = {0}#

Number of events generated in run.

unsigned int litigious_events = {0}#

Number of litigious events in run.

std::string process_name#

Unique name of the process generated in this run.

std::string process_parameters#

Serialised process parameters.

Public Static Attributes

static constexpr const char *TREE_NAME = "run"#

Output tree name.

class CepGenEvent#

All useful information about a generated event.

Public Functions

void clear()#

Reinitialise the event content.

inline TTree *tree()#

Retrieve the ROOT tree.

void create()#

Populate the tree and all associated branches.

void attach()#

Attach the event tree reader to a tree.

void attach(TFile *f, const std::string &events_tree = TREE_NAME)#

Attach the event tree reader to a file Attach the event tree reader to a file.

void fill(const cepgen::Event&, bool compress = false)#

Fill the tree with a new event.

bool next(cepgen::Event&)#

Read the next event in the file.

Public Members

float gen_time = {-1.}#

Event generation time.

float tot_time = {-1.}#

Total event generation time.

float weight = {-1.}#

Event weight.

int np = {0}#

Number of particles in the event.

double pt[MAX_PART]#

Particles transverse momentum.

double eta[MAX_PART]#

Particles pseudo-rapidity.

double phi[MAX_PART]#

Particles azimuthal angle.

double rapidity[MAX_PART]#

Particles rapidity.

double E[MAX_PART]#

Particles energy, in GeV.

double m[MAX_PART]#

Particles mass, in GeV/c \({}^2\).

double charge[MAX_PART]#

Particles charges, in e.

int pdg_id[MAX_PART]#

Integer particles PDG id.

int parent1[MAX_PART]#

First particles mother.

int parent2[MAX_PART]#

Last particles mother.

int stable[MAX_PART]#

Whether the particle must decay or not.

int role[MAX_PART]#

Particles role in the event.

int status[MAX_PART]#

Integer status code.

Public Static Attributes

static constexpr size_t MAX_PART = 5000#

Maximal number of particles in event.

static constexpr const char *TREE_NAME = "events"#

Output tree name.

delphes#

Added in version 0.9.7.

An interface to the Delphes [dFDD+14] fast simulation framework is provided through the CepGenDelphes add-on implemented here.

Beside the usual filename flag specifying the file name Delphes will use for its output, a path to the Tcl configuration card is also required to steer the output module through the inputCard string parameter.

Please refer to the Delphes manual and comprehensive list of examples for more information on the steering of the detector simulation.