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

Main HEJ 2 Namespace. More...

Namespaces

namespace  currents
 
namespace  detail
 
namespace  detail_HepMC
 
namespace  event_type
 Namespace for event types.
 
namespace  helicity
 
namespace  pid
 particle ids according to PDG
 
namespace  Version
 

Classes

struct  Analysis
 Analysis base class. More...
 
struct  Beam
 Beam parameters. More...
 
class  BufferedEventReader
 Event reader with internal buffer. More...
 
class  CombinedEventWriter
 Write event output to zero or more output files. More...
 
struct  Config
 
class  CrossSectionAccumulator
 Sum of Cross Section for different subproccess. More...
 
struct  DefaultRNG
 Helper struct with default implementations. More...
 
struct  EmptyAnalysis
 
class  Event
 An event with clustered jets. More...
 
struct  EventParameters
 Event parameters. More...
 
struct  EventReader
 Abstract base class for reading events from files. More...
 
class  EventReweighter
 Main class for reweighting events in HEJ. More...
 
struct  EventReweighterConfig
 Configuration options for the EventReweighter class. More...
 
class  EventWriter
 Pure abstract base class for event writers. More...
 
class  EWConstants
 Collection of electro-weak constants. More...
 
class  FixedScale
 Functor that returns a fixed scale regardless of the input event. More...
 
class  Fraction
 Class for floating point numbers between 0 and 1. More...
 
class  HDF5Reader
 Class for reading events from a file in the HDF5 file format. More...
 
class  HDF5Writer
 This is an event writer specifically for HDF5 output. More...
 
class  HepMC2Interface
 This class converts the Events into HepMC::GenEvents. More...
 
class  HepMC2Writer
 This is an event writer specifically for HepMC output. More...
 
class  HepMC3Interface
 This class converts the Events into HepMC3::GenEvents. More...
 
class  HepMC3Writer
 This is an event writer specifically for HepMC3 output. More...
 
struct  HiggsCouplingSettings
 Settings for Higgs boson coupling to gluons. More...
 
struct  invalid_type
 Exception indicating wrong option type. More...
 
class  istream
 Small wrapper around boost's filtering_istream. More...
 
struct  JetParameters
 Jet parameters. More...
 
class  JetSplitter
 Class to split jets into their constituents. More...
 
class  LesHouchesReader
 Class for reading events from a file in the Les Houches Event File format. More...
 
class  LesHouchesWriter
 Class for writing events to a file in the Les Houches Event File format. More...
 
class  MatrixElement
 Class to calculate the squares of matrix elements. More...
 
struct  MatrixElementConfig
 Configuration options for the MatrixElement class. More...
 
struct  missing_option
 Exception indicating missing option setting. More...
 
class  Mixmax
 MIXMAX random number generator. More...
 
struct  not_implemented
 Exception indicating functionality that has not been implemented yet. More...
 
struct  OutputFile
 Output file specification. More...
 
struct  ParameterDescription
 Description of event parameters, see also EventParameters. More...
 
struct  Parameters
 Collection of parameters, e.g. Weights, assigned to a single event. More...
 
struct  PartialUnweightConfig
 Settings for partial unweighting. More...
 
struct  Particle
 Class representing a particle. More...
 
struct  ParticleProperties
 collection of basic particle properties More...
 
class  PDF
 Class for interaction with a PDF set. More...
 
class  PhaseSpacePoint
 Generated point in resummation phase space. More...
 
struct  PhaseSpacePointConfig
 Configuration options for the PhaseSpacePoint class. More...
 
class  ProgressBar
 Class representing (and printing) a progress bar. More...
 
struct  pz_less
 Functor to compare momenta in z direction. More...
 
class  Ranlux64
 Ranlux64 random number generator. More...
 
struct  rapidity_less
 Functor to compare rapidities. More...
 
class  RivetAnalysis
 Class representing a Rivet analysis. More...
 
struct  RNG
 Interface for random number generator. More...
 
struct  RNGConfig
 Settings for random number generator. More...
 
struct  ScaleConfig
 Settings for scale variation. More...
 
class  ScaleFunction
 Class to calculate the scale associated with an event. More...
 
class  ScaleGenerator
 Generate combinations of renormalisation and factorisation scales. More...
 
class  SherpaLHEReader
 Les Houches Event file reader for LHE files created by Sherpa. More...
 
struct  UnclusteredEvent
 
struct  unknown_option
 Exception indicating unknown option. More...
 
class  Unweighter
 Unweight events. More...
 
struct  XSWithError
 Collection of Cross Section with its uncertainty. More...
 

Typedefs

using EventTreatMap = std::map< event_type::EventType, EventTreatment >
 Container to store the treatments for various event types. More...
 
template<typename T >
using optional = boost::optional< T >
 
using Weights = Parameters< double >
 Alias for weight container, e.g. used by the MatrixElement. More...
 
using Colour = std::pair< int, int >
 
using ParticleID = pid::ParticleID
 
template<typename T , std::size_t N, std::size_t... Ns>
using MultiArray = typename detail::ArrayTag< T, N, Ns... >::type
 

Enumerations

enum class  EventTreatment { reweight , keep , discard }
 
enum class  WeightType { weighted , unweighted_resum , partially_unweighted }
 Possible setting for the event weight. More...
 
enum class  FileFormat {
  Les_Houches , HepMC3 , HepMC2 , HDF5 ,
  HepMC =HepMC3
}
 Supported event file formats. More...
 
