Les Houches Accord

The Les Houches Accord (LHA) for user processes Boo01 is the standard way to input parton-level information from a matrix-elements-based generator into PYTHIA. The conventions for which information should be stored has been defined in a Fortran context, as two commonblocks. Here a C++ equivalent is defined, as a single class.

The LHAup class is a base class, containing reading and printout functions, plus two pure virtual functions, one to set initialization information and one to set information on each new event. Derived classes have to provide these two virtual functions to do the actual work. The existing derived classes are for reading information from a Les Houches Event File (LHEF), from the respective Fortran commonblocks, or from PYTHIA 8 itself.

You are free to write your own derived classes, using the rules and methods to be described below. Normally, pointers to objects of such derived classes should be handed in with the Pythia::init( LHAup*) method. However, with the LHEF format a filename can replace the pointer, see further below.

Let us now describe the methods at your disposal to do the job. the base class constructor takes the choice of mixing/weighting strategy as optional input argument, and calls setStrategy, see below. It also reserves some space for processes and particles. the destructor does not need to do anything. this method only sets the pointer that allows some information to be accessed, and is automatically called by Pythia::init(...).

Initialization

The LHAup class stores information equivalent to the /HEPRUP/ commonblock, as required to initialize the event generation chain. The main difference is that the vector container now allows a flexible number of subprocesses to be defined. For the rest, names have been modified, since the 6-character-limit does not apply, and variables have been regrouped for clarity, but nothing fundamental is changed. this pure virtual method has to be implemented in the derived class, to set relevant information when called. It should return false if it fails to set the info.

Inside setInit(), such information can be set by the following methods: sets the properties of the first and second incoming beam, respectively (cf. the Fortran IDBMUP(1), EBMUP(i), PDFGUP(i), PDFSUP(i), with i 1 or 2). The parton distribution information defaults to zero. These numbers can be used to tell which PDF sets were used when the hard process was generated, while the normal PDF Selection is used for the further event generation in PYTHIA. sets the event weighting and cross section strategy. The default, provided in the class constructor, is 3, which is the natural value e.g. for an LHEF. chosen strategy (cf. IDWTUP; see Sjo06 section 9.9.1 for extensive comments). events come with non-negative weight, given in units of pb, with an average that converges towards the cross section of the process. PYTHIA is in charge of the event mixing, i.e. for each new try decides which process should be generated, and then decides whether is should be kept, based on a comparison with xMax. Accepted events therefore have unit weight. as option 1, except that cross sections can now be negative and events after unweighting have weight +-1. You can use Info::weight() to find the weight of the current event. A correct event mixing requires that a process that can take both signs should be split in two, one limited to positive or zero and the other to negative or zero values, with xMax chosen appropriately for the two. events come with non-negative weight, in unspecified units, but such that xMax can be used to unweight the events to unit weight. Again PYTHIA is in charge of the event mixing. The total cross section of a process is stored in xSec. as option 2, except that cross sections can now be negative and events after unweighting have weight +-1. As for option -1 processes with indeterminate sign should be split in two. events come with unit weight, and are thus accepted as is. The total cross section of the process is stored in xSec. as option 3, except that events now come with weight +-1. Unlike options -1 and -2 processes with indeterminate sign need not be split in two, unless you intend to mix with internal PYTHIA processes (see below). events come with non-negative weight, given in units of pb, with an average that converges towards the cross section of the process, like for option 1. No attempt is made to unweight the events, however, but all are generated in full, and retain their original weight. For consistency with normal PYTHIA units, the weight stored in Info::weight() has been converted to mb, however. as option 4, except that events now can come either with positive or negative weights. Note 1: if several processes have already been mixed and stored in a common event file, either LHEF or some private format, it would be problematical to read back events in a different order. Since it is then not feasible to let PYTHIA pick the next process type, strategies +-1 and +-2 would not work. Instead strategy 3 would be the recommended choice, or -3 if negative-weight events are required. Note 2: it is possible to switch on internally implemented processes and have PYTHIA mix these with LHA ones according to their relative cross sections for strategies +-1, +-2 and 3. It does not work for strategy -3 unless the positive and negative sectors of the cross sections are in separate subprocesses (as must always be the case for -1 and -2), since otherwise the overall mixture of PYTHIA and LHA processes will be off. Mixing is not possible for strategies +-4, since the weighting procedure is not specified by the standard. (For instance, the intention may be to have events biased towards larger pT values in some particular functional form.) sets info on an allowed process (cf. LPRUP, XSECUP, XERRUP, XMAXUP). Each new call will append one more entry to the list of processes. The choice of strategy determines which quantities are mandatory: xSec for strategies +-2 and +-3, xErr never, and xMax for strategies +-1 and +-2. Note: PYTHIA does not make active use of the (optional) xErr values, but calculates a statistical cross section error based on the spread of event-to-event weights. This should work fine for strategy options +-1, but not for the others. Specifically, for options +-2 and +-3 the weight spread may well vanish, and anyway is likely to be an underestimate of the true error. If the author of the LHA input information does provide error information you may use that - this information is displayed at initialization. If not, then a relative error decreasing like 1/sqrt(n_acc), where n_acc is the number of accepted events, should offer a reasonable estimate. update the xSec value of the i'th process added with addProcess method (i.e. i runs from 0 through sizeProc() - 1, see below). update the xErr value of the i'th process added with addProcess method. update the xMax value of the i'th process added with addProcess method. set the header key to have value val. This is a wrapper function to the Info::setHeader function that should be used in any classes derived from LHAup.

