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... | |
Event (UnclusteredEvent const &ev, fastjet::JetDefinition const &jet_def, double min_jet_pt) | |
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, double min_pt=0.) const |
Check if given event could have been produced by HEJ. 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... | |
EventParameters & | central () |
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... | |
EventParameters & | variations (std::size_t i) |
Parameter (scale) variation. More... | |
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.
using HEJ::Event::ConstPartonIterator = boost::filter_iterator< bool (*)(Particle const &), std::vector<Particle>::const_iterator > |
Iterator over partons.
using HEJ::Event::ConstReversePartonIterator = std::reverse_iterator< ConstPartonIterator> |
Reverse Iterator over partons.
|
delete |
No default Constructor.
HEJ::Event::Event | ( | UnclusteredEvent const & | ev, |
fastjet::JetDefinition const & | jet_def, | ||
double | min_jet_pt | ||
) |
Event Constructor adding jet clustering to an unclustered event
ConstPartonIterator HEJ::Event::begin_partons | ( | ) | const |
Iterator to the first outgoing parton.
ConstPartonIterator HEJ::Event::cbegin_partons | ( | ) | const |
Iterator to the first outgoing parton.
ConstPartonIterator HEJ::Event::cend_partons | ( | ) | const |
Iterator to the end of the outgoing partons.
|
inline |
Central parameter choice.
|
inline |
Central parameter choice (const version)
ConstReversePartonIterator HEJ::Event::crbegin_partons | ( | ) | const |
Reverse Iterator to the first outgoing parton.
ConstReversePartonIterator HEJ::Event::crend_partons | ( | ) | const |
Reverse Iterator to the first outgoing parton.
|
inline |
Particle decays.
The key in the returned map corresponds to the index in the vector returned by outgoing()
ConstPartonIterator HEJ::Event::end_partons | ( | ) | const |
Iterator to the end of the outgoing partons.
bool HEJ::Event::generate_colours | ( | RNG & | ) |
Give colours to each particle.
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.
|
inline |
Incoming particles.
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]
|
inline |
Jet definition used for clustering.
|
inline |
The jets formed by the outgoing partons, sorted in rapidity.
|
inline |
Minimum jet transverse momentum.
|
inline |
Outgoing particles.
|
inline |
All chosen parameter, i.e. scale choices.
|
inline |
All chosen parameter, i.e. scale choices (const version)
|
inline |
particle_jet_indices() of the Event jets()
|
inline |
Indices of the jets the outgoing partons belong to.
jets | Jets to be tested |
ConstReversePartonIterator HEJ::Event::rbegin_partons | ( | ) | const |
Reverse Iterator to the first outgoing parton.
ConstReversePartonIterator HEJ::Event::rend_partons | ( | ) | const |
Reverse Iterator to the first outgoing parton.
|
inline |
Event type.
bool HEJ::Event::valid_hej_state | ( | double | soft_pt_regulator = DEFAULT_SOFT_PT_REGULATOR , |
double | min_pt = 0. |
||
) | const |
Check if given event could have been produced by HEJ.
A HEJ state has to fulfil:
soft_pt_regulator
of the total jet \( p_\perp \)soft_pt_regulator | Maximum transverse momentum fraction from soft radiation in tagging jets |
min_pt | Absolute minimal \( p_\perp \), deprecated use soft_pt_regulator instead |
|
inline |
Parameter (scale) variations.
|
inline |
Parameter (scale) variations (const version)
|
inline |
Parameter (scale) variation.
i | Index of the requested variation |
|
inline |
Parameter (scale) variation (const version)
i | Index of the requested variation |