enum  StatusCode {
  good , discard , empty_jets , failed_reshuffle ,
  failed_resummation_cuts , failed_split , too_much_energy , gluon_in_qqbar ,
  wrong_jets , unspecified
}
 Possible status codes from the event generation. More...
 

Functions

PhaseSpacePointConfig to_PhaseSpacePointConfig (Config const &conf)
 
MatrixElementConfig to_MatrixElementConfig (Config const &conf)
 
EventReweighterConfig to_EventReweighterConfig (Config const &conf)
 
std::ostream & operator<< (std::ostream &os, const CrossSectionAccumulator &xs)
 Print CrossSectionAccumulator to stream. More...
 
std::ostream & operator<< (std::ostream &os, Event const &ev)
 Print Event. More...
 
double shat (Event const &ev)
 Square of the partonic centre-of-mass energy \(\hat{s}\). More...
 
LHEF::HEPEUP to_HEPEUP (Event const &event, LHEF::HEPRUP *)
 Convert an event to a LHEF::HEPEUP. More...
 
std::unique_ptr< EventReadermake_reader (std::string const &filename)
 Factory function for event readers. More...
 
std::unique_ptr< Analysisget_analysis (YAML::Node const &parameters, LHEF::HEPRUP const &heprup)
 Load an analysis. More...
 
std::vector< std::unique_ptr< Analysis > > get_analyses (std::vector< YAML::Node > const &parameters, LHEF::HEPRUP const &heprup)
 Loads multiple analysis, vector version of get_analysis() More...
 
constexpr Helicity flip (Helicity h)
 Flip helicity. More...
 
std::tuple< fastjet::PseudoJet, fastjet::PseudoJet > incoming_momenta (std::vector< Particle > const &outgoing)
 Compute the incoming momentum from momentum conservation. More...
 
auto dot (CLHEP::HepLorentzVector const &pi, CLHEP::HepLorentzVector const &pj)
 "dot" product More...
 
std::complex< double > angle (CLHEP::HepLorentzVector const &pi, CLHEP::HepLorentzVector const &pj)
 "angle" product angle(pi, pj) = <i j> More...
 
std::complex< double > square (CLHEP::HepLorentzVector const &pi, CLHEP::HepLorentzVector const &pj)
 "square" product square(pi, pj) = [i j] More...
 
auto m2 (CLHEP::HepLorentzVector const &h1)
 Invariant mass. More...
 
auto plus (CLHEP::HepLorentzVector const &h1)
 Plus component. More...
 
auto minus (CLHEP::HepLorentzVector const &h1)
 Minus component. More...
 
auto perphat (CLHEP::HepLorentzVector const &h1)
 
std::pair< CLHEP::HepLorentzVector, CLHEP::HepLorentzVector > split_into_lightlike (CLHEP::HepLorentzVector const &P)
 Split a single Lorentz vector into two lightlike Lorentz vectors. More...
 
CLHEP::HepLorentzVector to_HepLorentzVector (fastjet::PseudoJet const &mom)
 
CLHEP::HepLorentzVector to_HepLorentzVector (Particle const &particle)
 
fastjet::PseudoJet to_PseudoJet (CLHEP::HepLorentzVector const &mom)
 
std::unique_ptr< RNGmake_RNG (std::string const &name, optional< std::string > const &seed)
 Factory function for random number generators. More...
 
std::unique_ptr< EventWritermake_format_writer (FileFormat format, std::string const &outfile, LHEF::HEPRUP const &heprup)
 Factory function for event writers. More...
 
std::string to_string (FileFormat f)
 Convert a file format to a string. More...
 
template<class T1 , class T2 >
Parameters< T1 > operator* (Parameters< T1 > a, Parameters< T2 > const &b)
 
template<class T >
Parameters< T > operator* (Parameters< T > a, double b)
 
template<class T >
Parameters< T > operator* (double b, Parameters< T > a)
 
template<class T1 , class T2 >
Parameters< T1 > operator/ (Parameters< T1 > a, Parameters< T2 > const &b)
 
template<class T >
Parameters< T > operator/ (Parameters< T > a, double b)
 
template<class T1 , class T2 >
Parameters< T1 > operator+ (Parameters< T1 > a, Parameters< T2 > const &b)
 
template<class T1 , class T2 >
Parameters< T1 > operator- (Parameters< T1 > a, Parameters< T2 > const &b)
 
std::string to_string (ParameterDescription const &p)
 generate human readable string name More...
 
std::string to_simple_string (ParameterDescription const &p)
 
EventParameters operator* (EventParameters a, double b)
 
EventParameters operator* (double b, EventParameters a)
 
EventParameters operator/ (EventParameters a, double b)
 
std::vector< fastjet::PseudoJet > to_PseudoJet (std::vector< Particle > const &v)
 Convert a vector of Particles to a vector of particle momenta. More...
 
constexpr bool is_parton (Particle const &p)
 Check if a particle is a parton, i.e. quark, antiquark, or gluon. More...
 
constexpr bool is_quark (Particle const &p)
 Check if a particle is a quark. More...
 
constexpr bool is_antiquark (Particle const &p)
 Check if a particle is an anti-quark. More...
 
constexpr bool is_anyquark (Particle const &p)
 Check if a particle is a quark or anit-quark. More...
 
constexpr bool is_lepton (Particle const &p)
 Function to determine if particle is a lepton. More...
 
constexpr bool is_antilepton (Particle const &p)
 Function to determine if particle is an antilepton. More...
 
constexpr bool is_anylepton (Particle const &p)
 Function to determine if particle is an (anti-)lepton. More...
 
constexpr bool is_neutrino (Particle const &p)
 Function to determine if particle is a neutrino. More...
 
constexpr bool is_antineutrino (Particle const &p)
 Function to determine if particle is an antineutrino. More...
 
constexpr bool is_anyneutrino (Particle const &p)
 Function to determine if particle is an (anti-)neutrino. More...
 
constexpr bool is_AWZ_boson (Particle const &particle)
 Check if a particle is a photon, W or Z boson. More...
 
constexpr bool is_AWZH_boson (Particle const &particle)
 Check if a particle is a photon, W, Z, or Higgs boson. More...
 
std::vector< Particlefilter_partons (std::vector< Particle > const &v)
 Extract all partons from a vector of particles. More...
 
std::vector< Particlefilter_AWZH_bosons (std::vector< Particle > const &v)
 Extract all AWZH bosons from a vector of particles. More...
 
ParticleID to_ParticleID (std::string const &name)
 Convert a particle name to the corresponding PDG particle ID. More...
 
constexpr bool is_quark (ParticleID id)
 Function to determine if particle is a quark. More...
 
constexpr bool is_antiquark (ParticleID id)
 Function to determine if particle is an antiquark. More...
 
constexpr bool is_anyquark (ParticleID id)
 Function to determine if particle is an (anti-)quark. More...
 
constexpr bool is_gluon (ParticleID id)
 Function to determine if particle is a gluon. More...
 
constexpr bool is_parton (ParticleID id)
 Function to determine if particle is a parton. More...
 
constexpr bool is_AWZ_boson (ParticleID id)
 function to determine if the particle is a photon, W or Z More...
 
constexpr bool is_AWZH_boson (ParticleID id)
 function to determine if the particle is a photon, W, Z, or Higgs boson More...
 
constexpr bool is_lepton (ParticleID id)
 Function to determine if particle is a lepton. More...
 
constexpr bool is_antilepton (ParticleID id)
 Function to determine if particle is an antilepton. More...
 
constexpr bool is_anylepton (ParticleID id)
 Function to determine if particle is an (anti-)lepton. More...
 
constexpr bool is_neutrino (ParticleID id)
 Function to determine if particle is a neutrino. More...
 
constexpr bool is_antineutrino (ParticleID id)
 Function to determine if particle is an antineutrino. More...
 
constexpr bool is_anyneutrino (ParticleID id)
 Function to determine if particle is an (anti-)neutrino. More...
 
Event::EventData to_EventData (PhaseSpacePoint psp)
 Extract Event::EventData from PhaseSpacePoint. More...
 
std::vector< fastjet::PseudoJet > resummation_jet_momenta (std::vector< fastjet::PseudoJet const * > const &p_born, fastjet::PseudoJet const &qperp)
 Calculate the resummation jet momenta. More...
 
double resummation_jet_weight (std::vector< fastjet::PseudoJet const * > const &p_born, fastjet::PseudoJet const &qperp)
 Calculate additional weight from changing the jet momenta. More...
 
ScaleFunction operator* (double factor, ScaleFunction base_scale)
 Multiply a scale choice by a constant factor. More...
 
ScaleFunction operator* (ScaleFunction factor, ScaleFunction base_scale)
 Multiply a scale choice by a second one. More...
 
ScaleFunction operator/ (ScaleFunction base_scale, double denom)
 Divide a scale choice by a constant factor. More...
 
ScaleFunction operator/ (ScaleFunction base_scale, ScaleFunction denom)
 Divide a scale choice by a second one. More...
 
double H_T (Event const &)
 Calculate \(H_T\) for the input event. More...
 
double max_jet_pt (Event const &)
 The maximum of all (scalar) jet transverse momentum. More...
 
double jet_invariant_mass (Event const &)
 The invariant mass of the sum of all jet momenta. More...
 
double m_j1j2 (Event const &)
 Invariant mass of the two hardest jets. More...
 
std::string to_string (StatusCode s)
 
std::string join (std::string const &)
 
std::string join (std::string const &, std::string const &str)
 
template<typename... Strings>
std::string join (std::string const &delim, std::string const &first, std::string const &second, Strings &&... rest)
 Join strings with a delimiter. More...
 
template<typename T >
std::string type_string (T &&)
 Return the name of the argument's type. More...
 
template<typename... T>
constexpr void ignore (T &&...)
 Eliminate compiler warnings for unused variables. More...
 
constexpr bool nearby_ep (double a, double b, double ep)
 Check whether two doubles are closer than ep > 0 to each other. More...
 
bool nearby_ep (fastjet::PseudoJet const &pa, fastjet::PseudoJet const &pb, double ep)
 Check whether all components of two PseudoJets are closer than ep to each other. More...
 
bool nearby (fastjet::PseudoJet const &pa, fastjet::PseudoJet const &pb, double const norm=1.)
 
Config load_config (std::string const &config_file)
 Load configuration from file. More...
 
template<typename T , typename... YamlNames>
void set_from_yaml (T &setting, YAML::Node const &yaml, YamlNames const &... names)
 Set option using the corresponding YAML entry. More...
 
