CKKW-L Merging

  1. Merging scale definitions
  2. Histogramming the events
  3. Merging with user-defined merging scale function
  4. Matrix element merging and HepMC output for RIVET
  5. Electroweak Merging
  6. Further variables
CKKW-L merging [Lon01] allows for a consistent combination of tree-level matrix elements containing multiple well-separated partons with each other and with parton showers. The result is a calculation that contains a mix of processes with different number of well-separated jets with fixed-order accuracy, improved by all-order resummation. The algorithm implemented in PYTHIA is described in [Lon11]. To perform matrix element merging, the user has to supply LHE files [Alw07] for the hard process and the corresponding process with up to N additional jets. This mix of processes is then internally disentangled to ensure that the inclusive fixed-order inputs can be converted to exclusive cross sections that no longer overlap. Please note that subtleties (and setting scheme) for the EW-improved way of disentangling processes presented in [Chr15a] is discussed in the section Electroweak Merging below.

The usage of the merging procedure is illustrated in a few example main programs (main81.cc, main82.cc, main83.cc, main84.cc and main85.cc, together with the input files main81.cmnd, main82.cmnd, main84.cmnd and main85.cmnd). These examples should of course only serve as an illustration, and as such will not make use of the merging in all possible ways. For full generality, the example programs link to LHAPDF, FastJet and HepMC. Of course the user is welcome to remove these dependencies. To remove the FastJet dependence, the functions calculating example observables have to be deleted. Removing the LHAPDF dependence requires changing the cmnd input files to choose an inbuilt PDF, as outlined in the PDF documentation. The HepMC dependence can be removed by erasing the code allowing for HepMC output.

The attached PYTHIA 8.3 Worksheet contains a Section 6 that walks you through the basics of setting up a merging in PYTHIA.

Three very short LHE files (w+_production_lhc_0.lhe, w+_production_lhc_1.lhe, w+_production_lhc_2.lhe) are included in the distribution. These files are not intended for physics studies, but only serve as input for the example main programs. For realistic studies, the user has to supply LHE files.

In the generation of LHE files, the value of the factorisation scale used in the PDFs is not important, since the cross section will be multiplied by ratios of PDFs to adjust to the PYTHIA starting scales. The same is true for the renormalisation scale (and starting value αs(MZ)) used to evaluate αs. Coupling and scale choices by the user will be transferred to the merging routines.

Multi-jet events can suffer from infrared divergences in the calculation. Sensible matrix element generator (MEG) outputs should not contain phase space points in which the calculation is ill-defined, meaning infrared regions need to be removed by cuts. This is most conveniently done by disallowing the MEG to produce partons below a minimal parton-parton separation in a certain jet algorithm. Using dedicated cuts to regularise MEG output is of course possible as well. Any regularisation criterion defines the matrix element region: The parts of phase space in which the fixed order calculation is considered valid and preferable to the parton shower. Matrix element merging is combining MEG events in the matrix element region with parton shower events in regions outside the regularisation cut (often called parton shower region). Because the regularisation cut defines a boundary between the matrix element and parton shower regions, i.e. the regions to be merged into one inclusive sample, it is usually called merging scale . Since many different cut choices may regularise the MEG calculation, many different merging scale definitions are possible. A few standard choices are listed below, as well as documentation on how to use a user-defined cut criterion. In combining matrix element and parton shower regions, the CKKW-L prescription tries to minimise the dependence on the merging scale. This can only be achieved if the combination of MEG events and parton shower populates the whole phase space. Additional cuts on the partons in the LHEF generation should hence be avoided as much as possible, meaning that the merging scale cut should always pose a more stringent cut than all other cuts on the partons. Of course, if the hard process itself is divergent, cuts need to be made. However, this should be chosen in such a way as to not exclude regions that will be available to the matrix elements with additional jets. An example is QCD di-jet production with additional jets: Say the 2 → 2 process is regularised with a pTmin cut of pTminCut = 100 GeV, and the 2 - >3 sample is regularised with a kTmin-cut of kTminCut = 50 GeV. This would mean that when reclustering the emission in the 2 → 3 sample, we could end up with a pT value pTminNow of the 2 → 2 configuration with pTminCut > pTminNow, which is excluded in the 2 → 2 sample. Thus, the 2 → 3 sample will include a Sudakov factor not included in the 2 → 2 sample, resulting in merging scale dependencies. Such dependencies can be avoided if the additional cuts on the hard process are minimal.

Of course, additional cuts on electroweak particles are allowed. These should be the same for all samples with 0 <= n <= N additional partons.

If it is not possible to generate LHE files with minimal cuts, the user can choose to use the MergingHooks structures in order to decide how much influence to attribute to parton shower histories in which the reclustered lowest multiplicity process does not pass the matrix element cuts. This is described below.

When generating LHE files, we advise against putting unstable particles (e.g. massive gauge bosons) in the final state. Rather, specify a resonance by its decay products, e.g. if Les Houches events for the pp → Z + jets → e+e- + jets process are desired, generate the matrix element events with the Z decay included. From a physical point of view, on-shell final massive gauge bosons should not be considered part of a hard process, since only the boson decay products will be detectable. Furthermore, non-narrow-width approximation contributions are not present if the ME generator only produces on-shell bosons. Interference effects between different production channels for the decay products would also be neglected. These points seem an unnecessary restriction on the accuracy of the ME calculation. In addition, there is a technical reason for this strategy. Since some matrix element generators choose to put additional information on intermediate bosons into Les Houches events, depending on if they pass a certain criterion (e.g. being close to the mass shell), without exact knowledge of this criterion, the only feasible way of bookkeeping the hard process is by identifying outgoing decay products.

Despite these considerations, (massive) gauge bosons in the final state are allowed in the hard process definition. This is useful particularly for Higgs physics, when different decays of the Higgs boson need to be simulated after the LHEF generation.

For all merging purposes, processes with different charge of outgoing leptons are considered different processes. That means e.g. that e+νe+ jets and e-ν̄e + jets are considered independent processes. If the user wishes to generate distributions including effects of more than one process, merged samples for all independent processes should be generated separately and added afterwards. Alternatively, to combine simple processes, combined LHE files can be used in conjunction with flexible containers (see below).

When the matrix element merging is used to produce HepMC [Dob01] files to be analysed with RIVET [Buc10], special care needs to taken in how the cross section is read by RIVET (see below).

To specify the merging conditions, additionally information on the merging scale value and the functional definition of the merging scale is needed. A few standard definitions of merging scales are available. We hope this makes the user interface less cumbersome.

Different choices intrinsic to the CKKW-L merging procedure might be relevant for the user as well. The base class MergingHooks gives the user the opportunity to define the functional form of the merging scale. In the following, the usage of the merging machinery to consistently include LHE files with additional jets into PYTHIA will be discussed.


Merging scale definitions

The quickest way to include processes with additional jets is to produce LHE files with one of the standard ways to define the merging scale. Three standard ways to define a merging scale (minimal kT, minimal evolution pT and by three cuts) are implemented. All of these prescriptions are equivalent - different definitions have only been introduced for the convenience of users, who might be limited by which cuts can be used in the generation of LHE files. Below, we describe how to switch on and use these different merging scale definitions.

Merging with merging scale defined in kT:

flag  Merging:doKTMerging   (default = off)
If the additional jets in the LHE files have been regulated by a kT cut, the user can supply the merging scale definition by setting this flag to on. kT here and below means cutting on Durham kT for e+e- collisions, and cutting on longitudinally invariant kT for hadronic collisions. Please note that this particular merging scale definition will check kT between all pairs of u,d,c,s,b,g partons.

Currently, the name longitudinally invariant kT is used for a few jet recombination algorithms with slightly different jet measures. A specific form can be chosen by setting the switch

mode  Merging:ktType   (default = 1; minimum = 1; maximum = 3)
Precise functional definition of longitudinally invariant kT. For e+e- collisions, Durham kT is always defined by the square root of min{ 2*min[ Ei2, Ej2] * [ 1 - cosθij] }, so that this switch will have no effect.
option 1 : Longitudinally invariant kT is defined by the square root of the minimum of minimal jet kinematic pT (pTkin,min2 = min{ pT,i2 } ) and pTlon,min2 = min{ min[ pT,i2, pT,j2] * [ (Δyij)2 + (Δφij)2 ] / D2 } , i.e. kT = min{ √pTkin,min2, √pTlon,min2 } for hadronic collisions. Note that the true rapidity of partons is used.
option 2 : Longitudinally invariant kT is defined by the square root of the minimum of minimal jet kinematic pT (pTkin,min2 = min{ pT,i2 } ) and pTlon,min2 = min{ min[ pT,i2, pT,j2] * [ (Δηij)2 + (Δφij)2 ] / D2 }, i.e. kT = min{ √pTkin,min2, √pTlon,min2 } for hadronic collisions. Note that the pseudorapidity of partons is used.
option 3 : Longitudinally invariant kT is defined by the square root of the minimum of minimal jet kinematic pT (pTkin,min2 = min{ pT,i2 } ) and pTlon,min2 = min{ min[ pT,i2, pT,j2] * [ cosh(Δηij) - cos(Δφij) ] / D2 } , i.e. kT = min{ √pTkin,min2 , √pTlon,min2 } for hadronic collisions.

parm  Merging:Dparameter   (default = 1.0)
The value of the D parameter needed in the definition of longitudinally invariant kT separation.

mode  Merging:nJetMax   (default = 0; minimum = 0)
Maximal number of additional jets in the matrix element. Note that the EW-improved "merging of mergings" strategy presented in [Chr15a] requires a different meaning of "additional", as explained in the "Electroweak Merging" section below.

parm  Merging:TMS   (default = 0.0)
The value of the merging scale. The name is inspired by the scale in evolution equations, which is often called 't', and the suffix 'MS' stands for merging scale. In the particular case of kT-merging, this would be the value of the kT-cut in GeV. For any merging scale definition, this input is considered the actual value of the merging scale.

word  Merging:Process   (default = void)
The string specifying the hard core process, in MG4/ME notation.

If e.g. W + jets merging should be performed, set this to pp>e+ve (without white spaces or quotation marks). This string may contain resonances in the MG/ME notation, e.g. for merging pp→Z W+→q q̄ e+νe + jets, the string pp>(z>jj)(w+>e+ve) would be applicable.

A lot more flexible hard process definitions are possible. To not dwell too much on these details here, we will come back to the process string at the end of this section.

flag  Merging:doMGMerging   (default = off)
Even easier, but highly non-general, is to perform the merging with MadGraph/MadEvent-produced LHE files, with a merging scale defined by a kT cut. For this, set this switch to on. The merging scale value will be read from the +1 jet LHE file by searching for the string ktdurham, and extracting the value from value = ktdurham. Also, the hard process will be read from the +0 jet LHE file, from the line containing the string @1 (the tag specifying the first process in the MadGraph process card). For this to work, PYTHIA should be initialised on LHE files called NameOfYourLesHouchesFile_0.lhe (+0 jet sample) and NameOfYourLesHouchesFile_1.lhe (+1 jet sample) and the same naming convention for LHE files with two or more additional jets. Since for this option, the merging scale value is read from the LHEF, no merging scale value needs to be supplied by setting Merging:TMS . Also, the hard process is read from LHEF, the input Merging::Process does not have to be defined. However, the maximal number of merged jets still has to be supplied by setting Merging:nJetMax.

flag  Merging:doMerging   (default = off)
Finally, it is also possible to switch on the merging machinery without switching on specific merging procedure, by turning this setting on. This is not sufficient for the internal merging of PYTHIA itself, but convenient when using an external shower plugin.

Merging with merging scale defined in Pythia evolution pT:

If the LHE file has been regularised by cutting on the minimal Pythia evolution pT in the state, this can also be used as a merging scale right away. For this, change the switch

flag  Merging:doPTLundMerging   (default = off)
The merging scale is then defined by finding the minimal Pythia evolution pT between sets of radiator, emitted and recoiler partons. For this particular merging scale definition, u,d,c,s,b,g are considered partons. The Pythia evolution pT of a single three-parton set is defined by

pTevol = zijk(1-zijk) Qij2 for FSR, where i is the radiating parton, j is the emitted parton and k is the recoiler, and Qij2 = (pi + pj)2 , and zijk = xi,jk / (xi,jk + xj,ik) with xi,jk = 2 pi (pi + pj + pk) / (pi + pj + pk)2

