hej
is hosted by
Hepforge
,
IPPP Durham
HEJ 2
2.0
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
include
HEJ
Event.hh
Go to the documentation of this file.
1
9
#pragma once
10
11
#include <array>
12
#include <memory>
13
#include <string>
14
#include <unordered_map>
15
#include <vector>
16
17
#include "
HEJ/event_types.hh
"
18
#include "
HEJ/Particle.hh
"
19
20
#include "fastjet/ClusterSequence.hh"
21
22
namespace
LHEF
{
23
class
HEPEUP;
24
class
HEPRUP;
25
}
26
27
namespace
fastjet
{
28
class
JetDefinition;
29
}
30
31
namespace
HEJ
{
32
33
34
struct
ParameterDescription;
35
37
struct
EventParameters
{
38
double
mur
;
39
double
muf
;
40
double
weight
;
42
std::shared_ptr<ParameterDescription>
description
=
nullptr
;
43
};
44
46
struct
ParameterDescription
{
48
std::string
scale_name
;
50
double
mur_factor
;
52
double
muf_factor
;
53
54
ParameterDescription
() =
default
;
55
ParameterDescription
(
56
std::string
scale_name
,
double
mur_factor
,
double
muf_factor
57
):
58
scale_name
{
scale_name
},
mur_factor
{
mur_factor
},
muf_factor
{
muf_factor
}
59
{};
60
};
61
63
struct
UnclusteredEvent
{
65
UnclusteredEvent
() =
default
;
67
UnclusteredEvent
(LHEF::HEPEUP
const
& hepeup);
68
69
std::array<Particle, 2>
incoming
;
70
std::vector<Particle>
outgoing
;
72
std::unordered_map<size_t, std::vector<Particle>>
decays
;
74
EventParameters
central
;
75
std::vector<EventParameters>
variations
;
76
};
77
84
class
Event
{
85
public
:
87
Event
() =
default
;
89
Event
(
90
UnclusteredEvent
ev,
91
fastjet::JetDefinition
const
&
jet_def
,
double
min_jet_pt
92
);
93
95
std::vector<fastjet::PseudoJet>
jets
()
const
;
96
98
UnclusteredEvent
const
&
unclustered
()
const
{
99
return
ev_;
100
}
101
103
EventParameters
const
&
central
()
const
{
104
return
ev_.
central
;
105
}
106
108
EventParameters
&
central
(){
109
return
ev_.
central
;
110
}
111
113
std::array<Particle, 2>
const
&
incoming
()
const
{
114
return
ev_.
incoming
;
115
}
116
118
std::vector<Particle>
const
&
outgoing
()
const
{
119
return
ev_.
outgoing
;
120
}
121
123
127
std::unordered_map<size_t, std::vector<Particle>>
const
&
decays
()
const
{
128
return
ev_.
decays
;
129
}
130
132
std::vector<EventParameters>
const
&
variations
()
const
{
133
return
ev_.
variations
;
134
}
135
137
std::vector<EventParameters> &
variations
(){
138
return
ev_.
variations
;
139
}
140
142
145
EventParameters
const
&
variations
(
size_t
i)
const
{
146
return
ev_.
variations
[i];
147
}
148
150
153
EventParameters
&
variations
(
size_t
i){
154
return
ev_.
variations
[i];
155
}
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
}
170
172
fastjet::JetDefinition
const
&
jet_def
()
const
{
173
return
cs_.jet_def();
174
}
175
177
double
min_jet_pt
()
const
{
178
return
min_jet_pt_;
179
}
180
182
event_type::EventType
type
()
const
{
183
return
type_;
184
}
185
186
private
:
187
UnclusteredEvent
ev_;
188
fastjet::ClusterSequence cs_;
189
double
min_jet_pt_;
190
event_type::EventType
type_;
191
};
192
194
double
shat
(
Event
const
& ev);
195
197
LHEF::HEPEUP
to_HEPEUP
(
Event
const
& event, LHEF::HEPRUP *);
198
199
}
Particle.hh
Contains the particle struct.
HEJ::Event
Definition:
Event.hh:84
HEJ::Event::particle_jet_indices
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
HEJ::Event::variations
EventParameters const & variations(size_t i) const
Parameter (scale) variation.
Definition:
Event.hh:145
HEJ::Event::Event
Event(UnclusteredEvent ev, fastjet::JetDefinition const &jet_def, double min_jet_pt)
Event Constructor adding jet clustering to an unclustered event.
HEJ::Event::central
EventParameters const & central() const
Central parameter choice.
Definition:
Event.hh:103
HEJ::Event::variations
std::vector< EventParameters > & variations()
Parameter (scale) variations.
Definition:
Event.hh:137
HEJ::Event::incoming
std::array< Particle, 2 > const & incoming() const
Incoming particles.
Definition:
Event.hh:113
HEJ::Event::variations
std::vector< EventParameters > const & variations() const
Parameter (scale) variations.
Definition:
Event.hh:132
HEJ::Event::min_jet_pt
double min_jet_pt() const
Minimum jet transverse momentum.
Definition:
Event.hh:177
HEJ::Event::variations
EventParameters & variations(size_t i)
Parameter (scale) variation.
Definition:
Event.hh:153
HEJ::Event::central
EventParameters & central()
Central parameter choice.
Definition:
Event.hh:108
HEJ::Event::decays
std::unordered_map< size_t, std::vector< Particle > > const & decays() const
Particle decays.
Definition:
Event.hh:127
HEJ::Event::jet_def
fastjet::JetDefinition const & jet_def() const
Jet definition used for clustering.
Definition:
Event.hh:172
HEJ::Event::unclustered
UnclusteredEvent const & unclustered() const
The corresponding event before jet clustering.
Definition:
Event.hh:98
HEJ::Event::type
event_type::EventType type() const
Event type.
Definition:
Event.hh:182
HEJ::Event::jets
std::vector< fastjet::PseudoJet > jets() const
The jets formed by the outgoing partons.
HEJ::Event::outgoing
std::vector< Particle > const & outgoing() const
Outgoing particles.
Definition:
Event.hh:118
HEJ::Event::Event
Event()=default
Default Event Constructor.
event_types.hh
Define different types of events.
HEJ::event_type::EventType
EventType
Possible event types.
Definition:
event_types.hh:18
HEJ
Main HEJ 2 Namespace.
Definition:
mainpage.dox:1
HEJ::to_HEPEUP
LHEF::HEPEUP to_HEPEUP(Event const &event, LHEF::HEPRUP *)
Convert an event to a LHEF::HEPEUP.
HEJ::shat
double shat(Event const &ev)
Square of the partonic centre-of-mass energy .
LHEF
Definition:
CombinedEventWriter.hh:17
fastjet
Definition:
Event.hh:27
HEJ::EventParameters
Event parameters.
Definition:
Event.hh:37
HEJ::EventParameters::mur
double mur
Definition:
Event.hh:38
HEJ::EventParameters::description
std::shared_ptr< ParameterDescription > description
Optional description.
Definition:
Event.hh:42
HEJ::EventParameters::muf
double muf
Definition:
Event.hh:39
HEJ::EventParameters::weight
double weight
Definition:
Event.hh:40
HEJ::ParameterDescription
Description of event parameters.
Definition:
Event.hh:46
HEJ::ParameterDescription::ParameterDescription
ParameterDescription(std::string scale_name, double mur_factor, double muf_factor)
Definition:
Event.hh:55
HEJ::ParameterDescription::muf_factor
double muf_factor
Actual factorisation scale divided by central scale.
Definition:
Event.hh:52
HEJ::ParameterDescription::scale_name
std::string scale_name
Name of central scale choice (e.g. "H_T/2")
Definition:
Event.hh:48
HEJ::ParameterDescription::ParameterDescription
ParameterDescription()=default
HEJ::ParameterDescription::mur_factor
double mur_factor
Actual renormalisation scale divided by central scale.
Definition:
Event.hh:50
HEJ::UnclusteredEvent
An event before jet clustering.
Definition:
Event.hh:63
HEJ::UnclusteredEvent::outgoing
std::vector< Particle > outgoing
Definition:
Event.hh:70
HEJ::UnclusteredEvent::UnclusteredEvent
UnclusteredEvent(LHEF::HEPEUP const &hepeup)
Constructor from LesHouches event information.
HEJ::UnclusteredEvent::decays
std::unordered_map< size_t, std::vector< Particle > > decays
Particle decays in the format {outgoing index, decay products}.
Definition:
Event.hh:72
HEJ::UnclusteredEvent::incoming
std::array< Particle, 2 > incoming
Definition:
Event.hh:69
HEJ::UnclusteredEvent::central
EventParameters central
Central parameter (e.g. scale) choice.
Definition:
Event.hh:74
HEJ::UnclusteredEvent::UnclusteredEvent
UnclusteredEvent()=default
Default Constructor.
HEJ::UnclusteredEvent::variations
std::vector< EventParameters > variations
Definition:
Event.hh:75
Generated by
1.9.5