Information is handed back by the following methods (that normally you would not need to touch): for the beam properties. for the strategy choice. for the number of subprocesses. for process i in the range 0 <= i < sizeProc(). the sum of the cross sections and errors (the latter added quadratically). Note that cross section errors are only meaningful for strategies +-3. prints the above initialization information. This method is automatically called from Pythia::init(...), so would normally not need to be called directly by the user.

Event input

The LHAup class also stores information equivalent to the /HEPEUP/ commonblock, as required to hand in the next parton-level configuration for complete event generation. The main difference is that the vector container now allows a flexible number of partons to be defined. For the rest, names have been modified, since the 6-character-limit does not apply, and variables have been regrouped for clarity, but nothing fundamental is changed.

The LHA standard is based on Fortran arrays beginning with index 1, and mother information is defined accordingly. In order to be compatible with this convention, the zeroth line of the C++ particle array is kept empty, so that index 1 also here corresponds to the first particle. One small incompatibility is that the sizePart() method returns the full size of the particle array, including the empty zeroth line, and thus is one larger than the true number of particles (NUP). this pure virtual method has to be implemented in the derived class, to set relevant information when called. For strategy options +-1 and +-2 the input idProcess value specifies which process that should be generated, while idProcess is irrelevant for strategies +-3 and +-4. The method should return false if it fails to set the info, i.e. normally that the supply of events in a file is exhausted. If so, no event is generated, and Pythia::next() returns false. You can then interrogate Info::atEndOfFile() to confirm that indeed the failure is caused in this method, and decide to break out of the event generation loop.

Inside a normal setEvent(...) call, information can be set by the following methods: tells which kind of process occured, with what weight, at what scale, and which alpha_EM and alpha_strong were used (cf. IDPRUP, XWTGUP, SCALUP, AQEDUP, AQCDUP). This method also resets the size of the particle list, and adds the empty zeroth line, so it has to be called before the addParticle method below. gives the properties of the next particle handed in (cf. IDUP, ISTUP, MOTHUP(1,..), MOTHUP(2,..), ICOLUP(1,..), ICOLUP(2,..), PUP(J,..), VTIMUP, SPINUP) .

Information is handed back by the following methods: process number. . Note that the weight stored in Info::weight() as a rule is not the same as the above weight(): the method here gives the value before unweighting while the one in info gives the one after unweighting and thus normally is 1 or -1. Only with strategy options +-3 and +-4 would the value in info be the same as here, except for a conversion from pb to mb for +-4. scale and couplings at that scale. the size of the particle array, which is one larger than the number of particles in the event, since the zeroth entry is kept empty (see above). for particle i in the range 0 <= i < sizePart(). (But again note that i = 0 is an empty line, so the true range begins at 1.)

