hej is hosted by Hepforge, IPPP Durham
HEJ 2 2.0
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
HEJ 2

Introduction

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.

Overview

The main functionality is contained in the HEJ namespace. Particles are defined via the Particle struct, which consists of the particle four-momentum and its identifier according to the PDG Monte Carlo numbering scheme . 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 UnclusteredEvent objects with the help of a jet definition according to the fastjet library. UnclusteredEvent objects can be assembled manually or converted from input events in the LesHouches standard, read from file with a LHEF::Reader.

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 HepMCWriter class.

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.

Example

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:

#include "HEJ/Event.hh"
Declares the Event class and helpers.
Contains the MatrixElement Class.

We then specify the incoming and outgoing particles. A particle has a type and four-momentum $(p_x, p_y, p_z, E)$. For instance, an incoming gluon could be defined as

fastjet::PseudoJet momentum{0, 0, 308., 308.};
HEJ::Particle gluon_in{HEJ::ParticleID::gluon, momentum};
Class representing a particle.
Definition: Particle.hh:19

We collect all incoming and outgoing particles in a partonic event. Here is an example for a partonic $gu \to gghu$ event:

HEJ::UnclusteredEvent partonic_event;
// incoming particles
partonic_event.incoming[0] = {
HEJ::ParticleID::gluon,
{ 0., 0., 308., 308.}
};
partonic_event.incoming[1] = {
HEJ::ParticleID::up,
{ 0., 0.,-164., 164.}
};
// outgoing particles
partonic_event.outgoing.push_back({
HEJ::ParticleID::higgs,
{ 98., 82., 14., 180.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::up,
{ 68.,-54., 36., 94.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::gluon,
{-72., 9., 48., 87.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::gluon,
{-94.,-37., 46., 111.}
});
An event before jet clustering.
Definition: Event.hh:63
std::vector< Particle > outgoing
Definition: Event.hh:70
std::array< Particle, 2 > incoming
Definition: Event.hh:69

Alternatively, we could read the event from a Les Houches event file, possibly compressed with gzip. For this, the additional header files HEJ/stream.hh and LHEF/LHEF.h have to be included.

HEJ::istream in{"events.lhe.gz"};
LHEF::Reader reader{in};
reader.readEvent();
HEJ::UnclusteredEvent partonic_event{reader.hepeup};
Small wrapper around boost's filtering_istream.
Definition: stream.hh:20

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:

partonic_event.central.mur = 50.;
double mur
Definition: Event.hh:38
EventParameters central
Central parameter (e.g. scale) choice.
Definition: Event.hh:74

It is possible to add more scales in order to perform scale variation:

partonic_event.variations.resize(2);
partonic_event.variations[0].mur = 25.;
partonic_event.variations[1].mur = 100.;
std::vector< EventParameters > variations
Definition: Event.hh:75

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.

const fastjet::JetDefinition jet_def{
fastjet::JetAlgorithm::antikt_algorithm, 0.4
};
const double min_jet_pt = 30.;
HEJ::Event event{partonic_event, jet_def, min_jet_pt};
Definition: Event.hh:84

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.

const auto alpha_s = [](double /* mu_r */) { return 0.118; };
// whether to include corrections from the
// evolution of \alpha_s in virtual corrections
ME_config.log_correction = false;
HEJ::MatrixElement ME{alpha_s, ME_config};
Class to calculate the squares of matrix elements.
Definition: MatrixElement.hh:28
Configuration options for the MatrixElement class.
Definition: config.hh:116
bool log_correction
Whether to include the logarithmic correction from running.
Definition: config.hh:118

If QCDLoop is installed, we can also take into account the full loop effects with finite top and bottom quark masses:

ME_config.Higgs_coupling.mt = 163;
ME_config.Higgs_coupling.include_bottom = true;
ME_config.Higgs_coupling.mb = 2.8;
double mt
Top quark mass.
Definition: HiggsCouplingSettings.hh:16
bool include_bottom
Whether to include bottom quark effects.
Definition: HiggsCouplingSettings.hh:22
double mb
Bottom quark mass.
Definition: HiggsCouplingSettings.hh:18
bool use_impact_factors
Whether to use impact factors.
Definition: HiggsCouplingSettings.hh:20
HiggsCouplingSettings Higgs_coupling
Settings for effective Higgs-gluon coupling.
Definition: config.hh:120

Finally, we can compute and print the square of the matrix element with

std::cout << "HEJ ME: " << ME(event).central << '\n';

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:

#include "HEJ/Event.hh"
int main(){
HEJ::UnclusteredEvent partonic_event;
// incoming particles
partonic_event.incoming[0] = {
HEJ::ParticleID::gluon,
{ 0., 0., 308., 308.}
};
partonic_event.incoming[1] = {
HEJ::ParticleID::up,
{ 0., 0.,-164., 164.}
};
// outgoing particles
partonic_event.outgoing.push_back({
HEJ::ParticleID::higgs,
{ 98., 82., 14., 180.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::up,
{ 68.,-54., 36., 94.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::gluon,
{-72., 9., 48., 87.}
});
partonic_event.outgoing.push_back({
HEJ::ParticleID::gluon,
{-94.,-37., 46., 111.}
});
const fastjet::JetDefinition jet_def{
fastjet::JetAlgorithm::antikt_algorithm, 0.4
};
const double min_jet_pt = 30.;
HEJ::Event event{partonic_event, jet_def, min_jet_pt};
const auto alpha_s = [](double /* mu_r */) { return 0.118; };
// whether to include corrections from the
// evolution of \alpha_s in virtual corrections
ME_config.log_correction = false;
HEJ::MatrixElement ME{alpha_s, ME_config};
std::cout
<< "HEJ ME: " << ME(event).central
<< " = tree * virtual = " << ME.tree(event).central
<< " * " << ME.virtual_corrections(event).central
<< '\n';
}

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

g++ -o matrix_element matrix_element.cc -lHEJ -lfastjet

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:

g++ $(fastjet-config --cxxflags) $(HEJ-config --cxxflags) -o matrix_element matrix_element.cc $(HEJ-config --libs) $(fastjet-config --libs)