pTevol = (1-zijk) Qij2 for ISR, where i is the radiating parton, j is the emitted parton and k is the second initial state parton, and Qij2 = -(pi - pj)2 , and zijk = (pi - pj + pk)2 / (pi + pk)2 .

When using this option, the merging scale is defined by the minimum pTevol for all combinations of three partons in the event, irrespective of flavour or colour-connections. The merging scale value will be read from the Merging:TMS parameter, so that this needs to be set just as in the case of the kT-merging prescription. Of course you will also need to set Merging:Process and the maximal number of additional matrix element jets Merging:nJetMax.

Merging with merging scale defined by a combination of cuts:

It is possible to regularise QCD divergences in a LHE file by applying cuts to the kinematical pT of jets (pTi), combined with a cut on ΔRij between jets and a cut on invariant mass Qij of jet pairs. The combination of these standard cuts can also serve as a merging scale. For this, use this setting

flag  Merging:doCutBasedMerging   (default = off)
This switch will use cuts on (pTi), ΔRij and Qij to define when parton shower emissions are allowed. Please note for this particular merging scale definition, only light jets (u,d,c,s,g) are checked.

The values of the cuts will then be read from

parm  Merging:QijMS   (default = 0.0)
The value of the invariant mass cut Qij of pairs of final state partons used in the matrix element generation.

parm  Merging:pTiMS   (default = 0.0)
The value of the minimal transverse momentum cut pTi on final state partons, as used in the matrix element generation.

parm  Merging:dRijMS   (default = 0.0)
The value of the minimal ΔRij separation between pairs of final state partons used in the matrix element generation, where ΔRij2 = (Δyij)2 + (Δφij)2.

With knowledge of these values, and Merging:doCutBasedMerging, Pythia will use these cuts as a separation between matrix element phase space and parton shower region. If e.g. the Les Houches Events have been generated with the cuts ΔRij = 0.1 , pTi= 20 GeV and Qij = 40 GeV, set Merging:QijMS=40., Merging:pTjMS=20., Merging:dRijMS=0.1 to perform a cut-based merging. Of course you will also need to set Merging:Process and the maximal number of additional matrix element jets Merging:nJetMax.

Les Houches events outside the matrix element region

Before continuing, we would like to point out that in general, the user should make sure that the events in the Les Houches file are actually calculated using the regularisation cut definition and value(s) supplied to Pythia as merging scale definition and value(s). However, if LH files with a large number of events and loose merging scale cuts are available, it might be useful to choose a higher merging scale value, e.g. for merging scale variations as part of uncertainty assessments. If CKKW-L merging is enabled, Pythia will by default check if events read from Les Houches file are in the matrix element region as defined by the merging scale definition and merging scale value. Events outside the matrix element region will be discarded. This will lead to warnings of the form "Les Houches Event fails merging scale cut. Cut by rejecting event". These warnings should, in this case, rather be regarded as information. To change the default behaviour, use the flag

flag  Merging:enforceCutOnLHE   (default = on)
This will check if the events read from LHE file are in the matrix element region as defined by the merging scale definition and value(s). If on, LHE input outside the matrix element region will be rejected. If off, every event is assumed to pass the merging scale cut.

Defining the hard process

To perform CKKW-L matrix element merging, the user has to decide on a hard process, and communicate this choice to Pythia. This is done by setting the input Merging:Process. Note that the EW-improved "merging of mergings" strategy presented in [Chr15a] requires a very loose process definition that is discussed in the section Electroweak Merging below.

For single processes in the Standard Model or the MSSM a notation similar to MG4/MEis applicable. However, for some purposes, using a single simple process string is not satisfactory. Mixed W+ and W- events in a single LHE file is a common example. For this case, it would of course be perfectly allowed to perform twice, once for W+ production and once for W- production, and then add the results. Nevertheless, it seems reasonable to alleviate difficulties by allowing for less restrictive hard process definitions. Some generalisations of the process tag are available: Containers, user-defined particle tags and a very rudimentary option to "guess" the hard process. The syntax of these settings is described below.

In case you want multiple processes in a LHE file to be treated on equal footing (e.g. W+ + jets and W- + jets), you should use flexible containers do specify the hard process. So far, we allow the use of the containers LEPTONS, NEUTRINOS, BQUARKS. If you use these containers, the hard process definition is relatively flexible, meaning that Pythia will attempt a merging of QCD jets for each event in the LHE file, and assume that all particles matching one of the containers are products of the hard process. This is best explained by examples. If you want to have both pp → e+ νe + jets and pp → e- ν̄e + jets events in a single file, you can set Merging:Process=pp>LEPTONS,NEUTRINOS as hard process (note that for readability, you are allowed to use commata to separate container names). Combining e.g. the processes pp → e+ νe and pp → μ+ νμ is possible with the hard process definition pp>LEPTONS,NEUTRINOS.

For maximal flexibility, the definition of the hard process by these containers does not mean that each Les Houches event needs to contain particles to match each container. It is enough if one container is matched. This means that with the string pp>LEPTONS,NEUTRINOS, you can immediately process pp → e+ e- events mixed with pp → e+ νe events, since particles matching at least one container can be found in both cases. Another example for the usage of containers is mixing pp → e+ νe and pp → tt̄ → e+ νe e- ν̄e bb̄. This can be accommodated by the hard process string Merging:Process=pp>LEPTONS,NEUTRINOS,BQUARKS.

There is however a conceptual limitation to containers: The hard process definition is necessary to ensure that when constructing lower multiplicity states (that will be used to calculate the correct merging weight), the structure of the hard process will be preserved. If e.g. we want the hard process to be pp → Z → bb̄ , we should ensure that the lowest multiplicity state contains a colour-singlet bb̄-pair. When reconstructing intermediate lower multiplicity states from multi-jet matrix elements, we should thus always be able to find at least one bb̄-pair. By mixing different processes in a LHE file, this requirement might already be violated at the level of Les Houches events. Flexible containers cannot give strong conditions which flavours should be preserved in the construction of the hard process. In order to avoid non-sensible results, it is hence assumed that all particles matching any of the containers will be part of the lowest multiplicity process. This implies that if you decide to use the BQUARKS container, all b-quarks in the LHE file will be interpreted as hard process particles, and never as additional radiation.

