MLM jet merging

  1. Event input source
  2. Jet Matching parameters
  3. Alpgen-specific parameters
  4. Madgraph-specific parameters
  5. Alpgen-style parton-jet merging
  6. Madgraph-style parton-jet merging
  7. A note on combining UserHooks
This manual page describes the parton-jet matching interfaces for PYTHIA8. In this approach, usually referred to as MLM merging [Man02, Man07], the final jets after parton-shower evolution and jet clustering are matched to the original partons. The event is accepted if a reasonable match is found, and rejected if not. The rejection step in an approximate way introduces a Sudakov form factor on to the hard processes. Notably the parton shower should not generate an emission that would doublecount hard activity already included in the matrix-element description. Within this general ansatz, different technical solutions can be adopted. We provide two alternatives, one based on the algorithm used in ALPGEN [Man03], and another on the one used in Madgraph [Alw11], both reimplemented from scratch here. The main points of these two algorithms are outlined further down on this page.

We also allow for two alternative sources of external events, one in the ALPGEN native format and one in the Madgraph LHEF-based one. All four combinations of input format and scheme are implemented. In the following it is therefore important to keep the two aspects apart, whenever the ALPGEN and Madgraph labels are used.

Currently all the files of interest are located in the include/Pythia8Plugins/ subdirectory:

Event input source

External sources of partons are used in the parton-jet matching process. The source of the partons has been separated from the implementation of the algorithm. By default, PYTHIA8 contains a machinery to process Les Houches Event Files (LHEFs) Madgraph5 adheres to this format, but also contains some further non-standardized information that can be used. The parsing of the native ALPGEN file format is described on the Alpgen Event Interface page.

Commonly, the source of external partons also contains information about how a particular type of algorithm should be employed. This information is handled by the AlpgenPar class for ALPGEN files, and MadgraphPar for LHEFs. The user can choose to set default merging parameters using the Alpgen:setMLM flag for ALPGEN files. For LHEFs, instead, the setting of default parameters is controlled with the JetMatching:setMad flag:

flag  JetMatching:setMad   (default = on)
When enabled, the merging parameters are set according to the values in the LHEF header. Specifically, the header must set the ickkw, xqcut, maxjetflavor and alpsfact values, and ickkw must be nonzero. Note that these labels are Madgraph-specific. For other programs with LHEF output, or for Madgraph files lacking this information, these parameters should be set by the user (or one can rely on the default values). The following parameters (described below) must then be specified:

With this flag on, the values from the LHEF for these parameters take precedence over other values.

Jet Matching parameters

A class JetMatching, derived from UserHooks, is used to define the basic structure of a parton-jet matching algorithm. Two versions are implemented here, based on the FORTRAN code provided by the ALPGEN and Madgraph packages, respectively: JetMatchingAlpgen and JetMatchingMadgraph. The merging parameters are defined with the JetMatching:* keyword.

Scheme and Usage

flag  JetMatching:merge   (default = off)
Master switch to activate parton-jet matching. When off, all external events are accepted (unless they are rejected due to weighting or event processing problems).

mode  JetMatching:scheme   (default = 1; minimum = 1; maximum = 2)
The parton-jet MLM-style merging scheme.
option 1 : The one inspired by the Madgraph scheme, here implemented in the JetMatchingMadgraph class.
option 2 : The one inspired by the ALPGEN scheme, here implemented in the JetMatchingAlpgen class.

Jet algorithm

The choice of jet algorithm and associated parameters can be adjusted with the settings below. The PYTHIA8 internal CellJet and SlowJet routines are used for jet finding. See the Event Analysis page for more details.

mode  JetMatching:jetAlgorithm   (default = 1; minimum = 1; maximum = 2)
The choice of jet algorithm to use when merging against hard partons. Currently, only SlowJet with the kT algorithm (and useStandardR = false) is supported for Madgraph-style merging, while there is full freedom for the ALPGEN-style merging.
option 1 : The CellJet cone algorithm.
option 2 : The SlowJet clustering algorithm.

mode  JetMatching: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.

mode  JetMatching:nEta   (default = 100; minimum = 50)
Specific to the CellJet algorithm, the number of bins in pseudorapidity.

mode  JetMatching:nPhi   (default = 60; minimum = 30)
Specific to the CellJet algorithm, the number of bins in phi.

parm  JetMatching: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  JetMatching:eTthreshold   (default = 0.1)
Specific to the CellJet algorithm, cells with eT < eTthreshold are completely neglected by the jet 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. These values can be set automatically based on the information in the ALPGEN file or LHEF.

parm  JetMatching: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. For the Madgraph scheme, this parameter should match the qCut parameter described later.

parm  JetMatching: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  JetMatching: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. For the Alpgen scheme, particles within etaJetMax + coneRadius are passed to the jet algorithm, with only jets within etaJetMax retained in the merging. For the Madgraph scheme, only particles within etaJetMax are used.

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  JetMatching:exclusive   (default = 2; minimum = 0; maximum = 2)
Exclusive or inclusive merging.
option 0 : The merging is run in inclusive mode. All partons must match jets, but additional jets are allowed, provided they are not harder than the matched jets.
option 1 : The merging is run in exclusive mode. All partons must match jets, and no additional jets are allowed.
option 2 : If nJet < nJetMax, then the merging is run in exclusive mode, otherwise it is run in inclusive mode. For the Madgraph scheme, this is checked on an event-by-event basis, which is useful when an LHEF contains a "soup" of partonic multiplicities. If nJetMax < 0 or nJet < 0, then the algorithm defaults to exclusive mode.

mode  JetMatching:nJet   (default = -1; minimum = -1)
When JetMatching:exclusive = 2, nJet indicates the minimum number of additional light jets in the incoming process. This value may be set automatically.

mode  JetMatching:nJetMax   (default = -1; minimum = -1)
When JetMatching:exclusive = 2, nJetMax is used to indicate the maximum number of jets that will be matched.

Jet matching

The following parameters control the criteria for matching a clustered jet to a hard parton.

mode  JetMatching: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.

Alpgen-specific parameters

mode  JetMatching: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 coneMatchLight * coneRadius for light jets or coneMatchHeavy * coneRadiusHeavy for heavy jets. Note that by default coneRadiusHeavy = 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  JetMatching:coneMatchLight   (default = 1.5; minimum = 0.1)
The coneMatchLight parameter used when JetMatching:jetMatch = 1.

parm  JetMatching:coneRadiusHeavy   (default = -1.0)
The coneRadiusHeavy parameter used when JetMatching:jetMatch = 1. When assigned a negative value, the value of JetMatching:coneRadius is used.

parm  JetMatching:coneMatchHeavy   (default = 1.0; minimum = 0.1)
The coneMatchHeavy parameter used when JetMatching:jetMatch = 1.

Madgraph-specific parameters

flag  JetMatching:doShowerKt   (default = off)
This switch changes the merging prescription to the shower-kT scheme outlined in [Alw08]. This scheme differs from "classical" MLM jet merging with respect to when the matching veto is checked. The shower-kT scheme considers already immediately after the first shower emission if an event should be discarded.

parm  JetMatching:qCut   (default = 10.0; minimum = 0.0)
kT scale for merging shower products into jets.

mode  JetMatching:nQmatch   (default = 5; minimum = 3; maximum = 6)
Controls the treatment of heavy quarks.
option 5 : All quarks (except top) are treated as light quarks for matching.
option 4 : Bottom quarks are treated separately. Currently, they are unmatched.

parm  JetMatching:clFact   (default = 1.0)
The clFact parameter determines how jet-to parton matching is done. A match is defined as a squared cluster scale that equals:
|clFact| * qCut for inclusive mode,
|clFact| * max(qCut,min(pT(parton))) for exclusive mode, clFact ≥ 0, or
|clFact| * min(kT(parton)) for exclusive mode, clFact < 0.

flag  JetMatching:doVeto   (default = on)
If turned off, then no jet matching veto will be applied internally in Pythia. Instead, it is assumed that the (expert) user enforces all necessary vetoes externally by some other means. Do not change the default value unless you are an expert in MLM jet matching and want to use your own code to perform the necessary vetoes.

An implementation of the FxFx prescription for combining multiple NLO calculations [Fre12] and its updated treatment of non-enhanced jets as described in [Fre21] is available. FxFx merging with aMC@NLO shares most parameters with the leading-order (MadGraph-style) MLM prescriptions and can be activated by using the three additional settings below.

flag  JetMatching:doFxFx   (default = off)
If turned on, then FxFx merging with aMC@NLO inputs is performed. Note that this requires event samples that are specifically generated for this task.

mode  JetMatching:nPartonsNow   (default = -1)
The number of partons in Born-like events for the current input LHEF. If the current sample e.g. contains pp → e+e- + 2 partons Born-like configurations, and pp → e+e- + 3 partons Real-emission-type events, then JetMatching:nPartonsNow = 2 applies.

parm  JetMatching:qCutME   (default = 10.0)
The cut applied to regulate multi-jet matrix elements. Note that this cut can differ from the matching scale.

Alpgen-style parton-jet merging

This section describes the Alpgen-style MLM merging algorithm for PYTHIA8. 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.

First, either the CellJet or SlowJet jet algorithm is chosen. Both of these algorithms have an R and an 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:

Some different options are provided, specified further above in the parameters section. 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 User Hooks machinery, and is therefore already taken into account in the cross section. In the output from Pythia::stat(), 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 JetMatching:nJetMax parameter.

Madgraph-style parton-jet merging

This section describes the Madgraph parton-jet merging algorithm for PYTHIA8.

First, the kT jet algorithm is applied using the PYTHIA8 SlowJet implementation. The useStandardR = false is used, ie. the (delta R)^2 separation is defined as 2 (cosh(delta y) - cos(delta phi)) rather than the more common (delta y)^2 + delta phi)^2. The R, etaMax, and a pTmin parameters are specified. By default, R = 1 and pTmin = qCut . It is not recommended to change these. These should match the algorithm parameters used in the Madgraph Matrix Element (ME) generation.

ME samples with different jet multiplicities are run through the event generator, and the generation is interrupted after parton showers have been applied, but before resonance decays and beam remnants have been processed. In what follows, the hardness measure of jets/partons is taken to be kT relative to qCut. The hard system (ignoring all MPI systems) is analyzed:

In exclusive mode, it is expected that ME samples with higher parton multiplicity are available to fill the phase space above qCut. The inclusive mode is when there are no such samples, and the parton shower is used to fill the phase space.

Some different options are provided, specified further above. These are set so that, by default, the algorithm closely follows the FORTRAN interface ME2Pythia provided in the Madgraph package.

All vetoing of events is done through the usual User Hooks machinery, and is therefore already taken into account in the cross section. In the output from Pythia::stat(), 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 the hard process parameters are used to guide the merging procedure, events will be matched in the exclusive mode.

Madgraph scheme with no internal vetoes (assuming an external veto implementation)

This section describes the facilities that allow expert users to use their own veto code to perform a Madgraph-style merging. This can e.g. be useful to assess parameter changes without having to process the same input events multiple times.

As a first step, any vetoes in the Pythia Jet Matching need to be disabled by using JetMatching:doVeto = off. In this mode, Pythia only stores all the information that is necessary to check (and apply) the shower-kT or kT-MLM vetoes externally by hand. This information can be accessed by calling the functions

Event JetMatchingMadgraph::getWorkEventJet()  
Return the event after parton showering, without resonance decay products and particles that should not be included in the jet matching, as necessary to implement the vetoes in the kT-MLM scheme.

Event JetMatchingMadgraph::getProcessSubset()  
Return the event record containing the hard process, without resonance decay products and particles that should not be included in the jet matching, as necessary to implement the vetoes in the shower-kT and kT-MLM schemes. In the former, this event is needed to find the lowest pT in the ME inputs. In the latter, the event record is used to count the number of hard-process partons, minimal hard process pT, and to perform the matching of hard-process particles to shower jets.

bool JetMatchingMadgraph::getExclusive()  
Return flag to identify if exclusive or inclusive vetoes should be applied to this event.

double JetMatchingMadgraph::getPTfirst()  
Return the transverse momentum (w.r.t. the beam) of the first parton shower emission, as needed for the shower-kT scheme.

vector <double> JetMatchingMadgraph::getDJR()  
Return a vector of jet clustering scales produced by running the jet algorithm used for jet matching on the event record without resonance decay products and particles that should not be included in the matching. In this vector, clustering scales for combining few jets appear before scales from combining many jets. This function is useful for the kT-MLM scheme, or to have quick access to this information for histogramming and sanity checks.

We do not currently supply example code for this very advanced functionality. Interested expert users should feel free to contact the Pythia authors for further explanations.

A note on combining UserHooks

As have been noted above, the matching is implemented using classes derived from the UserHooks class, thereby gaining access to the event generation process at the relevant locations. For native ALPGEN files, which do not adhere to the Les Houches standards, it is also necessary to intervene with a UserHooks-derived AlpgenHooks to handle the extraction and setting of relevant extra information.

One must then combine multiple UserHooks classes, 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 JetMatching : virtual public UserHooks 
Without this option, when combining two UserHooks-derived classes, two copies of the base UserHooks class would be created, leading to ambiguities.

The two first classes in CombineMatchingInput.h combine ALPGEN input with the two different matching schemes, e.g. for the first

 
class JetMatchingAlpgenInputAlpgen : public AlpgenHooks, 
  public JetMatchingAlpgen { 
public: 
  // Constructor and destructor. 
  JetMatchingAlpgenInputAlpgen(Pythia& pythia) : AlpgenHooks(pythia), 
    JetMatchingAlpgen() { } 
  ~JetMatchingAlpgenInputAlpgen() {} 
  // Initialisation. 
  virtual bool initAfterBeams() { 
    if (!AlpgenHooks::initAfterBeams()) return false; 
    if (!JetMatchingAlpgen::initAfterBeams()) return false; 
    return true; 
  } 
  // Process level vetos. 
  virtual bool canVetoProcessLevel() { 
    return JetMatchingAlpgen::canVetoProcessLevel(); 
  } 
  .... 
}; 
This class inherits from both AlpgenHooks and JetMatchingAlpgen. 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().