Loading [MathJax]/extensions/tex2jax.js
hej is hosted by Hepforge, IPPP Durham
HEJ  2.3.0
High energy resummation for hadron colliders
Event.hh
Go to the documentation of this file.
1 
8 #pragma once
9 
10 #include <array>
11 #include <cstddef>
12 #include <iosfwd>
13 #include <iterator>
14 #include <string>
15 #include <unordered_map>
16 #include <utility>
17 #include <vector>
18 
19 #include "boost/iterator/filter_iterator.hpp"
20 
21 #include "fastjet/ClusterSequence.hh"
22 #include "fastjet/PseudoJet.hh"
23 
24 #include "HEJ/Constants.hh"
25 #include "HEJ/Parameters.hh"
26 #include "HEJ/Particle.hh"
27 #include "HEJ/event_types.hh"
28 
29 namespace LHEF {
30  class HEPEUP;
31  class HEPRUP;
32 }
33 
34 namespace fastjet {
35  class JetDefinition;
36 }
37 
38 namespace HEJ {
39  class EWConstants;
40  struct RNG;
41 
43  size_t implemented_types(std::vector<Particle> const & bosons);
44 
51  class Event {
52  public:
53  class EventData;
54 
56  using ConstPartonIterator = boost::filter_iterator<
57  bool (*)(Particle const &),
58  std::vector<Particle>::const_iterator
59  >;
61  using ConstReversePartonIterator = std::reverse_iterator<
64  Event() = delete;
65 
68 
70  std::array<Particle, 2> const & incoming() const{
71  return incoming_;
72  }
74  std::vector<Particle> const & outgoing() const{
75  return outgoing_;
76  }
81 
86 
95 
97 
101  std::unordered_map<std::size_t, std::vector<Particle>> const & decays()
102  const {
103  return decays_;
104  }
106  std::vector<fastjet::PseudoJet> const & jets() const{
107  return jets_;
108  }
110 
113 
116  return parameters_;
117  }
120  return parameters_;
121  }
122 
124  EventParameters const & central() const{
125  return parameters_.central;
126  }
129  return parameters_.central;
130  }
131 
133  std::vector<EventParameters> const & variations() const{
134  return parameters_.variations;
135  }
137  std::vector<EventParameters> & variations(){
138  return parameters_.variations;
139  }
140 
142 
145  EventParameters const & variations(std::size_t i) const{
146  return parameters_.variations.at(i);
147  }
149 
152  EventParameters & variations(std::size_t i){
153  return parameters_.variations.at(i);
154  }
156 
158 
165  std::vector<int> particle_jet_indices(
166  std::vector<fastjet::PseudoJet> const & jets
167  ) const {
168  return cs_.particle_jet_indices(jets);
169  }
171  std::vector<int> particle_jet_indices() const {
172  return particle_jet_indices(jets());
173  }
174 
176  fastjet::JetDefinition const & jet_def() const{
177  return cs_.jet_def();
178  }
179 
181  double min_jet_pt() const{
182  return min_jet_pt_;
183  }
184 
187  return type_;
188  }
189 
191 
199  bool generate_colours(RNG & /*ran*/);
200 
202 
212  bool is_leading_colour() const;
213 
231  bool valid_hej_state(double soft_pt_regulator = DEFAULT_SOFT_PT_REGULATOR) const;
232 
234 
239  bool valid_incoming() const;
240 
241  private:
244 
249  Event(
250  std::array<Particle, 2> && incoming,
251  std::vector<Particle> && outgoing,
252  std::unordered_map<std::size_t, std::vector<Particle>> && decays,
254  fastjet::JetDefinition const & jet_def,
255  double min_jet_pt
256  );
257 
259  using PartonIterator = boost::filter_iterator<
260  bool (*)(Particle const &),
261  std::vector<Particle>::iterator
262  >;
264  using ReversePartonIterator = std::reverse_iterator<PartonIterator>;
265 
267  PartonIterator begin_partons();
269  PartonIterator end_partons();
270 
272  ReversePartonIterator rbegin_partons();
274  ReversePartonIterator rend_partons();
275 
276  std::array<Particle, 2> incoming_;
277  std::vector<Particle> outgoing_;
278  std::unordered_map<std::size_t, std::vector<Particle>> decays_;
279  std::vector<fastjet::PseudoJet> jets_;
280  Parameters<EventParameters> parameters_;
281  fastjet::ClusterSequence cs_;
282  double min_jet_pt_;
283  event_type::EventType type_;
284  }; // end class Event
285 
286 
288  inline
289  bool is_backward_g_to_h(Event const & ev) {
290  return ev.outgoing().front().type == pid::higgs
291  && ev.incoming().front().type == pid::gluon;
292  }
293 
295  inline
296  bool is_forward_g_to_h(Event const & ev) {
297  return ev.outgoing().back().type == pid::higgs
298  && ev.incoming().back().type == pid::gluon;
299  }
300 
303  public:
305  EventData() = default;
307  EventData(LHEF::HEPEUP const & hepeup);
310  std::array<Particle, 2> incoming,
311  std::vector<Particle> outgoing,
312  std::unordered_map<std::size_t, std::vector<Particle>> decays,
314  ):
315  incoming(std::move(incoming)), outgoing(std::move(outgoing)),
316  decays(std::move(decays)), parameters(std::move(parameters))
317  {}
318 
320 
332  fastjet::JetDefinition const & jet_def, double min_jet_pt);
333 
336  fastjet::JetDefinition const & jet_def, double const min_jet_pt){
337  return cluster(jet_def, min_jet_pt);
338  }
339 
341  void sort();
342 
344 
349  void reconstruct_intermediate(EWConstants const & /*ew_parameters*/);
350 
352 
359  void repair_momenta(double tolerance);
360 
362  std::array<Particle, 2> incoming;
364  std::vector<Particle> outgoing;
366  std::unordered_map<std::size_t, std::vector<Particle>> decays;
369  }; // end class EventData
370 
372  std::ostream& operator<<(std::ostream & os, Event const & ev);
373 
374  std::string to_string(Event const & ev);
375 
377  double shat(Event const & ev);
378 
380  static constexpr double TOL = 1e-6;
381 
383  LHEF::HEPEUP to_HEPEUP(Event const & event, LHEF::HEPRUP * /*heprup*/);
384 
385 } // namespace HEJ
Header file defining all global constants used for HEJ.
Containers for Parameter variations, e.g. different Weights.
Contains the particle struct.
Collection of electro-weak constants.
Definition: EWConstants.hh:24
Class to store general Event setup, i.e. Phase space and weights.
Definition: Event.hh:302
void reconstruct_intermediate(EWConstants const &)
Reconstruct intermediate particles from final-state leptons.
Parameters< EventParameters > parameters
Parameters, e.g. scale or inital weight.
Definition: Event.hh:368
Event operator()(fastjet::JetDefinition const &jet_def, double const min_jet_pt)
Alias for cluster()
Definition: Event.hh:335
void repair_momenta(double tolerance)
Repair momenta of massless particles.
EventData(std::array< Particle, 2 > incoming, std::vector< Particle > outgoing, std::unordered_map< std::size_t, std::vector< Particle >> decays, Parameters< EventParameters > parameters)
Constructor with all values given.
Definition: Event.hh:309
EventData()=default
Default Constructor.
std::vector< Particle > outgoing
Outcoing particles.
Definition: Event.hh:364
std::unordered_map< std::size_t, std::vector< Particle > > decays
Particle decays in the format {outgoing index, decay products}.
Definition: Event.hh:366
Event cluster(fastjet::JetDefinition const &jet_def, double min_jet_pt)
Generate an Event from the stored EventData.
std::array< Particle, 2 > incoming
Incoming particles.
Definition: Event.hh:362
void sort()
Sort particles in rapidity.
EventData(LHEF::HEPEUP const &hepeup)
Constructor from LesHouches event information.
An event with clustered jets.
Definition: Event.hh:51
Parameters< EventParameters > & parameters()
All chosen parameter, i.e. scale choices.
Definition: Event.hh:119
ConstPartonIterator end_partons() const
Iterator to the end of the outgoing partons.
ConstPartonIterator cbegin_partons() const
Iterator to the first outgoing parton.
fastjet::JetDefinition const & jet_def() const
Jet definition used for clustering.
Definition: Event.hh:176
EventParameters & variations(std::size_t i)
Parameter (scale) variation.
Definition: Event.hh:152
ConstPartonIterator begin_partons() const
Iterator to the first outgoing parton.
boost::filter_iterator< bool(*)(Particle const &), std::vector< Particle >::const_iterator > ConstPartonIterator
Iterator over partons.
Definition: Event.hh:59
Parameters< EventParameters > const & parameters() const
All chosen parameter, i.e. scale choices (const version)
Definition: Event.hh:115
bool is_leading_colour() const
Check that current colours are leading in the high energy limit.
std::unordered_map< std::size_t, std::vector< Particle > > const & decays() const
Particle decays.
Definition: Event.hh:101
std::vector< fastjet::PseudoJet > const & jets() const
The jets formed by the outgoing partons, sorted in rapidity.
Definition: Event.hh:106
std::vector< EventParameters > & variations()
Parameter (scale) variations.
Definition: Event.hh:137
double min_jet_pt() const
Minimum jet transverse momentum.
Definition: Event.hh:181
EventParameters & central()
Central parameter choice.
Definition: Event.hh:128
ConstPartonIterator cend_partons() const
Iterator to the end of the outgoing partons.
ConstReversePartonIterator crend_partons() const
Reverse Iterator to the first outgoing parton.
ConstReversePartonIterator rend_partons() const
Reverse Iterator to the first outgoing parton.
std::vector< Particle > const & outgoing() const
Outgoing particles.
Definition: Event.hh:74
ConstReversePartonIterator rbegin_partons() const
Reverse Iterator to the first outgoing parton.
std::array< Particle, 2 > const & incoming() const
Incoming particles.
Definition: Event.hh:70
std::vector< int > particle_jet_indices(std::vector< fastjet::PseudoJet > const &jets) const
Indices of the jets the outgoing partons belong to.
Definition: Event.hh:165
ConstReversePartonIterator crbegin_partons() const
Reverse Iterator to the first outgoing parton.
std::vector< EventParameters > const & variations() const
Parameter (scale) variations (const version)
Definition: Event.hh:133
EventParameters const & central() const
Central parameter choice (const version)
Definition: Event.hh:124
EventParameters const & variations(std::size_t i) const
Parameter (scale) variation (const version)
Definition: Event.hh:145
std::vector< int > particle_jet_indices() const
particle_jet_indices() of the Event jets()
Definition: Event.hh:171
event_type::EventType type() const
Event type.
Definition: Event.hh:186
std::reverse_iterator< ConstPartonIterator > ConstReversePartonIterator
Reverse Iterator over partons.
Definition: Event.hh:62
bool valid_hej_state(double soft_pt_regulator=DEFAULT_SOFT_PT_REGULATOR) const
Check if given event could have been produced by HEJ.
bool generate_colours(RNG &)
Give colours to each particle.
Event()=delete
No default Constructor.
bool valid_incoming() const
Check that the incoming momenta are valid.
Define different types of events.
EventType
Possible event types.
Definition: event_types.hh:19
@ higgs
Definition: PDG_codes.hh:87
@ e
Definition: PDG_codes.hh:40
@ gluon
Definition: PDG_codes.hh:76
Main HEJ 2 Namespace.
Definition: mainpage.dox:1
size_t implemented_types(std::vector< Particle > const &bosons)
Method for accessing implemented types.
LHEF::HEPEUP to_HEPEUP(Event const &event, LHEF::HEPRUP *)
Convert an event to a LHEF::HEPEUP.
std::string to_string(Event const &ev)
bool is_forward_g_to_h(Event const &ev)
Detect if a forward incoming gluon turns into a forward outgoing Higgs boson.
Definition: Event.hh:296
bool is_backward_g_to_h(Event const &ev)
Detect if a backward incoming gluon turns into a backward outgoing Higgs boson.
Definition: Event.hh:289
double shat(Event const &ev)
Square of the partonic centre-of-mass energy .
std::ostream & operator<<(std::ostream &os, const CrossSectionAccumulator &xs)
Print CrossSectionAccumulator to stream.
constexpr double DEFAULT_SOFT_PT_REGULATOR
default value for the maximal pt fraction of soft radiation in any tagging jets This cut is needed to...
Definition: Constants.hh:29
Definition: Analysis.hh:14
Definition: Event.hh:34
Event parameters.
Definition: Parameters.hh:107
Collection of parameters, e.g. Weights, assigned to a single event.
Definition: Parameters.hh:26
Class representing a particle.
Definition: Particle.hh:25
Interface for random number generator.
Definition: RNG.hh:22