HepMC Interface

An interface to the HepMC [Dob01] standard event record format has been provided by M. Kirsanov. To use it, the relevant libraries need to be linked, as explained in the README file. Only version 2 of HepMC is supported. (Version 1 requires a different interface structure, which was only supported up until Pythia 8.107.)

The (simple) procedure to translate PYTHIA 8 events into HepMC ones is illustrated in the main41.cc, main42.cc main61.cc and main62.cc main programs. At the core is a call to the

HepMC::I_Pythia8::fill_next_event( pythia, hepmcevt, ievnum = -1, convertGluonTo0 = false ) 
which takes a reference of the generator object and uses it, on the one hand, to read out and covert the event record in pythia.event and, on the other hand, to extract and store parton-density (PDF) information for the hard subprocess from pythia.info. The optional last argument, if true, allows you to store gluons as "PDG" code 0 rather than the normal 21; this only applies to the PDF information, not the event record.

The earlier version of this routine,

HepMC::I_Pythia8::fill_next_event( pythia.event, hepmcevt, ievnum = -1 ) 
is retained (for now) for backwards compatibility. It takes a PYTHIA event as input and returns a HepMC one, but without storing the PDF information. The latter could then instead be stored by a separate call
HepMC::I_Pythia8::pdf_put_info( hepmcevt, pythia, convertGluonTo0 = false ) 
or not, as wished.

The translation routine stores momenta, energies and masses in units of GeV, and distances and times in units of mm, with c = 1, exactly as used in PYTHIA. This only works seamlessly for HepMC 2.04 or later, where the choice of units can (and should) be specified in the HepMC::GenEvent call. To work also with older HepMC versions, where normally momenta are stored in MeV, the HEPMC_HAS_UNITS environment variable is tested. If not defined, i.e. if with an earlier HepMC version, conversion from GeV to MeV is done in fill_next_event.

The status code is now based on the proposed standard for HepMC 2.05, see the Event::statusHepMC(...) conversion routine for details. The earlier behaviour, where all final particles had status 1 and all initial or intermediate ones status 2, is available as a commented-out line in the I_Pythia8::fill_next_event(...) method, and so is simple to recover.

In HepMC 2.05 it also becomes possible to store the generated cross section and its error. The environment variable HEPMC_HAS_CROSS_SECTION is used to check whether this possibility exists and, if it does 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 it 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, conversion from mb to fb for that weight must be set by hand, see the last method below.)

The public methods

Here comes a complete list of all public methods of the I_Pythia8 class in the HepMC (not Pythia8!) namespace.

I_Pythia8::I_Pythia8()  
virtual I_Pythia8::~I_Pythia8()  
the constructor and destructor take no arguments.

bool I_Pythia8::fill_next_event( Pythia8::Pythia& pythia, GenEvent* evt, int ievnum = -1, bool convertGluonTo0 = false)  
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 the parton density 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 convertGluonTo0 : 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.

bool I_Pythia8::fill_next_event( Pythia8::Event& pyev, GenEvent* evt, int ievnum = -1 )  
convert a Pythia event to a HepMC one. Will return true if succeeded. Do not store parton-density information.
argument pyev : the input Pythia event, in its standard form.
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.

bool I_Pythia8::put_pdf_info( GenEvent* evt, Pythia8::Pythia& pythia, bool convertGluonTo0 = false )  
append parton-density information to an event already stored by the previous method.
argument evt : the output GenEvt event record, in its standard form.
argument pythia : the input Pythia generator object, from which both the event and the parton density information can be obtained.
argument convertGluonTo0 : 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.

The following methods can be used to set, and in some cases interrogate, the status of some switches that can be used to 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 I_Pythia8::set_trust_mothers_before_daughters( bool b = true )  
bool I_Pythia8::trust_mothers_before_daughters()  
if there is a conflict in the history information, then trust the information on mothers above that on daughters. Currently this is the only option implemented.

void I_Pythia8::set_trust_both_mothers_and_daughters( bool b = false )  
bool I_Pythia8::trust_both_mothers_and_daughters()  
currently dummy methods intended to resolve conflicts in the event history.

void I_Pythia8::set_print_inconsistency_errors( bool b = true )  
bool I_Pythia8::print_inconsistency_errors()  
print a warning line, on cerr, when inconsistent mother and daughter information is encountered.

void I_Pythia8::set_crash_on_problem( bool b = false )  
if problems are encountered then the run is interrupted by an exit(1) command. Default is not to crash.

void I_Pythia8::set_freepartonwarnings( bool b = true )  
interrupt the run by an exit(1) command if unhadronized gluons or quarks are encountered in the event record, unless hadronization is switched off. Default is to crash.

void I_Pythia8::set_convert_to_mev( bool b = false )  
convert the normal GeV energies, momenta and masses to MeV.