template<typename T , typename... YamlNames>
void set_from_yaml_if_defined (T &setting, YAML::Node const &yaml, YamlNames const &... names)
 Set option using the corresponding YAML entry, if present. More...
 
JetParameters get_jet_parameters (YAML::Node const &node, std::string const &entry)
 Extract jet parameters from YAML configuration. More...
 
HiggsCouplingSettings get_Higgs_coupling (YAML::Node const &node, std::string const &entry)
 Extract Higgs coupling settings from YAML configuration. More...
 
EWConstants get_ew_parameters (YAML::Node const &node)
 Extract the EW parameters from YAML configuration. More...
 
ScaleConfig to_ScaleConfig (YAML::Node const &yaml)
 Extract scale setting parameters from YAML configuration. More...
 
RNGConfig to_RNGConfig (YAML::Node const &node, std::string const &entry)
 Extract random number generator settings from YAML configuration. More...
 
void assert_all_options_known (YAML::Node const &conf, YAML::Node const &supported)
 Check whether all options in configuration are supported. More...
 

Variables

QCD parameters
constexpr double N_C = 3.
 number of Colours More...
 
constexpr double C_A = N_C
 \(C_A\) More...
 
constexpr double C_F = (N_C*N_C - 1.)/(2.*N_C)
 \(C_F\) More...
 
constexpr double T_F = 0.5
 \(t_f\) More...
 
constexpr double N_F = 5.
 number light flavours More...
 
constexpr double BETA0 = 11./3.*C_A - 4./3.*T_F*N_F
 
Generation Parameters
constexpr double CLAMBDA = 0.2
 Default scale for virtual correction \(\lambda\) cf. eq. (20) in [1]. More...
 
constexpr double CMINPT = 0.2
 
constexpr double DEFAULT_SOFT_PT_REGULATOR = 0.1
 default value for the maximal pt fraction of soft radiation in any tagging jets This cut is needed to regulate an otherwise uncancelled divergence. More...
 
Conventional Parameters
constexpr int COLOUR_OFFSET = 501
 

Detailed Description

Main HEJ 2 Namespace.

Authors
The HEJ collaboration (see AUTHORS for details)
Date
2019-2020

Typedef Documentation

◆ Colour

using HEJ::Colour = typedef std::pair<int,int>

◆ EventTreatMap

Container to store the treatments for various event types.

◆ MultiArray

template<typename T , std::size_t N, std::size_t... Ns>
using HEJ::MultiArray = typedef typename detail::ArrayTag<T, N, Ns...>::type

◆ optional

template<typename T >
using HEJ::optional = typedef boost::optional<T>

◆ ParticleID

◆ Weights

using HEJ::Weights = typedef Parameters<double>

Alias for weight container, e.g. used by the MatrixElement.

Enumeration Type Documentation

◆ EventTreatment

enum class HEJ::EventTreatment
strong

! Possible treatments for fixed-order input events.

The program will decide on how to treat an event based on the value of this enumeration.

Enumerator
reweight 

Perform resummation

keep 

Keep the event

discard 

Discard the event

◆ FileFormat

enum class HEJ::FileFormat
strong

Supported event file formats.

Enumerator
Les_Houches 

Les Houches Output

HepMC3 

HepMC3 Output

HepMC2 

HepMC2 Output

HDF5 

HDF5 Output

HepMC 

HepMC3 Output

◆ StatusCode

Possible status codes from the event generation.

Enumerator
good 
discard 
empty_jets 
failed_reshuffle 
failed_resummation_cuts 
failed_split 
too_much_energy 
gluon_in_qqbar 
wrong_jets 
unspecified 

◆ WeightType

enum class HEJ::WeightType
strong

Possible setting for the event weight.

Enumerator
weighted 

weighted events

unweighted_resum 

unweighted only resummation part

partially_unweighted 

mixed weighted and unweighted

Function Documentation

◆ angle()

std::complex< double > HEJ::angle ( CLHEP::HepLorentzVector const &  pi,
CLHEP::HepLorentzVector const &  pj 
)

"angle" product angle(pi, pj) = <i j>

◆ assert_all_options_known()

void HEJ::assert_all_options_known ( YAML::Node const &  conf,
YAML::Node const &  supported 
)

Check whether all options in configuration are supported.

Parameters
confConfiguration to be checked
supportedTree of supported options

If conf contains an entry that does not appear in supported an unknown_option exception is thrown. Sub-entries of "analysis" (if present) are not checked.

See also
unknown_option

◆ dot()

auto HEJ::dot ( CLHEP::HepLorentzVector const &  pi,
CLHEP::HepLorentzVector const &  pj 
)
inline

"dot" product

◆ filter_AWZH_bosons()

std::vector< Particle > HEJ::filter_AWZH_bosons ( std::vector< Particle > const &  v)
inline

Extract all AWZH bosons from a vector of particles.

◆ filter_partons()

std::vector< Particle > HEJ::filter_partons ( std::vector< Particle > const &  v)
inline

Extract all partons from a vector of particles.

◆ flip()

constexpr Helicity HEJ::flip ( Helicity  h)
inlineconstexpr

Flip helicity.

Parameters
hHelicity to flip
Returns
plus iff the input helicity was minus and vice versa

◆ get_analyses()

std::vector< std::unique_ptr< Analysis > > HEJ::get_analyses ( std::vector< YAML::Node > const &  parameters,
LHEF::HEPRUP const &  heprup 
)