Another way to specify the hard process particles is to explicitly define the particle names and identifiers. This is necessary if the matrix element merging in Pythia does not contain the particles of interest. To make sure that the hard process is still treated correctly, it is possible to define particles in the process string. If you e.g. want the hard process to contain a particle "zeta~" with PDG identifier "12345", produced in proton collisions, you have to include a user-defined particle tag by setting the process string to pp>{zeta~,12345}. The user-defined particle is enclosed in curly brackets, with syntax {particle_name,particle_identifier}, where "particle_name" and "particle_identifier" are the particle name and particle identifier used for this particle in the input LHE file. User-defined particles are only allowed in the final state. You are free to fix user-defined particles with more common ones, as long as user-defined particles are put before more common particles in the process string. This means that if you e.g. wanted the hard process to contain a graviton in association with a positron and an electron-neutrino, you have to define the hard process as pp>{G,39}e+ve.

Below you can find a list of particles predefined in the merging. If you wish to include a hard process with different final state particles, you may use the "curly bracket notation" outlined above.

The set of incoming particles us limited to: e- (electron), e+ (positron), mu- (muon), mu+ (antimuon), p (proton, container to hold all initial state coloured particles), p~ (identical to p container).

The following intermediate particles are allowed: a (photon), z (Z boson), w- (W- boson), w+ (W+ boson), h (scalar Higgs boson), W (container to hold both W- and W+ boson), t (top quark), t~ (anti-top), dl, dl~, ul, ul~, sl, sl~, cl, cl~, b1, b1~, t1, t1~, dr, dr~, ur, ur~, sr, sr~, cr, cr~, b2, b2~, t2, t2~ (all MSSM squarks).

We have pre-defined the outgoing particles: e+, e-, ve~, ve, mu+, mu-, vm~, vm, ta+, ta-, vt~, vt (all SM leptons and neutrinos), j~ (container to hold all final state coloured particles), j (container to hold all final state coloured particles), NEUTRINOS (container to hold all final state neutrinos and anti-neutrinos), LEPTONS (container to hold all final state leptons and anti-leptons), BQUARKS (container to hold final state b-quarks), d~, d, u~, u, s~, s, c~, c, b~, b, t~, t (all SM quarks), a, z, w-, w+ (all SM electro-weak bosons), h (scalar Higgs boson), W (container to hold both W- and W+ boson), n1 (MSSM neutralino), dl~, dl, ul~, ul, sl~, sl, cl~, cl, b1~, b1, t1~, t1, dr~, dr, ur~, ur, sr~, sr, cr~, cr, b2~, b2, t2~, t2 (all MSSM squarks). Other outgoing particles are possible if you use the "curly bracket notation" described earlier.

Finally, it is worth noting that including the header Pythia8Plugins/aMCatNLOHooks.h (as is e.g. the case in the main89.cc example) allows the user to instruct Pythia to "guess" the hard process by setting Merging:Process = guess. In this case, the code will inspect the hard process before starting the merging code, and set Merging:Process on an event-by-event basis according to the rules implemented in Pythia8Plugins/aMCatNLOHooks.h. It is instructive to understand these rules before using this feature. At present, the guess feature assumes that any non-colored particle, any particle with PDG ID larger than 21, and any quark not labelled as "additional parton" through the Merging:nQuarksMerge mode are part of the hard process definition. This e.g. enables usage for simple BSM processes. Please use this feature with caution.


Histogramming the events

After the event has been processed, histograms for observables of interest need to be filled. In order to achieve good statistical accuracy for all jet multiplicities and all subprocesses contributing to one jet multiplicity, generally a fixed number of unit-weighted events is read from each Les Houches Event file. To then arrive at the correct prediction, for each of these events, histogram bins should be filled with the corresponding cross section, or weighted with unit weight and normalised at the end to the generated cross section for each jet multiplicity separately.

Still another, even more important, event weight that has to applied on an event-by-event basis is the CKKW-L-weight. This corrective weight is the main outcome of the merging procedure and includes the correct no-emission probabilities, PDF weights and coupling (αs or αem) factors. This means that the merging implementation will generate weighted events. The CKKW-L-weight can be accessed by the following function:

double Info::mergingWeight()  
Returns the CKKW-L weight for the current event.

Note that to avoid confusion, this function does not include the the weight of a phase space point (given by Info::weight()). This weight will differ from unity when reading in weighted Les Houches events. In this case, the full weight with which to fill histogram bins is Info::mergingWeight() * Info::weight().

Finally, to arrive at a correct relative normalisation of the contributions from different number of additional jets in the matrix element, each histogram should be rescaled with the accepted cross section given by Info::sigmaGen(). The accepted cross section includes the effect of vetoes generating Sudakov form factors for the matrix elements, and is in general only known after the run.

This final step can of course be skipped if the accepted cross section had been estimated before the histogramming run, and histogram bins had instead been filled with the weight Info::mergingWeight() * σest(number of additional jets in current ME sample). This is the way HepMC events should be weighted to produce correct relative weights of events (see below, and particularly examine the example programs main84.cc and main85.cc).

Examples how to use these options are given in main81.cc (kT merging), main84.cc (automatic MG/ME merging for RIVET usage), and main85.cc (HepMC output for RIVET usage).


Merging with user-defined merging scale function

For all other merging scale definitions, the procedure is slightly more complicated, since the user has to write a small piece of code defining the merging scale. To allow for a user defined procedure, set the input

flag  Merging:doUserMerging   (default = off)
General user defined merging on/off.

Then, set the Merging:nJetMax, Merging:TMS and Merging:Process input as before.

Since during execution, PYTHIA needs to evaluate the merging scale with the definition of the user, the user interface is designed in a way similar to the UserHooks strategy. The class controlling the merging scale definition is called MergingHooks.

Initialisation

To initialise the merging with user-defined merging scale, we should construct a class derived from MergingHooks, with a constructor and destructor

MergingHooks::MergingHooks()  

virtual MergingHooks::~MergingHooks()  
The constructor and destructor do not need to do anything.

For the class to be called during execution, a pointer to an object of the class should be handed in with the
Pythia::setMergingHooksPtr( MergingHooks*) method. An examples of this procedure are given in main82.cc.

Defining a merging scale

Then, in the spirit of the UserHooks class, the user needs to supply the process to be merged by defining a methods to evaluate the merging scale variable.

virtual double MergingHooks::tmsDefinition(const Event& event)  
This method will have to calculate the value of the merging scale defined in some variable from the input event record. An example of such a function is given in main82.cc.

The base class MergingHooks contains many functions giving information on the hard process, to make the definition of the merging scale as easy as possible:

int MergingHooks::nMaxJets()  
Return the maximum number of additional jets to be merged.

int MergingHooks::nHardOutPartons()  
Returns the number of outgoing partons in the hard core process.

int MergingHooks::nHardOutLeptons()  
Returns the number of outgoing leptons in the hard core process.

int MergingHooks::nHardInPartons()  
Returns the number of incoming partons in the hard core process.

int MergingHooks::nHardInLeptons()  
Returns the number of incoming leptons in the hard core process.

int MergingHooks::nResInCurrent()  
The number of resonances in the hard process reconstructed from the current event. If e.g. the ME configuration was pp → (w+→e+ve)(z → mu+mu-)jj, and the ME generator put both intermediate bosons into the LHE file, this will return 2.

double MergingHooks::tms()  
Returns the value used as the merging scale.

Filling output histograms for the event then proceeds along the lines described above in "Histogramming the events".

The full procedure is outlined in main82.cc. Special care needs to be taken when the output is stored in the form of HepMC files for RIVET usage.

Defining a cut on lowest jet multiplicity events

It can sometimes happen that when generating LHE files, a fairly restrictive cut has been used when generating the lowest multiplicity matrix element configurations. Then, it can happen that states that are (in the generation of a parton shower history) constructed by reclustering from higher multiplicity configurations, do not pass this matrix element cut.

Consider as an example pure QCD dijet merging, when up to one additional jet should be merged. Three-jet matrix element configurations for which the reclustered two-jet state does not pass the cuts applied to the two-jet matrix element would never have been produced by showering the two-jet matrix element. This means that the three-jet matrix element includes regions of phase space that would never have been populated by the parton shower. Thus, since the matrix element phase space is larger than the shower phase space, merging scale dependencies are expected. A priori, this is not troublesome, since the aim of matrix element merging is to include regions of phase space outside the range of the parton shower approximation into the shower. An example is the inclusion of configurations with only unordered histories.

Clearly, if the parton shower phase space is very constrained by applying stringent cuts to the two-jet matrix element, merging scale dependencies can become sizable, as was e.g. seen in [Lon11] when forcing shower emissions to be ordered both in the evolution variable and in rapidity. To influence the effect of large phase space differences for shower emissions and matrix element configurations due to LHEF generation cuts, the user has to write a small piece of code overwriting method

virtual double MergingHooks::dampenIfFailCuts(const Event& event)  
multiplicity reclustered state as an input Event. From this input event, the user can then check if matrix element cuts are fulfilled. The return value will be internally multiplied to the CKKW-L weight of the current event. Thus, if the user wishes to suppress contributions not passing particular cuts, a number smaller than unity can be returned.

Note that this method gives the user access to the lowest multiplicity state, which ( e.g. in the case of incomplete histories) does not have to be a 2 → 2 configuration. Also, changing the weight of the current event by hand is of course a major intervention in the algorithm, and should be considered very carefully. Generally, if this facility would have to be used extensively, it is certainly preferable to be less restrictive when applying additional, non-merging-scale-related cuts to the matrix element.

An example how to force a cut on lowest multiplicity reclustered states for pure QCD matrix element configurations is given by main83.cc (to be used with e.g. main82.cmnd).

Influencing the construction of all possible histories

Even more powerful - and dangerous - is influencing the construction of histories directly. This should only be attempted by expert users. If you believe manipulations completely unavoidable, we advise you to take great care when redefining the following functions.

virtual bool MergingHooks::canCutOnRecState()  
In the base class this method returns false. If you redefine it to return true then the method doCutOnRecState(...) will be called for each reclustered state encountered in the generation of all possible histories of the matrix element state.

virtual bool MergingHooks::doCutOnRecState(const Event& event)  
This routine will be supplied internally with every possible reclustered event that can be reached by reclustering any number of partons in the matrix element input state. The new, reclustered, states can then be analysed. If the method returns false, the history to which the state belongs will be treated as if it were unordered, i.e. this path will only be chosen if no other histories are available. In this way, the number of histories not fulfilling the user criterion will be minimised.

Clearly, these methods are highly intrusive. It could e.g. happen that no history is allowed, which would make merging impossible. One example where this method could be useful is if cuts on the core 2 → 2 processes have to be checked, and the method MergingHooks::dampenIfFailCuts(const Event& event) is not sufficiently effective.

Defining the hard process matrix element

The MergingHooks class also allows the expert user to define the matrix element of the hard process, by defining the method

virtual double MergingHooks::hardProcessME(const Event& inEvent)  
This routine will be supplied internally with the reconstructed lowest-multiplicity event. From this, it is possible to calculate the squared matrix element of the hard process, by using the information stored in the event record. The function should return a double value that corresponds to the matrix element at the phase space point given by the input event record. This number will then be multiplied to the product of splitting functions that define the probability of the current path of the parton shower history. In this way, the hard process configuration can be taken into account when choosing the parton shower history, which is, internally, used to generate the "merging weight".

The inclusion of the hard process matrix element into the choice of histories becomes relevant when the hard process matrix element has very strong phase space dependencies. QCD dijet cross sections for example strongly depend on the transverse momentum of the jets. So far, the authors have not encountered any changes upon inclusion of the full hard process matrix element, even for the QCD dijet case.


Matrix element merging and HepMC output for RIVET

Examples how to produce matrix element merged events to be analysed with RIVET are given by main84.cc and main85.cc.

