hej is hosted by Hepforge, IPPP Durham
HEJ 2.2.2
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
HEJ::Event Class Reference

An event with clustered jets. More...

#include <Event.hh>

Classes

class  EventData
 Class to store general Event setup, i.e. Phase space and weights. More...
 

Public Types

using ConstPartonIterator = boost::filter_iterator< bool(*)(Particle const &), std::vector< Particle >::const_iterator >
 Iterator over partons. More...
 
using ConstReversePartonIterator = std::reverse_iterator< ConstPartonIterator >
 Reverse Iterator over partons. More...
 

Public Member Functions

 Event ()=delete
 No default Constructor. More...
 
std::vector< int > particle_jet_indices (std::vector< fastjet::PseudoJet > const &jets) const
 Indices of the jets the outgoing partons belong to. More...
 
std::vector< int > particle_jet_indices () const
 particle_jet_indices() of the Event jets() More...
 
fastjet::JetDefinition const & jet_def () const
 Jet definition used for clustering. More...
 
double min_jet_pt () const
 Minimum jet transverse momentum. More...
 
event_type::EventType type () const
 Event type. More...
 
bool generate_colours (RNG &)
 Give colours to each particle. More...
 
bool is_leading_colour () const
 Check that current colours are leading in the high energy limit. More...
 
bool valid_hej_state (double soft_pt_regulator=DEFAULT_SOFT_PT_REGULATOR) const
 Check if given event could have been produced by HEJ. More...
 
bool valid_incoming () const
 Check that the incoming momenta are valid. More...
 
Particle Access
std::array< Particle, 2 > const & incoming () const
 Incoming particles. More...
 
std::vector< Particle > const & outgoing () const
 Outgoing particles. More...
 
ConstPartonIterator begin_partons () const
 Iterator to the first outgoing parton. More...
 
ConstPartonIterator cbegin_partons () const
 Iterator to the first outgoing parton. More...
 
ConstPartonIterator end_partons () const
 Iterator to the end of the outgoing partons. More...
 
ConstPartonIterator cend_partons () const
 Iterator to the end of the outgoing partons. More...
 
ConstReversePartonIterator rbegin_partons () const
 Reverse Iterator to the first outgoing parton. More...
 
ConstReversePartonIterator crbegin_partons () const
 Reverse Iterator to the first outgoing parton. More...
 
ConstReversePartonIterator rend_partons () const
 Reverse Iterator to the first outgoing parton. More...
 
ConstReversePartonIterator crend_partons () const
 Reverse Iterator to the first outgoing parton. More...
 
std::unordered_map< std::size_t, std::vector< Particle > > const & decays () const
 Particle decays. More...
 
std::vector< fastjet::PseudoJet > const & jets () const
 The jets formed by the outgoing partons, sorted in rapidity. More...
 
Weight variations
Parameters< EventParameters > const & parameters () const
 All chosen parameter, i.e. scale choices (const version) More...
 
Parameters< EventParameters > & parameters ()
 All chosen parameter, i.e. scale choices. More...
 
EventParameters const & central () const
 Central parameter choice (const version) More...
 
EventParameterscentral ()
 Central parameter choice. More...
 
std::vector< EventParameters > const & variations () const
 Parameter (scale) variations (const version) More...
 
std::vector< EventParameters > & variations ()
 Parameter (scale) variations. More...
 
EventParameters const & variations (std::size_t i) const
 Parameter (scale) variation (const version) More...
 
EventParametersvariations (std::size_t i)
 Parameter (scale) variation. More...
 

Detailed Description

An event with clustered jets.

This is the main HEJ 2 event class. It contains kinematic information including jet clustering, parameter (e.g. scale) settings and the event weight.

Member Typedef Documentation

◆ ConstPartonIterator

using HEJ::Event::ConstPartonIterator = boost::filter_iterator< bool (*)(Particle const &), std::vector<Particle>::const_iterator >

Iterator over partons.

◆ ConstReversePartonIterator

Reverse Iterator over partons.

Constructor & Destructor Documentation

◆ Event()

HEJ::Event::Event ( )
delete

No default Constructor.

Member Function Documentation

◆ begin_partons()

ConstPartonIterator HEJ::Event::begin_partons ( ) const

Iterator to the first outgoing parton.

◆ cbegin_partons()

ConstPartonIterator HEJ::Event::cbegin_partons ( ) const

Iterator to the first outgoing parton.

◆ cend_partons()

ConstPartonIterator HEJ::Event::cend_partons ( ) const

Iterator to the end of the outgoing partons.

◆ central() [1/2]

EventParameters & HEJ::Event::central ( )
inline

Central parameter choice.

◆ central() [2/2]

EventParameters const & HEJ::Event::central ( ) const
inline

Central parameter choice (const version)

◆ crbegin_partons()

ConstReversePartonIterator HEJ::Event::crbegin_partons ( ) const