Loads multiple analysis, vector version of get_analysis()

Parameters
parametersVector of Analysis parameters
heprupGeneral run informations
Returns
Vector of pointers to an Analysis instance

◆ get_analysis()

std::unique_ptr< Analysis > HEJ::get_analysis ( YAML::Node const &  parameters,
LHEF::HEPRUP const &  heprup 
)

Load an analysis.

Parameters
parametersAnalysis parameters
heprupGeneral run informations
Returns
A pointer to an Analysis instance

If parameters["plugin"] exists, an analysis (deriving from the Analysis class) will be loaded from the library parameters["plugin"]. Otherwise, if parameters["rivet"] exists, the corresponding RivetAnalysis will be loaded. If none of these parameters are specified, a pointer to the default EmptyAnalysis is returned.

◆ get_ew_parameters()

EWConstants HEJ::get_ew_parameters ( YAML::Node const &  node)

Extract the EW parameters from YAML configuration.

◆ get_Higgs_coupling()

HiggsCouplingSettings HEJ::get_Higgs_coupling ( YAML::Node const &  node,
std::string const &  entry 
)

Extract Higgs coupling settings from YAML configuration.

◆ get_jet_parameters()

JetParameters HEJ::get_jet_parameters ( YAML::Node const &  node,
std::string const &  entry 
)

Extract jet parameters from YAML configuration.

◆ H_T()

double HEJ::H_T ( Event const &  )

Calculate \(H_T\) for the input event.

\(H_T\) is the sum of the (scalar) transverse momenta of all final-state particles

◆ ignore()

template<typename... T>
constexpr void HEJ::ignore ( T &&  ...)
constexpr

Eliminate compiler warnings for unused variables.

◆ incoming_momenta()

std::tuple< fastjet::PseudoJet, fastjet::PseudoJet > HEJ::incoming_momenta ( std::vector< Particle > const &  outgoing)

Compute the incoming momentum from momentum conservation.

Parameters
outgoingOutgoing particles

◆ is_antilepton() [1/2]

constexpr bool HEJ::is_antilepton ( Particle const &  p)
inlineconstexpr

Function to determine if particle is an antilepton.

Parameters
pthe particle
Returns
true if the particle is an antilepton, false otherwise

◆ is_antilepton() [2/2]

constexpr bool HEJ::is_antilepton ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an antilepton.

Parameters
idPDG ID of particle
Returns
true if the particle is an antilepton, false otherwise

◆ is_antineutrino() [1/2]

constexpr bool HEJ::is_antineutrino ( Particle const &  p)
inlineconstexpr

Function to determine if particle is an antineutrino.

Parameters
pthe particle
Returns
true if the particle is an antineutrino, false otherwise

◆ is_antineutrino() [2/2]

constexpr bool HEJ::is_antineutrino ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an antineutrino.

Parameters
idPDG ID of particle
Returns
true if the particle is an antineutrino, false otherwise

◆ is_antiquark() [1/2]

constexpr bool HEJ::is_antiquark ( Particle const &  p)
inlineconstexpr

Check if a particle is an anti-quark.

◆ is_antiquark() [2/2]

constexpr bool HEJ::is_antiquark ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an antiquark.

Parameters
idPDG ID of particle
Returns
true if the particle is an antiquark, false otherwise

◆ is_anylepton() [1/2]

constexpr bool HEJ::is_anylepton ( Particle const &  p)
inlineconstexpr

Function to determine if particle is an (anti-)lepton.

Parameters
pthe particle
Returns
true if the particle is a lepton or antilepton, false otherwise

◆ is_anylepton() [2/2]

constexpr bool HEJ::is_anylepton ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an (anti-)lepton.

Parameters
idPDG ID of particle
Returns
true if the particle is a lepton or antilepton, false otherwise

◆ is_anyneutrino() [1/2]

constexpr bool HEJ::is_anyneutrino ( Particle const &  p)
inlineconstexpr

Function to determine if particle is an (anti-)neutrino.

Parameters
pthe particle
Returns
true if the particle is a neutrino or antineutrino, false otherwise

◆ is_anyneutrino() [2/2]

constexpr bool HEJ::is_anyneutrino ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an (anti-)neutrino.

Parameters
idPDG ID of particle
Returns
true if the particle is a neutrino or antineutrino, false otherwise

◆ is_anyquark() [1/2]

constexpr bool HEJ::is_anyquark ( Particle const &  p)
inlineconstexpr

Check if a particle is a quark or anit-quark.

◆ is_anyquark() [2/2]

constexpr bool HEJ::is_anyquark ( ParticleID  id)
inlineconstexpr

Function to determine if particle is an (anti-)quark.

Parameters
idPDG ID of particle
Returns
true if the particle is a quark or antiquark, false otherwise

◆ is_AWZ_boson() [1/2]

constexpr bool HEJ::is_AWZ_boson ( Particle const &  particle)
inlineconstexpr

Check if a particle is a photon, W or Z boson.

◆ is_AWZ_boson() [2/2]

constexpr bool HEJ::is_AWZ_boson ( ParticleID  id)
inlineconstexpr

function to determine if the particle is a photon, W or Z

Parameters
idPDG ID of particle
Returns
true if the partice is an A,W,Z, or H, false otherwise

◆ is_AWZH_boson() [1/2]

constexpr bool HEJ::is_AWZH_boson ( Particle const &  particle)
inlineconstexpr