The main issue is that the output of separate RIVET runs can not in general be combined. To perform a matrix element merging, we however need to runs over different LHE files. The solution to this problem (so far) is to only perform one RIVET run for all matrix elements, i.e. print the events for all ME parton multiplicities, with the correct weights, to a single HepMC file. Since the correct weight includes the cross section of the different samples after Sudakov vetoes --- which is not a priori known --- the cross sections have to be estimated in a test run, before the actual production run is performed. Finally, the cross section of the last event in the HepMC file has to be taken as the full merged cross section sigma_merge = Sum_{i=0}^N Sum_{j=0}*^{nEvents} sigma_est(i)*wckkwl(j).

This procedure is outlined in main84.cc. Input LHE files with only very inclusive cuts pose further difficulties. For such files (which were already addressed under the heading Les Houches events outside the matrix element region), the cross section after the merging scale cut is not known before the cut is performed. Using Pythia's UserHooks facilities, it is possible to produce a valid estimate of the cross section after cuts. This however entails a careful cut definition by the user, which might become cumbersome for some in-built merging scale definitions. A reasonable alternative is using the switch

flag  Merging:doXSectionEstimate   (default = off)
If on, estimate cross section after merging scale cut. This switch has to be used in conjunction with a merging scale definition (e.g. Merging:doPTLundMerging = on). Then, this merging scale definition will be used as a cut on the input events. After the requested number of Monte Carlo events, the cross section after the cut can be extracted by inferring the Info::sigmaGen() method, and the number of accepted events by using Info::nAccepted()

This switch also relies on knowledge on how many partons a LHE file should contain. This is important for real-emission kinematics in the case of NLO merging. The number of (additional) partons in a LHE file can be set with

mode  Merging:nRequested   (default = -1; minimum = -1)
Exact number of additional jets requested for a particular LHE file. If a file should for example only contain W+ g g events, this switch should be set to "2" for this LHE file. For NLO merging schemes (see NLO Merging), this number has to be set.

The usage of these switches to obtain the necessary cross section estimate is illustrated in main85.cc. The example main85.cc program is intended as a "front-end" for CKKW-L merging in Pythia8, so we will discuss the program briefly. main85.cc should be used together with an input file (like main85.cmnd). The executable should be invoked with three arguments: the input file, the "name" of the input LHE files, and the name of the output HepMC event file. To use the LHE files that are shipped with the Pythia distribution, a valid usage would be

./main85.exe ./main85.cmnd ./w_production ./myhepmc.hepmc

If you want to use other input LHE files, note that main85.cc assumes the naming convention name_tree_#nAdditionalJets.lhe. All settings can be included though the input file, so that main85.cc does not have to be changed explicitly. main85.cc first switches off showers, multiparton interactions and hadronization, and estimates the cross sections (after enforcing the merging scale cut) of samples for different numbers of additional jets. Then, showers, MPI and hadronization are switched on again, and the CKKW-L merging procedure is performed. Events will be read in a decreasing sequence of jet multiplicities, which also means that e.g. events with two additional partons in the LHE file will be printed to the HepMC file before events with one additional parton.


Electroweak Merging

Merging strategies like CKKW-L usually assume that the description of a (relatively simple) underlying process should be improved by combining with states that contain additional well-separated partons - with "additional" measured with respect to the underlying process. As discussed in [Chr15a], this philosophy is not always sensible, and may lead to an unconvincing physics model. The bias can be greatly reduced by considering that in perturbation theory, corrections to seemingly very different underlying processes mix, so that there is no justification to classify some states as corrections to only one underlying process. Interactions that only contain vertices of only one theory (e.g. QCD) will mix with processes that contain only vertices of another interactions (e.g. QED). Underlying processes with very different coupling structures should thus be considered. This is the main aim of the electroweak merging scheme. The process p p → W jet jet provides a good example, since it can be interpreted either as double-real-QCD-emission correction to p p → W or as real-electroweak-emission correction to p p → jet jet. The distinction is artificial, but the all-order resummation is very different in either case, leading to distinctly different predictions. Thus, a minimally biased method for assigning an underlying process has to be found.

The method of [Chr15a] chooses the underlying process probabilistically for each phase space point based on a product of splitting kernels and full hard process matrix elements, and includes the correct all-order factors after this choice. In our previous example ( p p → W jet jet), this would mean that the hard scattering could be p p → W or p p → jet jet (or, depending on phase space considerations, also p p → W jet). The availability of electro-weak showers in Pythia 8 is crucial for a consistent scheme that removes of double-counting, since double-counting is (in part) removed by applying no-emission probabilities. In order to e.g. ensure that a state in the p p → jet jet does not, via W-boson emission, evolve into a state overlapping with p p → W + two QCD emissions, the former has to be reweighted with an all-order no-electroweak-emission probability. As always, a small merging scale dependence is facilitated by accounting for for dynamical PDF evaluation and running couplings. An electroweak merging thus includes a reweighting with αem ratios that are automatically included in the "merging weight".

As another consequence of probabilistically assigning the underlying process is the "merging of mergings": Since purely partonic final states can evolve into jets + electroweak bosons states, it is necessary to treat bosons and partons on equal footing, meaning that well-separated boson states should be corrected with fixed-order inputs, while soft/collinear bosons should be associated with parton showering. Shortly, bosons and partons are treated identically. For our previous example, this means that the "correct" set of fixed-order corrections that include up to three final state particles is p p → n Ws + m jets, where any combination of n and m subject to n + m <= 3 has to be included ( p p → W, p p → W jet, p p → W jet jet, p p → W W, p p → W W W, p p → W W jet, p p → jet jet, p p → jet jet jet) For practical purposes, it is sometimes permitted to not combine a complete set of processes. Only the single W-state has been explicitly validated and in addition the weak PS does not include all possible splittings for multiple W emissions, therefore caution has to be taken if using this for multiple W states.

The "merging of mergings" has important consequences for the (interested) user. Below, we give instructions on the usage of the electroweak merging.

Fixed-order inputs for electroweak merging

The electroweak merging leads to the idea of a "merging of mergings". This has to be enforced also at the fixed-order sample generation stage, with two main requirements.

