Les Houches Accord
- Initialization
- Event input
- A LHEF reader class
- A runtime Fortran interface
- A LHEF writer class
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 most common application is to read input from a Les Houches Event File
(LHEF) [Alw06], but it is also possible to have a runtime
interface to another program. This page will discuss runtime interface. For
LHEF handling, please refer to
Les Houches Event Files.
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, 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::setLHAupPtr( 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.
LHAup::LHAup( int strategy = 3)
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.
virtual LHAup::~LHAup()
the destructor does not need to do anything.
void LHAup::setPtr(Info* infoPtr)
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.
virtual bool LHAup::setInit()
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. In the no-beams extension this method need not
do anything, since by default strategy 3 is chosen and the rest is set
vanishing, but the method must exist.
Inside setInit()
, such information can be set by the following
methods:
void LHAup::setBeamA( int identity, double energy, int pdfGroup = 0, int pdfSet = 0)
void LHAup::setBeamB( int identity, double energy, int pdfGroup = 0, int pdfSet = 0)
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). 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.
void LHAup::setStrategy( int strategy)
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.
argument
strategy :
chosen strategy (cf. IDWTUP
; see [Sjo06]
section 9.9.1 for extensive comments).
argumentoption
1 : 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.
argumentoption
-1 : 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.
argumentoption
2 : 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
.
argumentoption
-2 : 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.
argumentoption
3 : events come with unit weight, and are thus accepted
as is. The total cross section of the process is stored in
xSec
.
argumentoption
-3 : 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).
argumentoption
4 : 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.
argumentoption
-4 : 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.)
void LHAup::addProcess( int idProcess, double xSec, double xErr, double xMax)
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.
void LHAup::setXSec( int i, double xSec)
update the xSec
value of the i
'th process
added with addProcess
method (i.e. i
runs
from 0 through sizeProc() - 1
, see below).
void LHAup::setXErr( int i, double xErr)
update the xErr
value of the i
'th process
added with addProcess
method.
void LHAup::setXMax( int i, double xMax)
update the xMax
value of the i
'th process
added with addProcess
method.
void LHAup::setInfoHeader(string &key, string &val)
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):
int LHAup::idBeamA()
int LHAup::idBeamB()
double LHAup::eBeamA()
double LHAup::eBeamB()
int LHAup::pdfGroupBeamA()
int LHAup::pdfGroupBeamB()
int LHAup::pdfSetBeamA()
int LHAup::pdfSetBeamB()
for the beam properties.
int LHAup::strategy()
for the strategy choice.
int LHAup::sizeProc()
for the number of subprocesses.
int LHAup::idProcess(int i)
double LHAup::xSec(int i)
double LHAup::xErr(int i)
double LHAup::xMax(int i)
for process i
in the range 0 <= i <
sizeProc()
.
double LHAup::xSecSum()
double LHAup::xErrSum()
the sum of the cross sections and errors (the latter added quadratically).
Note that cross section errors are only meaningful for strategies +-3.
void LHAup::listInit()
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
).
virtual bool LHAup::setEvent(int idProcess = 0)
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:
void LHAup::setProcess( int idProcess, double weight, double scale, double alphaQED, double alphaQCD)
tells which kind of process occurred, 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.
void LHAup::addParticle( int id, int status, int mother1, int mother2, int colourTag1, int colourTag2, double p_x, double p_y, double p_z, double e, double m, double tau, double spin, double scale)
gives the properties of the next particle handed in (cf. IDUP, ISTUP,
MOTHUP(1,..), MOTHUP(2,..), ICOLUP(1,..), ICOLUP(2,..), PUP(J,..),
VTIMUP, SPINUP
; while scale
is a new optional property,
see further below) .
Information is handed back by the following methods:
int LHAup::idProcess()
process number.
double LHAup::weight()
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.
double LHAup::scale()
double LHAup::alphaQED()
double LHAup::alphaQCD()
scale and couplings at that scale.
int LHAup::sizePart()
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).
int LHAup::id(int i)
int LHAup::status(int i)
int LHAup::mother1(int i)
int LHAup::mother2(int i)
int LHAup::col1(int i)
int LHAup::col2(int i)
double LHAup::px(int i)
double LHAup::py(int i)
double LHAup::pz(int i)
double LHAup::e(int i)
double LHAup::m(int i)
double LHAup::tau(int i)
double LHAup::spin(int i)
double LHAup::scale(int i)
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
void LHAup::setIdX(int id1, int id2, double x1, double x2)
This information is returned by the methods
int LHAup::id1()
int LHAup::id2()
double LHAup::x1()
double LHAup::x2()
the flavour and 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
void LHAup::setPdf( int id1pdf, int id2pdf, double x1pdf, double x2pdf, double scalePDF, double pdf1, double pdf2, bool pdfIsSet)
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
bool LHAup::pdfIsSet()
int LHAup::id1pdf()
int LHAup::id2pdf()
double LHAup::x1pdf()
double LHAup::x2pdf()
double LHAup::scalePDF()
double LHAup::pdf1()
double LHAup::pdf2()
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.
The maximum scale for parton-shower evolution of a Les Houches event is
regulated by the
TimeShower:pTmaxMatch
and
SpaceShower:pTmaxMatch
modes. If you want to guarantee that the input scale
value
is respected, as is often the case in matching/merging procedures, you
should set both of these modes to 1. That only affects the hard process,
while resonance decays are still processed using the resonance mass to
set the upper limit. However, the optional
Beams:strictLHEFscale = on
setting restricts also resonance-decay emissions to be below the input
scale
value.
As a further non-standard feature, it is also possible to read in the
separate scale values of all final particles. Such scale values could be used
e.g. to restrict the maximum scale for shower evolutions for each parton
separately. This reading will only be applied if the
Beams:setProductionScaleFromLHEF
switch is true (see
Beam Parameters
for details).
This information is returned by the method
double LHAup::scale(int i)
. When no such information
has been read from the LHEF, the scale defaults to -1.
void LHAup::listEvent()
prints the above information for the current event. In cases where the
LHAup
object is not available to the user, the
Pythia::LHAeventList()
method can
be used, which is a wrapper for the above.
virtual bool LHAup::skipEvent(int nSkip)
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.
The LHA standard only allows for one hard subcollision in an event.
Further multiparton interactions are supposed to be handled by the
internal MPI machinery. As a nonstandard feature, it is possible
to input two hard subcollisions in the same event, to match the internal
second hard process machinery.
In such cases two partons are extracted from each of the two incoming
hadrons. A restriction is that, unlike the single-subprocess case,
it is important that the partons are input in the order that PYTHIA
later would need. That is, the two subcollisions should follow each
other, with instate preceding outstate. Any resonance decay chain
should be put at the end, after both interactions. As illustration,
consider double W production. With 1 and 2
labelling the two subcollisions, and A and B the two
incoming hadron beams, the event record ordering should be
in_A1 - in_B1 - W_1 - in_A2 - in_B2 - W_2 - f_1 - fbar_1 - f_2 -
fbar_2, where f fbar is the fermion decay products of
the respective W. A limitation is that currently only one
input scale is used, that thereby limits all further partonic activity
in the same way for both processes.
When transferring events through the runtime interface, it is worth
noting that PYTHIA offers some settings to ensure the consistency of
momenta, e.g. reshuffling of particle momenta to guarantee on-shell
conditions, or matching the sum of incoming to the sum of outgoing momenta.
The documentation of settings related to this can be found under the
Transfer to the PYTHIA process record heading of
Les Houches Event Files .
A LHEF reader class
The LHEF standard ([Alw06], [But14]) 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, and is discussed in detail in the
Les Houches Event Files section.
Internally the file handling and reading is us handled by instance of the
derived class LHAupLHEF
.
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.
bool LHAup::setInitLHEF(ifstream& is, bool readHeaders = false)
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.
bool LHAup::setNewEventLHEF(ifstream& is)
read in event information from the specified stream into a staging area
where it can be reused by setOldEventLHEF
.
bool LHAup::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.
These three main methods build on a number of container classes and a
generic LHEF reader class (called Reader
) found in
LHEF3.h
and LHEF3.cc
. The Reader
handles all the parsing and storage necessary to adhere with
[But14]. (A matching Writer
class is also
available; see documentation in LHEF3.h
how it can be
used.) All parsing that is not strictly part of the LHEF format
(e.g. the reading of header information) is instead performed directly in
the LHAupLHEF
methods.
Some other small utility routines are:
bool LHAup::fileFound()
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.
bool LHAup::useExternal()
always returns false in the base class, but in LHAupLHEF
it returns false if the LHAupLHEF
instance is constructed to
work on an input LHE file, while it returns true if the LHAupLHEF
instance is constructed to use externally provided input streams instead.
For the latter, the LHAupLHEF
instance should have been
constructed with the class constructor
LHAupLHEF(Info* infoPtrIn, istream* isIn, istream* isHeadIn,
bool readHeadersIn, bool setScalesFromLHEFIn)
.
void LHAup::setInfoHeader(const string &key, const string &val)
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.
bool LHAup::newEventFile(const char* fileIn)
close current event input file/stream and open a new one, to
continue reading events of the same kind as before.
istream* LHAup::openFile(const char *fn, ifstream &ifs)
void LHAup::closeFile(istream *&is, ifstream &ifs)
open and close a file, also gzip files, where an intermediate
decompression layer is needed.
void LHAupLHEF::closeAllFiles()
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
include/Pythia8Plugins/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, and used to create a derived class that
contains 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
this support has been removed, but can still be recuperated from
previous code versions, in a reduced form up to version 8.176.
A LHEF writer class
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. The main documentation of these features is found in
the Les Houches Event Files section. Nevertheless,
we quickly summarize the main workhorse functions here. There are four
routines in the base class that should be called in sequence to build up the
proper file structure.
bool LHAup::openLHEF(string filename)
Opens a file with the filename indicated, and writes a header plus a brief
comment with date and time information.
bool LHAup::initLHEF()
Writes initialization information to the file above. Such information should
already have been set with the methods described in the "Initialization"
section above.
bool LHAup::eventLHEF(bool verbose = true)
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.
bool LHAup::closeLHEF(bool updateInit = false)
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.
string LHAup::getFileName()
Return the name of the LHE file above.