Check if a particle is a photon, W, Z, or Higgs boson.

◆ is_AWZH_boson() [2/2]

constexpr bool HEJ::is_AWZH_boson ( ParticleID  id)
inlineconstexpr

function to determine if the particle is a photon, W, Z, or Higgs boson

Parameters
idPDG ID of particle
Returns
true if the partice is an A,W,Z, or H, false otherwise

◆ is_gluon()

constexpr bool HEJ::is_gluon ( ParticleID  id)
inlineconstexpr

Function to determine if particle is a gluon.

Parameters
idPDG ID of particle
Returns
true if the particle is a gluon, false otherwise

◆ is_lepton() [1/2]

constexpr bool HEJ::is_lepton ( Particle const &  p)
inlineconstexpr

Function to determine if particle is a lepton.

Parameters
pthe particle
Returns
true if the particle is a lepton, false otherwise

◆ is_lepton() [2/2]

constexpr bool HEJ::is_lepton ( ParticleID  id)
inlineconstexpr

Function to determine if particle is a lepton.

Parameters
idPDG ID of particle
Returns
true if the particle is a lepton, false otherwise

◆ is_neutrino() [1/2]

constexpr bool HEJ::is_neutrino ( Particle const &  p)
inlineconstexpr

Function to determine if particle is a neutrino.

Parameters
pthe particle
Returns
true if the particle is a neutrino, false otherwise

◆ is_neutrino() [2/2]

constexpr bool HEJ::is_neutrino ( ParticleID  id)
inlineconstexpr

Function to determine if particle is a neutrino.

Parameters
idPDG ID of particle
Returns
true if the particle is a neutrino, false otherwise

◆ is_parton() [1/2]

constexpr bool HEJ::is_parton ( Particle const &  p)
inlineconstexpr

Check if a particle is a parton, i.e. quark, antiquark, or gluon.

◆ is_parton() [2/2]

constexpr bool HEJ::is_parton ( ParticleID  id)
inlineconstexpr

Function to determine if particle is a parton.

Parameters
idPDG ID of particle
Returns
true if the particle is a parton, false otherwise

◆ is_quark() [1/2]

constexpr bool HEJ::is_quark ( Particle const &  p)
inlineconstexpr

Check if a particle is a quark.

◆ is_quark() [2/2]

constexpr bool HEJ::is_quark ( ParticleID  id)
inlineconstexpr

Function to determine if particle is a quark.

Parameters
idPDG ID of particle
Returns
true if the particle is a quark, false otherwise

◆ jet_invariant_mass()

double HEJ::jet_invariant_mass ( Event const &  )

The invariant mass of the sum of all jet momenta.

◆ join() [1/3]

std::string HEJ::join ( std::string const &  )
inline

◆ join() [2/3]

std::string HEJ::join ( std::string const &  ,
std::string const &  str 
)
inline

◆ join() [3/3]

template<typename... Strings>
std::string HEJ::join ( std::string const &  delim,
std::string const &  first,
std::string const &  second,
Strings &&...  rest 
)

Join strings with a delimiter.

Parameters
delimDelimiter to be put between consecutive strings
firstFirst string
secondSecond string
restRemaining strings

◆ load_config()

Config HEJ::load_config ( std::string const &  config_file)

Load configuration from file.

Parameters
config_fileName of the YAML configuration file
Returns
The HEJ 2 configuration

◆ m2()

auto HEJ::m2 ( CLHEP::HepLorentzVector const &  h1)
inline

Invariant mass.

◆ m_j1j2()

double HEJ::m_j1j2 ( Event const &  )

Invariant mass of the two hardest jets.

◆ make_format_writer()

std::unique_ptr< EventWriter > HEJ::make_format_writer ( FileFormat  format,
std::string const &  outfile,
LHEF::HEPRUP const &  heprup 
)

Factory function for event writers.

Parameters
formatThe format of the output file
outfileThe name of the output file
heprupGeneral process information
Returns
A pointer to an instance of an EventWriter for the desired format

◆ make_reader()

std::unique_ptr< EventReader > HEJ::make_reader ( std::string const &  filename)

Factory function for event readers.

Parameters
filenameThe name of the input file
Returns
A pointer to an instance of an EventReader for the input file

◆ make_RNG()

std::unique_ptr< RNG > HEJ::make_RNG ( std::string const &  name,
optional< std::string > const &  seed 
)

Factory function for random number generators.

Parameters
nameName of the random number generator
seedOptional seed
Returns
A pointer to an instance of a random number generator

At present, name should be one of "ranlux64" or "mixmax" (case insensitive). The interpretation of the seed depends on the random number generator. For ranlux64, it is the name of a seed file. For mixmax it should be a string convertible to a long integer.

◆ max_jet_pt()

double HEJ::max_jet_pt ( Event const &  )

The maximum of all (scalar) jet transverse momentum.

◆ minus()

auto HEJ::minus ( CLHEP::HepLorentzVector const &  h1)
inline

Minus component.

◆ nearby()

bool HEJ::nearby ( fastjet::PseudoJet const &  pa,
fastjet::PseudoJet const &  pb,
double const  norm = 1. 
)
inline

◆ nearby_ep() [1/2]

constexpr bool HEJ::nearby_ep ( double  a,
double  b,
double  ep 
)
inlineconstexpr

Check whether two doubles are closer than ep > 0 to each other.

◆ nearby_ep() [2/2]

