hej is hosted by Hepforge, IPPP Durham
HEJ 2 2.0
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
EventReweighter.hh
Go to the documentation of this file.
1
11#pragma once
12
13#include <array>
14#include <functional>
15#include <utility>
16#include <vector>
17
18#include "HEJ/config.hh"
19#include "HEJ/event_types.hh"
20#include "HEJ/MatrixElement.hh"
21#include "HEJ/PDF.hh"
22#include "HEJ/PDG_codes.hh"
23#include "HEJ/RNG.hh"
24#include "HEJ/ScaleFunction.hh"
25
26namespace LHEF {
27 class HEPRUP;
28}
29
30namespace HEJ{
31 class Event;
32 class Weights;
33
35
39 struct Beam{
40 double E;
41 std::array<ParticleID, 2> type;
42 };
43
47 public:
48
50 Beam beam,
51 int pdf_id,
52 ScaleGenerator scale_gen,
54 HEJ::RNG & ran
55 );
56
58 LHEF::HEPRUP const & heprup,
59 ScaleGenerator scale_gen,
61 HEJ::RNG & ran
62 );
63
65 PDF const & pdf() const;
66
67
69
86 std::vector<Event> reweight(
87 Event const & ev,
88 int num_events
89 );
90
91 private:
92
93 template<typename... T>
94 PDF const & pdf(T&& ...);
95
100 std::vector<Event> gen_res_events(
101 Event const & ev, int num_events
102 );
103 std::vector<Event> rescale(
104 Event const & Born_ev, std::vector<Event> events
105 ) const;
106
113 bool jets_pass_resummation_cuts(Event const & ev) const;
114
124 Weights pdf_factors(Event const & ev) const;
125
135 Weights matrix_elements(Event const & ev) const;
136
150 Weights fixed_order_scale_ME(Event const & ev) const;
151
161 double tree_matrix_element(Event const & ev) const;
162
163
166
168 double E_beam_;
169
171 PDF pdf_;
172
174 MatrixElement MEt2_;
176 ScaleGenerator scale_gen_;
182 std::reference_wrapper<HEJ::RNG> ran_;
183 };
184
185 template<typename... T>
186 PDF const & EventReweighter::pdf(T&&... t){
187 return pdf_ = PDF{std::forward<T>(t)...};
188 }
189
190}
Contains the MatrixElement Class.
Contains all the necessary classes and functions for interaction with PDFs.
Contains the Particle IDs of all relevant SM particles.
Interface for pseudorandom number generators.
Functions to calculate the (renormalisation and factorisation) scales for an event.
Main class for reweighting events in HEJ.
Definition: EventReweighter.hh:45
std::vector< Event > reweight(Event const &ev, int num_events)
Generate resummation events for a given fixed-order event.
EventReweighter(LHEF::HEPRUP const &heprup, ScaleGenerator scale_gen, EventReweighterConfig conf, HEJ::RNG &ran)
EventReweighter(Beam beam, int pdf_id, ScaleGenerator scale_gen, EventReweighterConfig conf, HEJ::RNG &ran)
PDF const & pdf() const
Get the used pdf.
Definition: Event.hh:84
Class to calculate the squares of matrix elements.
Definition: MatrixElement.hh:28
Class for interaction with a PDF set.
Definition: PDF.hh:21
Generate combinations of renormalisation and factorisation scales.
Definition: ScaleFunction.hh:90
HEJ 2 configuration parameters.
Define different types of events.
EventType
Possible event types.
Definition: event_types.hh:18
Main HEJ 2 Namespace.
Definition: mainpage.dox:1
Definition: CombinedEventWriter.hh:17
Beam parameters.
Definition: EventReweighter.hh:39
double E
Definition: EventReweighter.hh:40
std::array< ParticleID, 2 > type
Definition: EventReweighter.hh:41
Configuration options for the EventReweighter class.
Definition: config.hh:124
Interface for random number generator.
Definition: RNG.hh:19
Collection of weights assigned to a single event.
Definition: Weights.hh:19