A fully consistent treatment requires the generation of samples containing all states with a number of emissions that is less than or equal to the "maximal possible number of emissions of any type" that should be corrected, cf. the W jet jet example above. In practise, it is often permitted to disregard some (set of) samples since their impact on an analysis is negligible. If e.g. an analyis always requires missing transverse momentum and a single lepton, and we assume perfect lepton acceptances, then it would be permitted to disregard the multi-W samples in the example. If the collision energy is in addition low (~up to fews of TeVs), then the probability for a pure QCD state to emit W-bosons is often low enough so that the pure jet samples can be neglected. However, you should think very carefully before settling on any shortcuts.

Any particle that could count as an emission has to be included in the calculation of the particle separations that define the merging scale. If e.g. W-bosons are considered emissions, then any state with W-bosons that are collinear with another emission should be removed from the fixed-order sample. Such configurations will instead be produced through parton showering. This requirement means that you might have to define your own cut including this condition in your favourite fixed-order matrix element generator. At the risk of losses in efficiency, you can also use samples with very loose cuts and have Pythia enforce the merging scale cut when reading your input events. The latter is only possible if you use the merging scale definition Merging:doPTLundMerging = on.

Enabling the electroweak merging

The electroweak merging is currently only tested for processes containing W-bosons and jets. For a consistent merging, it is necessary to enable W-boson emissions by using TimeShower:weakShower = on, TimeShower:weakShowerMode = 1, SpaceShower:weakShower = on, SpaceShower:weakShowerMode = 1, and WeakShower:externalSetup = on.

To enable the electroweak merging, use the following switch.

flag  Merging:allowWeakClustering   (default = off)
Allow clustering of weak bosons, as necessary if a merging of matrix elements with QCD and weak showering is attempted. Currently, only emissions of W-bosons are accounted for. This switch should only be used if weak showering is turned on (see Timelike showers, Spacelike showers and Weak showers for details on weak showering).

As explained above, the electroweak (EW) merging treats all shower-producible particles equally. This also means that the merging scale definition must include a separation of W-bosons and partons to define if a state is in the (well-separated) fixed-order region or if it is in the (soft/collinear) parton shower region. Such a cut can be implemented inside the (external) matrix element generator. On the other hand, Pythia 8 allows the usage of samples with very loose cuts and can enforce the correct merging scale cut by rejecting input events that do not pass the cut. This is also possible for the EW merging, albeit only for the merging scale definition that is enabled by using Merging:doPTLundMerging = on. We recommend using this strategy for users that do not wish to implement the cut directly into the ME generation, and who are prepared to accept a loss of efficiency because of Pythia's a-posteriori rejection.

Defining the "inclusive" merging process

Since the concept of a single "hard process" is not suitable for the EW merging, the process should be defined in a rather loose manner. This loose definition is still done by setting the input Merging:Process. Processes for EW merging should use the containers Jinc, Winc, Ainc and Zinc, which tell the code which particles could be possible "additional" emissions. No other particle defnitions are allowed, and none of the settings discussed in the "Defining the hard process" section are relevant here. Examples of allowed process definitions are

Merging:Process = pp > Jinc,Winc meaning that W-bosons and partons are treated on equal footing (i.e. this is the setting applicable to the example used earlier). The merging will then include pure QCD multijet events, W+jets events, multi-W+jets events and pure multi-W events;

Merging:Process = pp > Jinc,Zinc meaning that Z-bosons and partons are treated on equal footing;

Merging:Process = pp > Jinc,Winc,Zinc meaning that W-bosons, Z-bosons and partons are treated on equal footing.

Setting the number of additional particles

Since the EW merging probabilistically decides on the "underlying process", it is a priori not possible to set the maximal number of additional emissions on top of this underlying process. A p p → W jet jet state would e.g. contain two additional emissions if interpreted as correction to p p → W, and only one additional emission if interpreted as correction to p p → jet jet. Pythia 8 consequently decides dynamically how to set the additional number of emissions.

The maximal number of emissions, set by using the Merging:nJetMax setting, still has to be defined to allow a sensible treatment of the "highest-multiplicity" states. We thus redefine the meaning of Merging:nJetMax to "maximal possible number of emissions of any type". As an example, Merging:nJetMax = 3 if you want to perform a "merging of mergins" containing states with up to three partons, or up to two partons and one W-boson, or up to one parton and two W-bosons, or up to three W-bosons.


Further variables

For more advanced manipulations of the merging machinery, all parameter changes that were investigated in [Lon11] are supplied. Please check [Lon11] for a detailed discussion of the switches.

These switches allow enthusiastic users to perform a systematic assessment of the merging prescription. Apart from this, we advise the non-expert user to keep the default values.