Reverse Iterator to the first outgoing parton.

◆ crend_partons()

ConstReversePartonIterator HEJ::Event::crend_partons ( ) const

Reverse Iterator to the first outgoing parton.

◆ decays()

std::unordered_map< std::size_t, std::vector< Particle > > const & HEJ::Event::decays ( ) const
inline

Particle decays.

The key in the returned map corresponds to the index in the vector returned by outgoing()

◆ end_partons()

ConstPartonIterator HEJ::Event::end_partons ( ) const

Iterator to the end of the outgoing partons.

◆ generate_colours()

bool HEJ::Event::generate_colours ( RNG )

Give colours to each particle.

Returns
true if new colours are generated, i.e. same as is_resummable()

Colour ordering is done according to leading colour in the MRK limit, see [2]. This only affects HEJ configurations, all other EventTypes will be ignored.

Note
This overwrites all previously set colours.

◆ incoming()

std::array< Particle, 2 > const & HEJ::Event::incoming ( ) const
inline

Incoming particles.

◆ is_leading_colour()

bool HEJ::Event::is_leading_colour ( ) const

Check that current colours are leading in the high energy limit.

Checks that the colour configuration can be split up in multiple, rapidity ordered, non-overlapping ladders. Such configurations are leading in the MRK limit, see [2]

Note
This is not to be confused with is_resummable(), however for all resummable states it is possible to create a leading colour configuration, see generate_colours()

◆ jet_def()

fastjet::JetDefinition const & HEJ::Event::jet_def ( ) const
inline

Jet definition used for clustering.

◆ jets()

std::vector< fastjet::PseudoJet > const & HEJ::Event::jets ( ) const
inline

The jets formed by the outgoing partons, sorted in rapidity.

◆ min_jet_pt()

double HEJ::Event::min_jet_pt ( ) const
inline

Minimum jet transverse momentum.

◆ outgoing()

std::vector< Particle > const & HEJ::Event::outgoing ( ) const
inline

Outgoing particles.

◆ parameters() [1/2]

Parameters< EventParameters > & HEJ::Event::parameters ( )
inline

All chosen parameter, i.e. scale choices.

◆ parameters() [2/2]

Parameters< EventParameters > const & HEJ::Event::parameters ( ) const
inline

All chosen parameter, i.e. scale choices (const version)

◆ particle_jet_indices() [1/2]

std::vector< int > HEJ::Event::particle_jet_indices ( ) const
inline

◆ particle_jet_indices() [2/2]

std::vector< int > HEJ::Event::particle_jet_indices ( std::vector< fastjet::PseudoJet > const &  jets) const
inline

Indices of the jets the outgoing partons belong to.

Parameters
jetsJets to be tested
Returns
A vector containing, for each outgoing parton, the index in the vector of jets the considered parton belongs to. If the parton is not inside any of the passed jets, the corresponding index is set to -1.

◆ rbegin_partons()

ConstReversePartonIterator HEJ::Event::rbegin_partons ( ) const

Reverse Iterator to the first outgoing parton.

◆ rend_partons()

ConstReversePartonIterator HEJ::Event::rend_partons ( ) const

Reverse Iterator to the first outgoing parton.

◆ type()

event_type::EventType HEJ::Event::type ( ) const
inline

Event type.

◆ valid_hej_state()

bool HEJ::Event::valid_hej_state ( double  soft_pt_regulator = DEFAULT_SOFT_PT_REGULATOR) const

Check if given event could have been produced by HEJ.

A HEJ state has to fulfil:

  1. type() has to be resummable
  2. Soft radiation in the tagging jets contributes at most to soft_pt_regulator of the total jet \( p_\perp \)
  3. Partons related to subleading configurations (uno gluon, qqbar) must be in separate jets.
Note
This is true for any resummed stated produced by the EventReweighter or any resummable Leading Order state.
Parameters
soft_pt_regulatorMaximum transverse momentum fraction from soft radiation in tagging jets
Returns
True if this state could have been produced by HEJ

◆ valid_incoming()

bool HEJ::Event::valid_incoming ( ) const

Check that the incoming momenta are valid.

Checks that the incoming parton momenta are on-shell and have vanishing transverse components.

◆ variations() [1/4]

std::vector< EventParameters > & HEJ::Event::variations ( )
inline

Parameter (scale) variations.

◆ variations() [2/4]

std::vector< EventParameters > const & HEJ::Event::variations ( ) const
inline

Parameter (scale) variations (const version)

◆ variations() [3/4]

EventParameters & HEJ::Event::variations ( std::size_t  i)
inline

Parameter (scale) variation.

Parameters
iIndex of the requested variation

◆ variations() [4/4]

EventParameters const & HEJ::Event::variations ( std::size_t  i) const
inline

Parameter (scale) variation (const version)

Parameters
iIndex of the requested variation

The documentation for this class was generated from the following file: