HEJ 2.2.2
High energy resummation for hadron colliders
|
HEJ 2 is a library for all-order resummation of high-energy logarithms. It includes a program to add resummation to fixed-order events. User documentation for the program can be found here. This documentation is instead aimed at users of the library itself.
The main functionality is contained in the HEJ namespace. Particles are defined via the Particle struct, which consists of the particle four-momentum, its identifier according to the PDG Monte Carlo numbering scheme and an optional Colour charge. Given a number of incoming and outgoing particles, the square of the resummation matrix element can be calculated with the help of the MatrixElement class.
The EventReweighter class adds resummation to existing fixed-order events. Both fixed-order and resummation events are objects of the Event class, which are created from EventData objects with the help of a jet definition according to the fastjet library. EventData objects can be assembled manually or converted from input events in the LesHouches standard, read from file with a EventReader (e.g. LesHouchesReader or HDF5Reader).
Events can be saved with one of the EventWriter classes. Currently, there is support for the Les Houches event file format with the LesHouchesWriter class. If HEJ 2 was installed with HepMC 2 or 3 support, the respective format is available through the HepMC2Writer or HepMC3Writer class. The HDF5Writer class can be used if HEJ 2 was compiled with HDF5 support.
Further classes of interest are the interfaces to the Mixmax and Ranlux64 random number generators, the PDF class to interact with LHAPDF and the ScaleGenerator and ScaleConfig classes to calculate renormalisation and factorisation scales for a given Event.
As an example, we show a toy program that computes the square of a matrix element in the HEJ approximation for a single event. First, we include the necessary header files:
We then specify the incoming and outgoing particles. A particle has a type, a four-momentum \((p_x, p_y, p_z, E)\) and optionally a colour charge. For instance, an incoming gluon could be defined as
We collect all incoming and outgoing particles in a partonic event. Here is an example for a partonic \(gu \to gghu\) event (omitting colours):
Alternatively, we could read the event from e.g. a Les Houches event file, possibly compressed with gzip. For this, the additional header file HEJ/EventReader.hh
have to be included.
In this specific example we will later choose a constant value for the strong coupling, so that the HEJ matrix element does not depend on the renormalisation scale. However, in a more general scenario, we will want to set a central scale:
It is possible to add more scales in order to perform scale variation:
In the next step, we leverage FastJet to construct an event with clustered jets. Here, we use antikt jets with R=0.4 and transverse momenta of at least 30 GeV.
In order to calculate the Matrix element, we now have to fix the physics parameters. For the sake of simplicity, we assume an effective coupling of the Higgs boson to gluons in the limit of an infinite top-quark mass and a fixed value of \(\alpha_s = 0.118\) for the strong coupling.
If QCDLoop is installed, we can also take into account the full loop effects with finite top and bottom quark masses:
Finally, we can compute and print the square of the matrix element with
In the case of scale variation, the weight associated with the scale event.variations
[i].mur is ME(event)
.variations[i].
Collecting the above pieces, we have the following program:
After saving the above code to a file matrix_element.cc
, it can be compiled into an executable matrix_element
with a suitable compiler. For example, with g++
this can be done with the command
If HEJ or any of the required libraries was installed to a non-standard location, it may be necessary to explicitly specify the paths to the required header and library files. This can be done with the HEJ-config
executable and similar programs for the other dependencies: