hej is hosted by Hepforge, IPPP Durham
HEJ  2.1.4
High energy resummation for hadron colliders
EventReweighter.hh
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <array>
14 #include <cstddef>
15 #include <memory>
16 #include <utility>
17 #include <vector>
18 
19 #include "HEJ/Config.hh"
20 #include "HEJ/MatrixElement.hh"
21 #include "HEJ/PDF.hh"
22 #include "HEJ/PDG_codes.hh"
23 #include "HEJ/Parameters.hh"
24 #include "HEJ/ScaleFunction.hh"
25 #include "HEJ/StatusCode.hh"
26 #include "HEJ/event_types.hh"
27 
28 namespace LHEF {
29  class HEPRUP;
30 }
31 
32 namespace HEJ {
33  class Event;
34  struct RNG;
35 
37 
41  struct Beam{
42  double E;
43  std::array<ParticleID, 2> type;
44  };
45 
49 
50  public:
52  Beam const & beam,
53  int pdf_id,
54  ScaleGenerator scale_gen,
56  std::shared_ptr<RNG> ran
57  );
58 
60  LHEF::HEPRUP const & heprup,
61  ScaleGenerator scale_gen,
63  std::shared_ptr<RNG> ran
64  );
65 
67  PDF const & pdf() const;
68 
71 
73 
89  std::vector<Event> reweight(
90  Event const & ev,
91  std::size_t num_events
92  );
93 
95 
99  std::vector<StatusCode> const & status() const {
100  return status_;
101  }
102 
103  private:
104 
109  std::vector<Event> gen_res_events(
110  Event const & ev, std::size_t phase_space_points
111  );
112  std::vector<Event> rescale(
113  Event const & Born_ev, std::vector<Event> events
114  ) const;
115 
122  bool jets_pass_resummation_cuts(Event const & ev) const;
123 
133  Weights pdf_factors(Event const & ev) const;
134 
144  Weights matrix_elements(Event const & ev) const;
145 
159  Weights fixed_order_scale_ME(Event const & ev) const;
160 
170  double tree_matrix_element(Event const & ev) const;
171 
173  EventReweighterConfig param_;
174 
176  double E_beam_;
177 
179  PDF pdf_;
180 
182  MatrixElement MEt2_;
184  ScaleGenerator scale_gen_;
186  std::shared_ptr<RNG> ran_;
188  std::vector<StatusCode> status_;
189  };
190 
191 } // namespace HEJ
HEJ 2 configuration parameters.
Contains the MatrixElement Class.
Contains all the necessary classes and functions for interaction with PDFs.
Contains the Particle IDs of all relevant SM particles.
Containers for Parameter variations, e.g. different Weights.
Functions to calculate the (renormalisation and factorisation) scales for an event.
Header file for status codes of event generation.
Main class for reweighting events in HEJ.
Definition: EventReweighter.hh:47
std::vector< Event > reweight(Event const &ev, std::size_t num_events)
Generate resummation events for a given fixed-order event.
EventReweighter(Beam const &beam, int pdf_id, ScaleGenerator scale_gen, EventReweighterConfig conf, std::shared_ptr< RNG > ran)
std::vector< StatusCode > const & status() const
Gives all StatusCodes of the last reweight()
Definition: EventReweighter.hh:99
PDF const & pdf() const
Get the used pdf.
EventTreatment treatment(EventType type) const
Get event treatment.
EventReweighter(LHEF::HEPRUP const &heprup, ScaleGenerator scale_gen, EventReweighterConfig conf, std::shared_ptr< RNG > ran)
An event with clustered jets.
Definition: Event.hh:47
Class to calculate the squares of matrix elements.
Definition: MatrixElement.hh:28
Class for interaction with a PDF set.
Definition: PDF.hh:23
Generate combinations of renormalisation and factorisation scales.
Definition: ScaleFunction.hh:104
Define different types of events.
@ beam
Beam.
Definition: HepMCInterface_common.hh:26
EventType
Possible event types.
Definition: event_types.hh:19
Main HEJ 2 Namespace.
Definition: mainpage.dox:1
EventTreatment
Definition: Config.hh:66
Definition: Analysis.hh:14
Beam parameters.
Definition: EventReweighter.hh:41
double E
Definition: EventReweighter.hh:42
std::array< ParticleID, 2 > type
Definition: EventReweighter.hh:43
Configuration options for the EventReweighter class.
Definition: Config.hh:202
Collection of parameters, e.g. Weights, assigned to a single event.
Definition: Parameters.hh:26