HepMC Interface
An interface to the HepMC [Dob01] standard event record
format has been provided by M. Kirsanov. The code is stored in
include/Pythia8Plugins/HepMC2.h
. To use it,
the relevant libraries need to be linked, as explained in the
README
file.
Only version 2.06 (or later) of HepMC is supported, by agreement
with the LHC experimental community.
The (simple) procedure to translate PYTHIA 8 events into HepMC ones
is illustrated in the main41.cc
, main42.cc
and main43.cc
main programs. At the core is a call to the
HepMC::Pythia8ToHepMC::fill_next_event( pythia, hepmcevt, ievnum = -1)
which takes a reference of the generator object and uses it, on the one
hand, to read out and convert the event record in pythia.event
and, on the other hand, to extract and store parton-density (PDF),
cross section and other information for the hard subprocess from
pythia.info
. There is also an alternative form that
does not requires access to the full pythia
object,
but only the event record, at the expense of a reduced information
storage, see below.
While PYTHIA always stores momenta in units of GeV, with c = 1,
HepMC nowadays can be built either for MeV or GeV as default, a choice
that can then be overridden on an event-by-event basis, see e.g. the
main41.cc
code. When filling the HepMC event record, PYTHIA
will convert to the unit specified for the current HepMC event record.
Also for length units there are choices, and again the PYTHIA interface
will convert to the units set for the HepMC event record. Here the mm
choice of PYTHIA seems to be shared by most other programs, and is
HepMC default.
The status code is now based on the new standard introduced for HepMC 2.05,
see the Event::statusHepMC(...)
conversion routine for details.
The current values from pythia.info.sigmaGen()
and
pythia.info.sigmaErr()
are stored for each event,
multiplied by 10^9 to convert from mb to pb. Note that
PYTHIA improves its accuracy by Monte Carlo integration in the course
of the run, so the values associated with the last generated event
should be the most accurate ones. If events also come with a dimensional
weight, like in some Les Houches strategies, this weight is in units of pb.
The public methods
Here comes a complete list of all public methods of the
Pythia8ToHepMC
class in the HepMC
(not Pythia8
!) namespace.
Pythia8ToHepMC::Pythia8ToHepMC()
virtual Pythia8ToHepMC::~Pythia8ToHepMC()
the constructor and destructor take no arguments.
bool Pythia8ToHepMC::fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt, int ievnum = -1, bool append = false, GenParticle* rootParticle = 0, int iBarcode = -1)
convert a Pythia
event to a HepMC
one.
Will return true if succeeded.
argument
pythia :
the input Pythia
generator object, from which both the
event and other information can be obtained.
argument
evt :
the output GenEvt
event, in its standard form.
argument
iev :
set the event number of the current event. If negative then the
internal event number is used, which is incremented by one for
each new event.
argument
append :
if true
then the input event is appended to the current
HepMC
event record, rather than starting a new one.
argument
rootParticle :
the root particle that defines the new production vertex for the
particles to be added in the append = true
option.
argument
iBarcode :
used to set the bar code when append = true
.
If positive then start from iBarcode
itself, if negative
then start from the current size of the HepMC
event record,
and if 0 then set all bar codes to vanish.
bool Pythia8ToHepMC::fill_next_event( Pythia8::Event& pyev, GenEvent* evt, int ievnum = -1, Pythia8::Info* pyinfo = 0, Pythia8::Settings* pyset = 0, bool append = false, GenParticle* rootParticle = 0, int iBarcode = -1)
convert a Pythia
event to a HepMC
one.
Will return true if succeeded.
argument
pyev :
the input Pythia
event that is to be converted to HepMC
format.
argument
evt :
the output GenEvt
event, in its standard form.
argument
iev :
set the event number of the current event. If negative then the
internal event number is used, which is incremented by one for
each new event.
argument
pyinfo :
pointer to the Pythia Info
object, which is used to
extract PFD values, and process and cross section information.
Without such a pointer this information therefore cannot be stored,
i.e. it is equivalent to the three set_store
methods
below being set false.
argument
pyset :
pointer to the Pythia Settings
object, which is used to
decide whether hadronization is carried out, and therefore whether
to warn about unhadronized partons. Without such a pointer the
set_free_parton_warnings
method below uniquely controls
the behaviour.
argument
append :
if true
then the input event is appended to the current
HepMC
event record, rather than starting a new one.
argument
rootParticle :
the root particle that defines the new production vertex for the
particles to be added in the append = true
option.
argument
iBarcode :
used to set the bar code when append = true
.
If positive then start from iBarcode
itself, if negative
then start from the current size of the HepMC
event record,
and if 0 then set all bar codes to vanish.
The following paired methods can be used to set and get the status of
some switches that modify the behaviour of the conversion routine.
The set
methods have the same default input values as
the internal initialization ones, so that a call without an argument
(re)stores the default.
void Pythia8ToHepMC::set_print_inconsistency(bool b = true)
bool Pythia8ToHepMC::print_inconsistency()
print a warning line, on cerr
, when inconsistent mother
and daughter information is encountered.
void Pythia8ToHepMC::set_free_parton_exception(bool b = true)
bool Pythia8ToHepMC::free_parton_exception()
check and throw an exception when unhadronized gluons or quarks are
encountered in the event record. Does not apply when Pythia hadronization
is switched off. Default is to do this check. If an exception is thrown the
PartonEndVertexException
class will return a warning message.
The calling code can choose action to take, also having access to the
location (index()
) and species (pdg_id()
) of a
bad parton.
void Pythia8ToHepMC::set_convert_gluon_to_0(bool b = false)
bool Pythia8ToHepMC::convert_gluon_to_0()
the normal gluon identity code 21 is used also when parton density
information is stored, unless this optional argument is set true to
have gluons represented by a 0. This choice does not affect the
normal event record, where a gluon is always 21.
void Pythia8ToHepMC::set_store_pdf(bool b = true)
bool Pythia8ToHepMC::store_pdf()
for each event store information on the two incoming flavours, their
x values and common factorization scale, and the values of the two
parton distributions, xf(x,Q).
void Pythia8ToHepMC::set_store_proc(bool b = true)
bool Pythia8ToHepMC::store_proc()
for each event store information on the Pythia process code, the
renormalization scale, and alpha_em and alpha_s
values used for the hard process.
void Pythia8ToHepMC::set_store_xsec(bool b = true)
bool Pythia8ToHepMC::store_xsec()
for each event store information on the Pythia cross section and its error,
in pb, and the event weight. If events also come with a dimensional weight,
like in some Les Houches strategies, this weight is in units of pb.