bool HEJ::nearby_ep ( fastjet::PseudoJet const &  pa,
fastjet::PseudoJet const &  pb,
double  ep 
)
inline

Check whether all components of two PseudoJets are closer than ep to each other.

◆ operator*() [1/7]

EventParameters HEJ::operator* ( double  b,
EventParameters  a 
)
inline

◆ operator*() [2/7]

template<class T >
Parameters< T > HEJ::operator* ( double  b,
Parameters< T >  a 
)
inline

◆ operator*() [3/7]

ScaleFunction HEJ::operator* ( double  factor,
ScaleFunction  base_scale 
)

Multiply a scale choice by a constant factor.

For example, multiplying 0.5 and a scale function for H_T will result in a scale function for H_T/2.

◆ operator*() [4/7]

EventParameters HEJ::operator* ( EventParameters  a,
double  b 
)
inline

◆ operator*() [5/7]

template<class T >
Parameters< T > HEJ::operator* ( Parameters< T >  a,
double  b 
)
inline

◆ operator*() [6/7]

template<class T1 , class T2 >
Parameters< T1 > HEJ::operator* ( Parameters< T1 >  a,
Parameters< T2 > const &  b 
)
inline

◆ operator*() [7/7]

ScaleFunction HEJ::operator* ( ScaleFunction  factor,
ScaleFunction  base_scale 
)

Multiply a scale choice by a second one.

For example, multiplying H_T and m_j1j2 will result in a scale function for H_T*m_j1j2.

◆ operator+()

template<class T1 , class T2 >
Parameters< T1 > HEJ::operator+ ( Parameters< T1 >  a,
Parameters< T2 > const &  b 
)
inline

◆ operator-()

template<class T1 , class T2 >
Parameters< T1 > HEJ::operator- ( Parameters< T1 >  a,
Parameters< T2 > const &  b 
)
inline

◆ operator/() [1/5]

EventParameters HEJ::operator/ ( EventParameters  a,
double  b 
)
inline

◆ operator/() [2/5]

template<class T >
Parameters< T > HEJ::operator/ ( Parameters< T >  a,
double  b 
)
inline

◆ operator/() [3/5]

template<class T1 , class T2 >
Parameters< T1 > HEJ::operator/ ( Parameters< T1 >  a,
Parameters< T2 > const &  b 
)
inline

◆ operator/() [4/5]

ScaleFunction HEJ::operator/ ( ScaleFunction  base_scale,
double  denom 
)

Divide a scale choice by a constant factor.

For example, dividing a scale function for H_T by 2 will result in a scale function for H_T/2.

◆ operator/() [5/5]

ScaleFunction HEJ::operator/ ( ScaleFunction  base_scale,
ScaleFunction  denom 
)

Divide a scale choice by a second one.

For example, dividing a scale function for H_T by m_j1j2 will result in a scale function for H_T/m_j1j2.

◆ operator<<() [1/2]

std::ostream & HEJ::operator<< ( std::ostream &  os,
const CrossSectionAccumulator xs 
)

Print CrossSectionAccumulator to stream.

◆ operator<<() [2/2]

std::ostream & HEJ::operator<< ( std::ostream &  os,
Event const &  ev 
)

Print Event.

◆ perphat()

auto HEJ::perphat ( CLHEP::HepLorentzVector const &  h1)
inline

◆ plus()

auto HEJ::plus ( CLHEP::HepLorentzVector const &  h1)
inline

Plus component.

◆ resummation_jet_momenta()

std::vector< fastjet::PseudoJet > HEJ::resummation_jet_momenta ( std::vector< fastjet::PseudoJet const * > const &  p_born,
fastjet::PseudoJet const &  qperp 
)

Calculate the resummation jet momenta.

Parameters
p_bornborn Jet Momenta
qperpSum of non-jet Parton Transverse Momenta
Returns
Resummation Jet Momenta

◆ resummation_jet_weight()

double HEJ::resummation_jet_weight ( std::vector< fastjet::PseudoJet const * > const &  p_born,
fastjet::PseudoJet const &  qperp 
)

Calculate additional weight from changing the jet momenta.

Parameters
p_bornborn Jet Momenta
qperpSum of non-jet Parton Transverse Momenta

Computes the Jacobian for changing the original delta functions expressed in terms of jet momenta to delta functions of the parton momenta in the resummation phase space

◆ set_from_yaml()

template<typename T , typename... YamlNames>
void HEJ::set_from_yaml ( T &  setting,
YAML::Node const &  yaml,
YamlNames const &...  names 
)

Set option using the corresponding YAML entry.

Parameters
settingOption variable to be set
yamlRoot of the YAML configuration
namesName of the entry

If the entry does not exist or has the wrong type or format an exception is thrown.

For example

set_from_yaml(foobar, yaml, "foo", "bar")
void set_from_yaml(T &setting, YAML::Node const &yaml, YamlNames const &... names)
Set option using the corresponding YAML entry.
Definition: YAMLreader.hh:202

is equivalent to

foobar = yaml["foo"]["bar"].as<decltype(foobar)>()

with improved diagnostics on errors.

See also
set_from_yaml_if_defined

◆ set_from_yaml_if_defined()

template<typename T , typename... YamlNames>
void HEJ::set_from_yaml_if_defined ( T &  setting,
YAML::Node const &  yaml,
YamlNames const &...  names 
)