From the information in the event record it is possible to set the flavour and x values of the initiators

This information is returned by the methods x values of the two initiators.

In the LHEF description Alw06 an extension to include information on the parton densities of the colliding partons is suggested. This optional further information can be set by which gives the flavours , the x and the Q scale (in GeV) at which the parton densities x*f_i(x, Q) have been evaluated. The last argument is normally true.

This information is returned by the methods where the first one tells whether this optional information has been set for the current event. (setPdf(...) must be called after the setProcess(...) call of the event for this to work.) Note that the flavour and x values usually but not always agree with those obtained by the same methods without pdf in their names, see explanation in the Event Information description.

prints the above information for the current event. In cases where the LHAup object is not available to the user, the Pythia::LHAeventList(ostream& os = cout) method can be used, which is a wrapper for the above. skip ahead nSkip events in the Les Houches generation sequence, without doing anything further with them. Mainly intended for debug purposes, e.g. when an event at a known location in a Les Houches Event File is causing problems. Will return false if operation fails, specifically if the end of an LHEF has been reached. The implementation in the base class simply executes setEvent() the requested number of times. The derived LHAupLHEF class (see below) only uses the setNewEventLHEF(...) part of its setEvent() method, and other derived classes could choose other shortcuts.

The LHA expects the decay of resonances to be included as part of the hard process, i.e. if unstable particles are produced in a process then their decays are also described. This includes Z^0, W^+-, H^0 and other short-lived particles in models beyond the Standard Model. Should this not be the case then PYTHIA will perform the decays of all resonances it knows how to do, in the same way as for internal processes. Note that you will be on slippery ground if you then restrict the decay of these resonances to specific allowed channels since, if this is not what was intended, you will obtain the wrong cross section and potentially the wrong mix of different event types. (Since the original intention is unknown, the cross section will not be corrected for the fraction of open channels, i.e. the procedure used for internal processes is not applied in this case.)

Even if PYTHIA can select resonance decay modes according to its internal tables, there is normally no way for it to know which decay angular correlations should exist in the simulated process. Therefore almost all decays are isotropic. The exceptions are Higgs and top decays, in the decay chains H -> WW/ZZ -> f fbar f' fbar' and t -> b W -> b f fbar, where the process-independent correlations implemented for internal processes are used. If part of the decay chain has already been set, however (e.g. H -> WW/ZZ or t -> b W), then decay is still isotropic.

An interface to Les Houches Event Files

The LHEF standard Alw06 specifies a format where a single file packs initialization and event information. This has become the most frequently used procedure to process external parton-level events in Pythia. Therefore a special Pythia::init(fileName) initialization option exists, where the LHEF name is provided as input. Internally this name is then used to create an instance of the derived class LHAupLHEF, which can do the job of reading an LHEF.

The LHEF reader can also read in and store header blocks. By default this option is switched on, but may be controlled through the Beams:readLHEFheaders flag if necessary. The information can later be read out through the Info class for further processing. Due to the non-standard nature of the information in these blocks they are stored whole, and PYTHIA itself makes no further attempt to process their meaning.

Because Les Houches Event files tend not to adhere strictly to XML conventions, to consistently read in header information, certain choices must be made. The primary goal is to make as much information available as possible. First, information sitting directly in the <header> block is stored under the key "base". Second, the tags starting and ending each sub block must be on their own line. Finally, the contents of comment blocks, <!-- -->, are still stored. The header keys are formed hierarchically from the names of the header blocks. This behaviour is illustrated in the following example:

  <header>
    BaseA
    <hblock1>
      1A
      <hblock11>
        11A <hblock111>
        </hblock111> 11B
      </hblock11>
      1B
    </hblock1>
    <hblock2>
      2A
      <!-- 2B -->
    </hblock2>
    BaseB
  </header>
