hej
is hosted by
Hepforge
,
IPPP Durham
HEJ
2.2.2
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
include
HEJ
MatrixElement.hh
Go to the documentation of this file.
1
8
#pragma once
9
10
#include <functional>
11
#include <vector>
12
13
#include "fastjet/PseudoJet.hh"
14
15
#include "
HEJ/Config.hh
"
16
#include "
HEJ/PDG_codes.hh
"
17
#include "
HEJ/Parameters.hh
"
18
19
namespace
CLHEP
{
20
class
HepLorentzVector;
21
}
22
23
namespace
HEJ
{
24
class
Event;
25
struct
Particle;
26
28
class
MatrixElement
{
29
public
:
35
MatrixElement
(
36
std::function<
double
(
double
)> alpha_s,
37
MatrixElementConfig
conf
38
);
39
54
Weights
operator()
(
Event
const
& event)
const
;
55
57
63
Weights
tree
(
Event
const
& event)
const
;
64
82
std::vector<Weights>
virtual_corrections
(
Event
const
& event)
const
;
83
97
Weights
tree_param
(
Event
const
& event)
const
;
98
118
std::vector<double>
tree_kin
(
Event
const
& event)
const
;
119
120
private
:
121
double
tree_param
(
122
Event
const
& event,
123
double
mur
124
)
const
;
125
126
double
virtual_corrections_W(
127
Event
const
& event,
128
double
mur,
129
Particle
const
& WBoson
130
)
const
;
131
std::vector <double> virtual_corrections_WW(
132
Event
const
& event,
133
double
mur
134
)
const
;
135
std::vector <double> virtual_corrections_Z_qq(
136
Event
const
& event,
137
double
mur,
138
Particle
const
& ZBoson
139
)
const
;
140
double
virtual_corrections_Z_qg(
141
Event
const
& event,
142
double
mur,
143
Particle
const
& ZBoson,
144
bool
is_gq_event
145
)
const
;
146
147
template
<
class
InputIterator>
148
std::vector <double> virtual_corrections_interference(
149
InputIterator begin_parton, InputIterator end_parton,
150
fastjet::PseudoJet
const
& q0_t,
151
fastjet::PseudoJet
const
& q0_b,
152
const
double
mur
153
)
const
;
154
155
std::vector<double>
virtual_corrections
(
156
Event
const
& event,
157
double
mur
158
)
const
;
159
161
double
omega0(
162
double
alpha_s,
double
mur,
163
fastjet::PseudoJet
const
& q_j
164
)
const
;
165
166
double
tree_kin_jets(
167
Event
const
& ev
168
)
const
;
169
double
tree_kin_W(
170
Event
const
& ev
171
)
const
;
172
std::vector <double> tree_kin_WW(
173
Event
const
& ev
174
)
const
;
175
std::vector <double> tree_kin_Z(
176
Event
const
& ev
177
)
const
;
178
double
tree_kin_Higgs(
179
Event
const
& ev
180
)
const
;
181
double
tree_kin_Higgs_first(
182
Event
const
& ev
183
)
const
;
184
double
tree_kin_Higgs_last(
185
Event
const
& ev
186
)
const
;
187
188
double
tree_kin_Higgs_between(
189
Event
const
& ev
190
)
const
;
191
192
double
tree_param_partons(
193
double
alpha_s,
double
mur,
194
std::vector<Particle>
const
& partons
195
)
const
;
196
197
std::vector<int> in_extremal_jet_indices(
198
std::vector<fastjet::PseudoJet>
const
& partons
199
)
const
;
200
201
double
MH2_backwardH(
202
ParticleID
type_forward,
203
CLHEP::HepLorentzVector
const
& pa,
204
CLHEP::HepLorentzVector
const
& pb,
205
CLHEP::HepLorentzVector
const
& pH,
206
CLHEP::HepLorentzVector
const
& pn
207
)
const
;
208
209
double
MH2_unob_forwardH(
210
CLHEP::HepLorentzVector
const
& pa,
211
CLHEP::HepLorentzVector
const
& pb,
212
CLHEP::HepLorentzVector
const
& pg,
213
CLHEP::HepLorentzVector
const
& p1,
214
CLHEP::HepLorentzVector
const
& pH
215
)
const
;
216
217
std::function<double (
double
)> alpha_s_;
218
219
MatrixElementConfig
param_;
220
};
221
222
}
// namespace HEJ
Config.hh
HEJ 2 configuration parameters.
PDG_codes.hh
Contains the Particle IDs of all relevant SM particles.
Parameters.hh
Containers for Parameter variations, e.g. different Weights.
HEJ::Event
An event with clustered jets.
Definition:
Event.hh:51
HEJ::MatrixElement
Class to calculate the squares of matrix elements.
Definition:
MatrixElement.hh:28
HEJ::MatrixElement::tree_param
Weights tree_param(Event const &event) const
Scale-dependent part of tree-level matrix element squares.
HEJ::MatrixElement::operator()
Weights operator()(Event const &event) const
squares of regulated HEJ matrix elements
HEJ::MatrixElement::tree
Weights tree(Event const &event) const
Squares of HEJ tree-level matrix elements.
HEJ::MatrixElement::tree_kin
std::vector< double > tree_kin(Event const &event) const
Kinematic part of tree-level matrix element squares.
HEJ::MatrixElement::MatrixElement
MatrixElement(std::function< double(double)> alpha_s, MatrixElementConfig conf)
MatrixElement Constructor.
HEJ::MatrixElement::virtual_corrections
std::vector< Weights > virtual_corrections(Event const &event) const
Virtual corrections to matrix element squares.
CLHEP
Definition:
MatrixElement.hh:19
HEJ::pid::ParticleID
ParticleID
The possible particle identities. We use PDG IDs as standard.
Definition:
PDG_codes.hh:25
HEJ
Main HEJ 2 Namespace.
Definition:
mainpage.dox:1
HEJ::MatrixElementConfig
Configuration options for the MatrixElement class.
Definition:
Config.hh:172
HEJ::Parameters
Collection of parameters, e.g. Weights, assigned to a single event.
Definition:
Parameters.hh:26
HEJ::Particle
Class representing a particle.
Definition:
Particle.hh:25
Generated by
1.9.5