mode  Merging:nQuarksMerge   (default = 5; minimum = 2; maximum = 5)
This switch controls which quarks flavours (labelled by PDG id's) are considered additional partons. If e.g. set to 4, then u-, d-, c- and s-quarks will be merged, while b-quarks will not be considered in the merging (corresponding to a 4-flavour merging scheme). We advise caution when changing this number. In particular, please ensure that the allowed flavour for additional partons in the input LHE file does not exceed this value, since unnecessary double-counting might occur otherwise.

flag  Merging:includeMassive   (default = on)
If on, use the correct massive evolution variable and massive splitting kernels in the reconstruction and picking of parton shower histories of the matrix element. If off, reconstruct evolution scales, kinematics and splitting kernels as if all partons were massless.

flag  Merging:enforceStrongOrdering   (default = off)
If on, preferably pick parton shower histories of the matrix element which have strongly ordered consecutive splittings, i.e. paths in which consecutive reclustered evolution scales are separated by a user-defined factor.

parm  Merging:scaleSeparationFactor   (default = 1.0; minimum = 1.0; maximum = 10.0)
The factor by which scales should differ to be classified as strongly ordered.

flag  Merging:orderInRapidity   (default = off)
If on, preferably pick parton shower histories of the matrix element with consecutive splittings ordered in rapidity and pT.

flag  Merging:pickByFullP   (default = on)
If on, pick parton shower histories of the matrix element by the full shower splitting kernels, including potential ME corrections and Jacobians from joined evolution measures.

flag  Merging:pickByPoPT2   (default = off)
If on, pick parton shower histories of the matrix element by the shower splitting kernels divided by the evolution pT.

flag  Merging:pickBySumPT   (default = off)
If on, exclusively pick parton shower histories of the matrix element for which have the smallest sum of scalar evolution pT for consecutive splittings has been calculated.

flag  Merging:includeRedundant   (default = off)
If on, then also include PDF ratios and αs factors in the splitting probabilities used for picking a parton shower history of the matrix element, when picking histories by the full shower splitting probability. As argued in [Lon11], this should not be done since a reweighting with PDF ratios and αs factors will be performed. However, it can give useful insight in how sensitive the results are to the prescription on how to choose PS histories.

parm  Merging:nonJoinedNorm   (default = 1.0; minimum = 0.0; maximum = 10.0)
Normalisation factor with which to multiply splitting probability for splittings without joined evolution equation.

parm  Merging:fsrInRecNorm   (default = 1.0; minimum = 0.0; maximum = 10.0)
Normalisation factor with which to multiply splitting probability for final state splittings with an initial state recoiler.

parm  Merging:aCollFSR   (default = 1.0; minimum = 0.0; maximum = 10.0)
Factor with which to multiply the scalar pT of a final state splitting, when choosing the history by the smallest sum of scalar pT. Default value taken from Herwig++ [Tul09].

parm  Merging:aCollISR   (default = 0.9; minimum = 0.0; maximum = 10.0)
Factor with which to multiply the scalar pT of an initial state splitting, when choosing the history by the smallest sum of scalar pT. Default value taken from Herwig++ [Tul09].

mode  Merging:unorderedScalePrescrip   (default = 0; minimum = 0; maximum = 1)
When the parton shower history of the matrix element contains a sequence of splittings which are not ordered in evolution pT (called an unordered history), this sequence is interpreted as a combined emission. Then, a decision on which starting scale for trial emissions off reconstructed states in this sequence of unordered splittings has to be made. Two options are available:
option 0 : Use larger of the two reconstructed (unordered) scales as starting scale.
option 1 : Use smaller of the two reconstructed (unordered) scales as starting scale.

mode  Merging:unorderedASscalePrescrip   (default = 1; minimum = 0; maximum = 1)
Prescription which scale to use to evaluate αs weight for splittings in a sequence of splittings which are not ordered in evolution pT.
option 0 : Use the combined splitting scale as argument in αs, for both splittings.
option 1 : Use the true reconstructed scale as as argument in αs, for each splitting separately.

mode  Merging:unorderedPDFscalePrescrip   (default = 0; minimum = 0; maximum = 1)
Prescription which scale to use to evaluate ratios of parton distributions for splittings in a sequence of splittings which are not ordered in evolution pT.
option 0 : Use the combined splitting scale as argument in PDF ratios, for both splittings.
option 1 : Use the true reconstructed scale as argument in PDF ratios, for each splitting separately.

mode  Merging:incompleteScalePrescrip   (default = 0; minimum = 0; maximum = 2)
When no complete parton shower history (i.e. starting from a 2 → 2 process) for a matrix element with additional jets can be found, such a configuration is said to have an incomplete history. Since in incomplete histories, not all shower starting scales are determined by clusterings, a prescription for setting the starting scale of trial showers in incomplete histories is needed. Three options are provided.
option 0 : Use factorisation scale as shower starting scale for incomplete histories.
option 1 : Use sHat as shower starting scale for incomplete histories.
option 2 : Use s as shower starting scale for incomplete histories.

flag  Merging:allowColourShuffling   (default = off)
If on, this will allow the algorithm to swap one colour index in the state, when trying to find all possible clusterings, if no clustering has been found, but more clusterings had been requested. In this way, some incomplete histories can be avoided. Generally, we advise the non-expert user to not touch this switch, because a slight change in the colour structure can change the radiation pattern. To however study the sensitivity of the predictions on these effects, allowing for colour reshuffling can be useful.

flag  Merging:usePythiaQRenHard   (default = on)
If on, this will allow the algorithm to use a dynamical renormalisation scale to evaluate the strong couplings of the core hard process in dijet and prompt photon events. This means that the value of αs used as coupling of the hard process in the matrix element generation will be replaced with a running coupling evaluated at the geometric mean of the squared transverse masses of the two outgoing particles, as is the default prescription in Pythia.

flag  Merging:usePythiaQFacHard   (default = on)
If on, this will allow the algorithm to use a dynamical factorisation scale to evaluate parton distributions associated with the hadronic cross section of the core hard process in dijet and prompt photon events. In the calculation of PDF ratios as part of the CKKW-L weight of an event, parton distributions that should be evaluated at the scale of the core 2 - >2 process will be evaluated using the dynamical factorisation scale Pythia would attribute to this process. This means that the hard process factorisation scale is set to the smaller of the squared transverse masses of the two outgoing particles.

flag  Merging:mayRemoveDecayProducts   (default = off)
Remove products of resonances in the hard process, in case Pythia generates decay products before merging. This makes merging possible even for an indeterminate final state, if Pythia itself has produced the decay products. The merging methods will instead be invoked on the "non-decayed" event, thus removing the limitation to only one decay channel when performing the merging. This switch is necessary e.g. for slepton pair production in association with additional QCD jets, if the input LHE file contains the resonant sleptons, and Pythia decides on a decay according to the branching fractions read from SLHA input.

flag  Merging:allowSQCDClustering   (default = off)
Allow clustering of gluon emission off squarks.

flag  Merging:useShowerPlugin   (default = off)
Use the splitting probabilities, evolution variables and phase space mappings of an external shower plugin. This will become possible as soon as new showers containing the necessary ingredients are available in Pythia.

flag  Merging:applyVeto   (default = on)
If off, no event veto based on the merging scale is applied in CKKW-L merging. This means that the user has to implement the veto by hand in the Pythia main program. It can be useful to postpone event vetoes for the purpose of merging scale variations.

flag  Merging:includeWeightInXsection   (default = on)
If on, then the reweighting of events in the CKKW-L scheme is included in the event weight Info::weight(), the merging weight Info:mergingWeight() is unity, and the cross section printed by Info::sigmaGen() includes the effect of CKKW-L merging.

flag  Merging:runtimeAMCATNLOInterface   (default = off)
This setting should only be touched by experts. If on, then the merging machinery will only be used to generate and extract information necessary for the runtime interface of Pythia and aMC@NLO. This is usually done via runtime calls by aMC@NLO. Read further details about the method in [Fre20].