which would lead to the following information being stored in the Info class:
Key Value
base BaseA
BaseB
hblock1 1A
1B
hblock1.hblock11 11A <hblock111>
</hblock111> 11B
hblock2 2A
<!-- 2B -->

Normally the LHEF would be in uncompressed format, and thus human-readable if opened in a text editor. A possibility to read gzipped files has been added, based on the Boost and zlib libraries, which therefore have to be linked appropriately in order for this option to work. See the README file in the main directory for details on how to do this.

An example how to generate events from an LHEF is found in main11.cc. Note the use of Info::atEndOfFile() to find out when the whole LHEF has been processed.

To allow the sequential use of several event files the Pythia::init(...) method has an optional second argument: Pythia::init(fileName, bool skipInit = false). If called with this argument true then there will be no initialization, except that the existing LHAupLHEF class instance will be deleted and replaced by ones pointing to the new file. It is assumed (but never checked) that the initialization information is identical, and that the new file simply contains further events of exactly the same kind as the previous one. An example of this possibility, and the option to mix with internal processes, is found in main12.cc. A variant, based on input in a command file, is given in main13.cc.

In C++, real numbers are printed with an 'E' to denote the exponent part, e.g. 1.23E+04, and are read in accordingly. Other languges may use other letters, e.g. Fortran allows either 'E' or 'D'. A file using the latter convention would not be readable by the standard routines. In case you have such an "incorrectly formatted" file, a conversion to a new corrected file could be done e.g. using sed, as a one-line command

  sed -e 's/\([0-9]\.\{0,1\}\)[dD]\([+-]\{0,1\}[0-9]\)/\1E\2/g' old.lhe > new.lhe
This replaces a 'd' or 'D' with an 'E' only when it occurs in the combination
(digit) ('.' or absent) ('d' or 'D') ('+', '-' or absent) (digit)
It will work on all parts of the file, also inside a <header>...</header> block. For conversion only inside the <init>...</init> and <event>...</event> blocks, create a file convert.sed containing
  /<init>/,/<\/init>/bconv
  /<event>/,/<\/event>/bconv
  b
  :conv
  s/\([0-9]\.\{0,1\}\)[dD]\([+-]\{0,1\}[0-9]\)/\1E\2/g
and run it with
  sed -f convert.sed old.lhe > new.lhe

The workhorses of the LHAupLHEF class are three methods found in the base class, so as to allow them to be reused in other contexts. Specifically, it allows derived classes where one parton-level configuration can be reused several times, e.g. in the context of matrix-element-to-parton-shower matching (example in preparation). Also two small utility routines. read in and set all required initialization information from the specified stream. With second argument true it will also read and store header information, as described above. Return false if it fails. read in event information from the specified stream into a staging area where it can be reused by setOldEventLHEF. store the event information from the staging area into the normal location. Thus a single setNewEventLHEF call can be followed by several setOldEventLHEF ones, so as to process the same configuration several times. This method currently only returns true, i.e. any errors should be caught by the preceding setNewEventLHEF call. always returns true in the base class, but in LHAupLHEF it returns false if the LHEF provided in the constructor is not found and opened correctly. is used to send header information on to the Info class.

A few other methods, most of them derived from the base class, streamlines file opening and closing, e.g. if several LHE files are to be read consecutively, without the need for a complete reinitialization. This presupposes that the events are of the same kind, only split e.g. to limit file sizes. close current event input file/stream and open a new one, to continue reading events of the same kind as before. open and close a file, also gzip files, where an intermediate decompression layer is needed. close main event file (LHEF) and, if present, separate header file.

A runtime Fortran interface

The runtime Fortran interface requires linking to an external Fortran code. In order to avoid problems with unresolved external references when this interface is not used, the code has been put in a separate LHAFortran.h file, that is not included in any of the other library files. Instead it should be included in the user-supplied main program, together with the implementation of two methods below that call the Fortran program to do its part of the job.

