ALPGEN and MLM merging
This manual page describes the ALPGEN [Man03] and MLM merging
[Man02] interfaces for PYTHIA 8. While future versions of
ALPGEN will be able to write out events in LHEF format, previous
versions always output events in an ALPGEN native format (a combination
of a ".unw" and a "_unw.par" file). The ALPGEN component of this code
contains a reader for this native format (for unweighted events), as
well as parameter reading for both ALPGEN native and LHE file formats.
Although designed to work together with the MLM component, they are
implemented entirely separately and it is possible to use one without
the other.
It should be noted that all the functionality described here is provided
through external routines, and therefore the presence of these features is
dependent on the main program being used. This structure allows for the
easy extensibility of the merging scheme. The files of interest are located
in the examples/
subdirectory:
-
MLMhooks.h
: provides a
UserHooks
derived
class for performing the MLM merging procedure. All MLM:*
options, described below, are implemented in this file. Further
technical details of this class are given at the end of this manual
page.
-
LHAupAlpgen.h
: provides three classes for the reading of
ALPGEN event and parameter files. LHAupAlpgen
is an
LHAup
derived
class for reading in ALPGEN native format event files.
AlpgenPar
is a class for the parsing of ALPGEN parameter
files, making the information available through a simple interface.
AlpgenHooks
is a
UserHooks
derived class that
provides the Alpgen:*
options, described below. Further
technical details of these classes are given at the end of this manual
page.
-
main32.cc, main32.cmnd
: a sample main program and card
file showing the usage of the previous two files. It combines the Alpgen
and MLM UserHooks classes together, such that the functionality of both
is available, and reads in a sample ALPGEN event file while performing
the MLM merging procedure. Some commented out sets of options are
provided in the card file, which can be activated to try different
merging setups.
-
main32.unw, main32_unw.par
: an ALPGEN format event and
parameter file containing 100 W + 3 jet events. It is not feasible
to package large event files with the PYTHIA distribution, but this
sample is enough to show the different components in action.
ALPGEN
NB: these options are provided by the AlpgenHooks class,
which must be loaded for this functionality to be present
ALPGEN event files that have been written out in LHEF format should be
read in through the normal LHEF machinery
(see beam parameters). Files in
ALPGEN's native format, instead, may be processed using the
Alpgen:file
option below. When using this option, the
ALPGEN parameter file is stored in the PYTHIA Info object under the key
AlpgenPar
, see the "Header information" section of the
Event Information manual page for
more details. Processes not implemented by the PYTHIA 6 interface
supplied with ALPGEN are also not implemented here.
When reading in ALPGEN native event files, some momenta are shifted by
the file reader to ensure energy-momentum conservation. The magnitude of
these shifts should be small (around the MeV level in the worst case)
and warnings will be produced if they are above a set threshold. A large
number of warnings may signify unexpected behaviour and should
potentially be investigated. It is also known that certain event
classes, for example an event with both light and heavy b
quarks may give rise to these warnings.
The ALPGEN file reader supports the reading of the event and parameter
files in gzip format with file extensions ".unw.gz" and "_unw.par.gz"
respectively. This requires the use of external libraries, however, and
the README
file in the main directory contains instructions
on how to enable this.
All other Alpgen:*
options apply to both LHE and native
file formats, and include options to guide the MLM merging procedure
based on the parameters that are read in with the events file.
word
Alpgen:file
(default = void
)
This option is used to read in ALPGEN format event files. Using this option
overrides any previously set beam options inside PYTHIA. The path to the
files, not including any file extension, should be provided e.g. for input
files input_unw.par and input.unw, the value
input should be used.
flag
Alpgen:setMasses
(default = on
)
When switched on, any particle masses provided by ALPGEN are set in
the PYTHIA particle database.
flag
Alpgen:setMLM
(default = on
)
When switched on, the merging parameters (see below) are set according to
the ALPGEN hard process cuts:
- MLM:eTjetMin = min(ptjmin + 5., 1.2 * ptjmin),
- MLM:coneRadius = drjmin,
- MLM:etaJetMax = etajmax,
where the ptjmin
, drjmin
and
etajmax
are the incoming ALPGEN parameters. Note that any
existing values of these parameters are overwritten.
flag
Alpgen:setNjet
(default = on
)
When switched on, the MLM:nJet
parameter (see below) is set
to the incoming njet
ALPGEN parameter. Note that any
existing value of this parameter is overwritten.
MLM Merging
NB: these options are provided by the MLMhooks class,
which must be loaded for this functionality to be present
This section describes the MLM merging interface for PYTHIA 8. The most
common reference to the algorithm is [Man02]. In many respects,
however, the implementation provided in the ALPGEN package should be
considered the official description of the MLM merging procedure.
Although designed primarily to work with events generated with ALPGEN, it
can in principle also be used with events from a different source. This
should not be done without thought, however, and it is up to the user to
understand the details of the algorithm and the implications of using a
different hard process generator.
A brief description of the MLM procedure, as implemented, is given here.
First, either the CellJet or SlowJet jet algorithm is chosen, see the
Event Analysis page for more details.
Both of these algorithms have an R and etaMax
parameter. In addition, CellJet has an eTmin and SlowJet has a
pTmin parameter. These are the primary three parameters of the
merging procedure, and in practice are set dependent on the cuts applied
to the matrix element (ME) generation. We stress that the merging
procedure is not tied to the geometry of a specific physical detector,
but only to the match between the original partons and the resulting
jets, using standard jet algorithms in the phase space region where
partons have been generated.
ME samples with different jet multiplicities are run through the event
generator, and the generation interrupted after parton showers have been
applied, but before resonance decays and beam remnants have been
processed. Note in particular that top quarks will not yet
be decayed, which may lead to slight differences with the PYTHIA 6
interface included with the ALPGEN package. In what follows, the
hardness measure of jets/partons is taken to be eT when CellJet
is used and pT when SlowJet is used. The hard system (ignoring
all MPI systems) is then analysed:
-
The particles in the original matrix element process are sorted into
light partons, heavy partons and other particles. For backwards
compatibility, a light parton is defined as the set (d, u, s, c,
b, g) with zero mass. A heavy parton is defined as the set
(c, b, t) with non-zero mass.
-
All particles not originating from the heavy partons or other
particles are passed to the jet algorithm and clustered.
-
Clustered jets are matched to the light partons in the original ME
process. There are two different methods which can be used:
-
Method 1: The following is done for each parton, in order
of decreasing hardness. The delta R between the parton
and all jets is calculated and the smallest value taken. If
this is less than the jet R parameter, possibly
multiplied by a constant, the jet and parton are considered to
match, and the jet is removed from further consideration.
Note that for CellJet the delta R measure is in
(eta, phi), while for SlowJet, it is in
(y, phi).
-
Method 2: This method is only possible when using the SlowJet
algorithm. Before the clustering is performed, extremely soft
"ghost" particles are added to the event at the
(y, phi) coordinates of the original matrix element
partons. If such a particle is clustered into a jet, the parton
and jet are considered to match. The idea of "ghost" particles
was originally introduced by FastJet as a way to measure jet
areas [Cac06] and should not affect clustering with an
infrared-safe jet algorithm.
-
If there is a light ME parton remaining which has not been matched
to a jet, then the event is vetoed. If all ME partons have been
matched to a jet, but there are still some extra jets remaining,
then two options are possible:
-
Exclusive mode: the event is vetoed. This is typically used when
there are ME samples with higher jet multiplicities, which would
fill in the extra jets.
-
Inclusive mode: the event is retained if the extra jets are softer
than the softest matched jet. This is typically used when
there is no ME sample with higher jet multiplicity, so the parton
shower should be allowed to give extra jets.
-
All particles originating from the heavy partons are passed to the
jet algorithm and clustered.
-
The clustered jets are again matched to the original partons, but
there is no requirement for a match to be present; all matched jets
are immediately discarded. The matching procedure is much the same
as for light partons, but with two differences when delta R
matching is used. First, a different R parameter than that
used by the jet algorithm may optionally be given. Second, all jets
that are within the given radius of the parton are matched, not
just the one with the smallest delta R measure. If there
are still extra jets remaining then in exclusive mode the event is
immediately vetoed, while in inclusive mode the event is retained if
the extra jets are softer than the softest light matched jet.
Some different options are provided, specified further below. These
are set so that, by default, the algorithm closely follows the official
MLM interface provided in the ALPGEN package. All vetoing of events
is done through the usual
UserHooks
machinery, and is
therefore already taken into account in the cross section.
In the output from
Pythia::statistics()
,
the difference between the "Selected" and "Accepted" columns gives the
number of events that have not survived the vetoing procedure. It is
still the responsibility of the user to add together the results from
runs with different jet multiplicities. In the simplest case, when
ALPGEN input is used and the hard process parameters are used to guide
the merging procedure, it is enough to set the MLM:nJetMax
parameter (see below).
Merging
flag
MLM:merge
(default = off
)
Master switch to activate MLM merging.
Exclusive mode
The following settings determine whether clustered jets which do not
match an original hard parton are allowed. They are typically permitted
in the highest jet multiplicity sample, where the parton shower may
produce extra hard jets, without risk of double counting. Any
extra jet produced by the shower must be softer than any matched light
jet, or else the event is vetoed.
mode
MLM:exclusive
(default = 2
; minimum = 0
; maximum = 2
)
Exclusive or inclusive merging.
option
0 :
The merging is run in inclusive mode.
option
1 :
The merging is run in exclusive mode.
option
2 :
If MLM:nJet < MLM:nJetMax, then the merging is run in
exclusive mode, otherwise it is run in inclusive mode. If
MLM:nJetMax < 0 or MLM:nJet < 0, then the
algorithm defaults to exclusive mode.
mode
MLM:nJet
(default = -1
; minimum = -1
)
The number of additional light jets in the incoming process.
This is used when MLM:exclusive = 2.
Note that for ALPGEN input, this value may be automatically set.
mode
MLM:nJetMax
(default = -1
; minimum = -1
)
This parameter is used to indicate the maximum number of jets that will be
matched. This is used when MLM:exclusive = 2.
Jet algorithm
The choice of jet algorithm and associated parameters can be adjusted with
the settings below. The PYTHIA 8 internal CellJet
and
SlowJet
routines are used, see the
Event Analysis page for more details.
mode
MLM:jetAlgorithm
(default = 1
; minimum = 1
; maximum = 2
)
The choice of jet algorithm to use when merging against hard partons.
option
1 : The CellJet algorithm.
option
2 : The SlowJet algorithm.
mode
MLM:nEta
(default = 100
; minimum = 50
)
Specific to the CellJet algorithm, the number of bins in pseudorapidity.
mode
MLM:nPhi
(default = 60
; minimum = 30
)
Specific to the CellJet algorithm, the number of bins in phi.
parm
MLM:eTseed
(default = 1.5
; minimum = 0.0
)
Specific to the CellJet algorithm, the minimum eT for a cell
to be acceptable as the trial center of a jet.
parm
MLM:eTthreshold
(default = 0.1
)
Specific to the CellJet algorithm, cells with eT < eTthreshold
are completely neglected by the jet algorithm.
mode
MLM:slowJetPower
(default = -1
; minimum = -1
; maximum = 1
)
The power to use in the SlowJet algorithm.
option
-1 : The anti-kT algorithm.
option
0 : The Cambridge/Aachen algorithm.
option
1 : The kT algorithm.
Merging parameters
The following options are the three main parameters for the merging
procedure. Although here they are in principle free parameters, they should
be heavily influenced by the hard process generation cuts. For ALPGEN
input, these values may be set automatically.
parm
MLM:eTjetMin
(default = 20.0
; minimum = 5.0
)
For the CellJet algorithm, this gives the minimum transverse energy
inside a cone for a jet to be accepted. For the SlowJet algorithm, this
is instead the minimum transverse momentum required for a cluster to be
accepted as a jet.
parm
MLM:coneRadius
(default = 0.7
; minimum = 0.1
)
For the CellJet algorithm, this gives the size of the cone in (eta, phi)
space drawn around the geometric center of the jet. For the SlowJet
algorithm, this gives the R parameter.
parm
MLM:etaJetMax
(default = 2.5
; minimum = 0.1
)
For both jet algorithms, this defines the maximum pseudorapidity that
the detector is assumed to cover. In this context, however, it is tied
to the phase space region in which partons have been generated.
Specifically, particles within etaJetMax + coneRadius are
passed to the jet algorithm, with only jets within etaJetMax
retained in the merging.
Jet matching
The following parameters control the criteria for matching a clustered jet
to a hard parton.
mode
MLM:jetAllow
(default = 1
; minimum = 1
; maximum = 2
)
Controls which particles are clustered by the jet algorithm.
option
1 :
This option explicitly disallows top quarks, leptons and photons. All other
particle types are passed to the jet algorithm.
option
2 :
No extra particles are disallowed.
mode
MLM:jetMatch
(default = 1
; minimum = 1
; maximum = 2
)
Criteria for matching a clustered jet to a parton.
option
1 :
This option can be used with both the CellJet and SlowJet
algorithms. The delta R between each parton and
jet is taken, and the minimal value compared against
MLM:coneMatchLight * MLM:coneRadius for light jets or
MLM:coneMatchHeavy * MLM:coneRadiusHeavy for heavy jets.
Note that by default MLM:coneRadiusHeavy = MLM:coneRadius,
see below. If below this value, the parton and jet are considered
to match. With CellJet, the delta R measure is in
(eta, phi), while with SlowJet it is in (y, phi).
option
2 :
This option can only be used with the SlowJet algorithm. The hard
partons are inserted into the parton level event as "ghost" particles,
but at the correct (y, phi) position. If this particle is then
clustered into a jet, it is considered a match.
parm
MLM:coneMatchLight
(default = 1.5
; minimum = 0.1
)
The coneMatchLight
parameter used when
MLM:jetMatch = 1.
parm
MLM:coneRadiusHeavy
(default = -1.0
)
The coneRadiusHeavy
parameter used when
MLM:jetMatch = 1. When assigned a negative value,
the value of MLM:coneRadius
is used.
parm
MLM:coneMatchHeavy
(default = 1.0
; minimum = 0.1
)
The coneMatchHeavy
parameter used when
MLM:jetMatch = 1.
Class information
Some more technical information about the different classes is given
below. For clarity, some limited information on certain private methods
is provided.
LHAupAlpgen
This class is derived from the
LHAup
base class, and
uses the standard machinery to pass initialisation and event data to
PYTHIA. These standard functions are not documented here. The complete
parameter file is stored in the PYTHIA Info object, if given, under the
key AlpgenPar
.
LHAupAlpgen::LHAupAlpgen(const char *baseFNin, Info *infoPtrIn = NULL)
The constructor for the class takes the base filename for the ALPGEN
format files (without file extensions) and optionally a pointer to a
PYTHIA Info class, used for warning/error message printing and for
storing the ALPGEN parameter file. The event and
parameter files are opened immediately, with the AlpgenPar
class, described below, used to parse the parameter file.
bool LHAupAlpgen::addResonances()
This is a private method used when an event is read in. The information
read from the event file does not always contain a complete listing of
all particles and four-momenta, and so various details must be
reconstructed. Exactly which details are filled in can vary based on the
ALPGEN process in question.
bool LHAupAlpgen::rescaleMomenta()
This is another private method used when an event is read in.
It shuffles and rescales momenta in an event to ensure energy-momentum
conservation. First, pT is made to balance by splitting any
imbalance between all outgoing particles with their energies also
scaled. Second, the e/pZ of the two incoming particles are
scaled to balance the outgoing particles. Finally, any intermediate
resonances are recalculated from their decay products.
AlpgenPar
This class parses an ALPGEN parameter file and makes the information
available through a simple interface. The information is stored
internally in key/value (string/double) format. All lines prior to:
************** run parameters
are ignored, and in the general case, a line e.g.
10 3.00000000000000 ! njets
would be stored with key "njets" and value "3.0". The following lines
are special cases where the line may be split or the key translated:
3 ! hard process code
0.000 4.700 174.300 80.419 91.188 120.000 ! mc,mb,mt,mw,mz,mh
912.905 0.0914176 ! Crosssection +- error (pb)
100 29787.4 ! unwtd events, lum (pb-1) Njob= 2
In the first line, the key "hard process code" is translated to
"hpc". In the second, the mass values are split and each given an entry
in the internal store. In the third, the cross section and cross section
error are stored under the keys "xsecup" and "xerrup" respectively.
Finally, the number of events and luminosity are stored under the keys
"nevent" and "lum" respectively. In the event that a duplicate key is
present, with differing values, the stored value is overwritten and a
warning given.
AlpgenPar::AlpgenPar(Info *infoPtrIn = NULL)
The constructor does nothing except for store the PYTHIA Info
pointer, if given. This is used for warning/error message printing.
bool AlpgenPar::parse(const string paramStr)
This method parses an ALPGEN parameter file. The parameter file is
passed as a single string, mainly intended to be read out from the
PYTHIA Info object using the header information methods.
bool AlpgenPar::haveParam(const string ¶mIn)
Method to check if a parameter with key paramIn
is present.
Returns true if present, else false.
double AlpgenPar::getParam(const string ¶mIn)
int AlpgenPar::getParamAsInt(const string ¶mIn)
Return the parameter with key paramIn
as a double or
integer. The presence of a parameter should have already been checked
using the haveParam()
function above. If the parameter is
not present, 0 is returned.
void AlpgenPar::void printParams()
Method to print a list of stored parameters.
AlpgenHooks
This UserHooks
derived class
provides all the Alpgen:*
options. It is provided as a
UserHooks class such that the code works regardless of whether ALPGEN
native or LHE file formats are used. It is declared with virtual
inheritance so that it may be combine with other UserHooks classes, see
the "Combining UserHooks" section below.
AlpgenHooks(Pythia &pythia)
The constructor takes a PYTHIA object as input, so that the beam
parameter settings can be overridden if the Alpgen:file
option is given. If this is the case, an LHAupAlpgen
instance is automatically created and passed to PYTHIA.
bool initAfterBeams()
This is the only UserHooks method that is overridden. It is called
directly after PYTHIA has initialised the beams, and therefore the
header information should be present in the PYTHIA Info object. The
AlpgenPar
class is used to parse ALPGEN parameters, if
present, which are then used to set further PYTHIA settings.
MLMhooks
This UserHooks
derived class
provides all the MLM:*
options. It is also declared with
virtual inheritance for combining with other UserHooks classes. It
uses standard UserHooks methods to perform the merging procedure
outlined previously in this manual page, and therefore detailed method
information is not given.
Combining UserHooks
It is possible to combine multiple UserHooks classes together, such that
the functionality of both is present. A prerequisite is that the
different UserHooks classes should be declared with virtual inheritance,
e.g.
class MLMhooks : virtual public UserHooks
Without this option, when combining two UserHooks derived classes
together, two copies of the base UserHooks class would be created
leading to ambiguity.
An example combining the AlpgenHooks and MLMhooks classes together is
given in main32.cc
:
class AlpgenAndMLMhooks : public AlpgenHooks, public MLMhooks {
public:
AlpgenAndMLMhooks(Pythia &pythia)
: AlpgenHooks(pythia), MLMhooks() {}
virtual bool initAfterBeams() {
// Call init of both AlpgenHooks and MLMhooks (order important)
if (!AlpgenHooks::initAfterBeams()) return false;
if (!MLMhooks::initAfterBeams()) return false;
return true;
}
};
This class inherits from both AlpgenHooks and MLMhooks. Any functions
which are present in both classes should be overridden with a function
that calls the different parent methods in the desired order. In the
above example, the only shared methods are the constructor and
initAfterBeams()
.