hej
is hosted by
Hepforge
,
IPPP Durham
HEJ
2.1.4
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
52
Weights
operator()
(
Event
const
& event)
const
;
53
55
61
Weights
tree
(
Event
const
& event)
const
;
62
80
std::vector<Weights>
virtual_corrections
(
Event
const
& event)
const
;
81
95
Weights
tree_param
(
Event
const
& event)
const
;
96
114
std::vector<double>
tree_kin
(
Event
const
& event)
const
;
115
116
private
:
117
double
tree_param
(
118
Event
const
& event,
119
double
mur
120
)
const
;
121
122
double
virtual_corrections_W(
123
Event
const
& event,
124
double
mur,
125
Particle
const
& WBoson
126
)
const
;
127
std::vector <double> virtual_corrections_Z_qq(
128
Event
const
& event,
129
double
mur,
130
Particle
const
& ZBoson
131
)
const
;
132
double
virtual_corrections_Z_qg(
133
Event
const
& event,
134
double
mur,
135
Particle
const
& ZBoson,
136
bool
is_gq_event
137
)
const
;
138
std::vector<double>
virtual_corrections
(
139
Event
const
& event,
140
double
mur
141
)
const
;
142
144
double
omega0(
145
double
alpha_s,
double
mur,
146
fastjet::PseudoJet
const
& q_j
147
)
const
;
148
149
double
tree_kin_jets(
150
Event
const
& ev
151
)
const
;
152
double
tree_kin_W(
153
Event
const
& ev
154
)
const
;
155
std::vector <double> tree_kin_Z(
156
Event
const
& ev
157
)
const
;
158
double
tree_kin_Higgs(
159
Event
const
& ev
160
)
const
;
161
double
tree_kin_Higgs_first(
162
Event
const
& ev
163
)
const
;
164
double
tree_kin_Higgs_last(
165
Event
const
& ev
166
)
const
;
167
176
double
tree_kin_Higgs_between(
177
Event
const
& ev
178
)
const
;
179
180
double
tree_param_partons(
181
double
alpha_s,
double
mur,
182
std::vector<Particle>
const
& partons
183
)
const
;
184
185
std::vector<int> in_extremal_jet_indices(
186
std::vector<fastjet::PseudoJet>
const
& partons
187
)
const
;
188
189
double
MH2_forwardH(
190
CLHEP::HepLorentzVector
const
& p1out,
191
CLHEP::HepLorentzVector
const
& p1in,
192
pid::ParticleID
type2,
193
CLHEP::HepLorentzVector
const
& p2out,
194
CLHEP::HepLorentzVector
const
& p2in,
195
CLHEP::HepLorentzVector
const
& pH,
196
double
t1,
double
t2
197
)
const
;
198
199
std::function<double (
double
)> alpha_s_;
200
201
MatrixElementConfig
param_;
202
};
203
204
}
// 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:47
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:23
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:24
Generated by
1.9.5