The LHAupFortran class derives from LHAup. It reads initialization and event information from the LHA standard Fortran commonblocks, assuming these commonblocks behave like two extern "C" struct named heprup_ and hepeup_. (Note the final underscore, to match how the gcc compiler internally names Fortran files.)

The instantiation does not require any arguments.

The user has to supply implementations of the fillHepRup() and fillHepEup() methods, that is to do the actual calling of the external Fortran routines that fill the HEPRUP and HEPEUP commonblocks. The translation of this information to the C++ structure is provided by the existing setInit() and setEvent() code.

Up to and including version 8.125 the LHAupFortran class was used to construct a runtime interface to PYTHIA 6.4. This was convenient in the early days of PYTHIA 8 evolution, when this program did not yet contain hard-process generation, and the LHEF standard did not yet exist. Nowadays it is more of a bother, since a full cross-platform support leads to many possible combinations. Therefore the support has been reduced in the current version. Only the main91.cc example remains as an illustration, where the previously separate interface code (include/Pythia6Interface.h) has been inserted in the beginning. You also need to modify the examples/Makefile to link main91.cc properly also to a PYTHIA 6.4 library version, see commented-out section for ideas how to to this.

Methods for LHEF output

The main objective of the LHAup class is to feed information from an external program into PYTHIA. It can be used to export information as well, however. Specifically, there are four routines in the base class that can be called to write a Les Houches Event File. These should be called in sequence in order to build up the proper file structure. Opens a file with the filename indicated, and writes a header plus a brief comment with date and time information. Writes initialization information to the file above. Such information should already have been set with the methods described in the "Initialization" section above. Writes event information to the file above. Such information should already have been set with the methods described in the "Event input" section above. This call should be repeated once for each event to be stored. By default the event information is lined up in columns. To save space, the alternative verbose = false only leaves a single blank between the information fields. Writes the closing tag and closes the file. Optionally, if updateInit = true, this routine will reopen the file from the beginning, rewrite the same header as openLHEF() did, and then call initLHEF() again to overwrite the old information. This is especially geared towards programs, such as PYTHIA itself, where the cross section information is not available at the beginning of the run, but only is obtained by Monte Carlo integration in parallel with the event generation itself. Then the setXSec( i, xSec), setXErr( i, xSec) and setXMax( i, xSec) can be used to update the relevant information before closeLHEF is called. Warning: overwriting the beginning of a file without upsetting anything is a delicate operation. It only works when the new lines require exactly as much space as the old ones did. Thus, if you add another process in between, the file will be corrupted.

PYTHIA 8 output to an LHEF

The above methods could be used by any program to write an LHEF. For PYTHIA 8 to do this, a derived class already exists, LHAupFromPYTHIA8. In order for it to do its job, it must gain access to the information produced by PYTHIA, specifically the process event record and the generic information stored in info. Therefore, if you are working with an instance pythia of the Pythia class, you have to instantiate LHAupFromPYTHIA8 with pointers to the process and info objects of pythia:
LHAupFromPYTHIA8 myLHA(&pythia.process, &pythia.info);

The method setInit() should be called to store the pythia initialization information in the LHA object, and setEvent() to store event information. Furthermore, updateSigma() can be used at the end of the run to update cross-section information, cf. closeLHEF(true) above. An example how the generation, translation and writing methods should be ordered is found in main20.cc.

Currently there are some limitations, that could be overcome if necessary. Firstly, you may mix many processes in the same run, but the cross-section information stored in info only refers to the sum of them all, and therefore they are all classified as a common process 9999. Secondly, you should generate your events in the CM frame of the collision, since this is the assumed frame of stored Les Houches events, and no boosts have been implemented for the case that Pythia::process is not in this frame.

The LHEF standard is the agreed format to store the particles of a hard process, as input to generators, whereas output of final states is normally handled using the HepMC standard. It is possible to use LHEF also here, however. It requires that the above initialization is replaced by
LHAupFromPYTHIA8 myLHA(&pythia.event, &pythia.info);
i.e. that process is replaced by event. In addition, the PartonLevel:all = off command found in main20.cc obviously must be removed if one wants to obtain complete events.