Set option using the corresponding YAML entry, if present.

Parameters
settingOption variable to be set
yamlRoot of the YAML configuration
namesName of the entry

This function works similar to set_from_yaml, but does not throw any exception if the requested YAML entry does not exist.

See also
set_from_yaml

◆ shat()

double HEJ::shat ( Event const &  ev)

Square of the partonic centre-of-mass energy \(\hat{s}\).

◆ split_into_lightlike()

std::pair< CLHEP::HepLorentzVector, CLHEP::HepLorentzVector > HEJ::split_into_lightlike ( CLHEP::HepLorentzVector const &  P)

Split a single Lorentz vector into two lightlike Lorentz vectors.

Parameters
PLorentz vector to be split
Returns
A pair (p, q) of Lorentz vectors with P = p + q and p^2 = q^2 = 0

P.perp() has to be positive.

p.e() is guaranteed to be positive. In addition, if either of P.plus() or P.minus() is positive, q.e() has the same sign as P.m2()

◆ square()

std::complex< double > HEJ::square ( CLHEP::HepLorentzVector const &  pi,
CLHEP::HepLorentzVector const &  pj 
)

"square" product square(pi, pj) = [i j]

◆ to_EventData()

Event::EventData HEJ::to_EventData ( PhaseSpacePoint  psp)

◆ to_EventReweighterConfig()

EventReweighterConfig HEJ::to_EventReweighterConfig ( Config const &  conf)
inline

◆ to_HEPEUP()

LHEF::HEPEUP HEJ::to_HEPEUP ( Event const &  event,
LHEF::HEPRUP *   
)

Convert an event to a LHEF::HEPEUP.

◆ to_HepLorentzVector() [1/2]

CLHEP::HepLorentzVector HEJ::to_HepLorentzVector ( fastjet::PseudoJet const &  mom)
inline

◆ to_HepLorentzVector() [2/2]

CLHEP::HepLorentzVector HEJ::to_HepLorentzVector ( Particle const &  particle)
inline

◆ to_MatrixElementConfig()

MatrixElementConfig HEJ::to_MatrixElementConfig ( Config const &  conf)
inline

◆ to_ParticleID()

ParticleID HEJ::to_ParticleID ( std::string const &  name)

Convert a particle name to the corresponding PDG particle ID.

◆ to_PhaseSpacePointConfig()

PhaseSpacePointConfig HEJ::to_PhaseSpacePointConfig ( Config const &  conf)
inline

! Extract PhaseSpacePointConfig from Config

◆ to_PseudoJet() [1/2]

fastjet::PseudoJet HEJ::to_PseudoJet ( CLHEP::HepLorentzVector const &  mom)
inline

◆ to_PseudoJet() [2/2]

std::vector< fastjet::PseudoJet > HEJ::to_PseudoJet ( std::vector< Particle > const &  v)
inline

Convert a vector of Particles to a vector of particle momenta.

◆ to_RNGConfig()

RNGConfig HEJ::to_RNGConfig ( YAML::Node const &  node,
std::string const &  entry 
)

Extract random number generator settings from YAML configuration.

◆ to_ScaleConfig()

ScaleConfig HEJ::to_ScaleConfig ( YAML::Node const &  yaml)

Extract scale setting parameters from YAML configuration.

◆ to_simple_string()

std::string HEJ::to_simple_string ( ParameterDescription const &  p)

generate simplified string, intended for easy parsing Format: Scale_SCALENAME_MuRxx_MuFyy

◆ to_string() [1/3]

std::string HEJ::to_string ( FileFormat  f)
inline

Convert a file format to a string.

◆ to_string() [2/3]

std::string HEJ::to_string ( ParameterDescription const &  p)

generate human readable string name

◆ to_string() [3/3]

std::string HEJ::to_string ( StatusCode  s)
inline

Get name of StatusCode

Todo:
better names

◆ type_string()

template<typename T >
std::string HEJ::type_string ( T &&  )

Return the name of the argument's type.

Variable Documentation

◆ BETA0

constexpr double HEJ::BETA0 = 11./3.*C_A - 4./3.*T_F*N_F
constexpr

\(\beta_0\)

◆ C_A

constexpr double HEJ::C_A = N_C
constexpr

\(C_A\)

◆ C_F

constexpr double HEJ::C_F = (N_C*N_C - 1.)/(2.*N_C)
constexpr

\(C_F\)

◆ CLAMBDA

constexpr double HEJ::CLAMBDA = 0.2
constexpr

Default scale for virtual correction \(\lambda\) cf. eq. (20) in [1].

◆ CMINPT

constexpr double HEJ::CMINPT = 0.2
constexpr

minimal \(p_t\) of all partons

◆ COLOUR_OFFSET

constexpr int HEJ::COLOUR_OFFSET = 501
constexpr

Value of first colour for colour dressing, according to LHE convention [3]

◆ DEFAULT_SOFT_PT_REGULATOR

constexpr double HEJ::DEFAULT_SOFT_PT_REGULATOR = 0.1
constexpr

default value for the maximal pt fraction of soft radiation in any tagging jets This cut is needed to regulate an otherwise uncancelled divergence.

◆ N_C

constexpr double HEJ::N_C = 3.
constexpr

number of Colours

◆ N_F

constexpr double HEJ::N_F = 5.
constexpr

number light flavours

◆ T_F

constexpr double HEJ::T_F = 0.5
constexpr

\(t_f\)