Loading [MathJax]/extensions/tex2jax.js
hej is hosted by Hepforge, IPPP Durham
HEJ  2.3.0
High energy resummation for hadron colliders
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 
29  public:
36  std::function<double (double)> alpha_s,
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_Zphoton(
136  Event const & event,
137  double mur,
138  Particle const & boson
139  ) const;
140  double virtual_corrections_Zphoton_single(
141  Event const & event,
142  double mur,
143  Particle const & boson,
144  bool emit_fwd
145  ) const;
146  std::vector<double> virtual_corrections_Zphoton_double(
147  Event const & event,
148  double mur,
149  Particle const & boson
150  ) const;
151  std::vector<double> virtual_corrections_Zphoton_triple(
152  Event const & event,
153  double mur,
154  Particle const & boson
155  ) const;
156 
157  template<class InputIterator>
158  double virtual_corrections_no_interference(
159  InputIterator begin_parton, InputIterator end_parton,
160  fastjet::PseudoJet & q,
161  const double mur
162  ) const;
163 
164  template<class InputIterator>
165  std::vector<double> virtual_corrections_interference(
166  InputIterator begin_parton, InputIterator end_parton,
167  fastjet::PseudoJet & q_t,
168  fastjet::PseudoJet & q_b,
169  const double mur
170  ) const;
171 
172  std::vector<double> virtual_corrections(
173  Event const & event,
174  double mur
175  ) const;
176 
178  double omega0(
179  double alpha_s, double mur,
180  fastjet::PseudoJet const & q_j
181  ) const;
182 
183  double tree_kin_jets(
184  Event const & ev
185  ) const;
186  double tree_kin_W(
187  Event const & ev
188  ) const;
189  std::vector <double> tree_kin_WW(
190  Event const & ev
191  ) const;
192  std::vector <double> tree_kin_Z(
193  Event const & ev
194  ) const;
195  std::vector <double> tree_kin_photon(
196  Event const & ev
197  ) const;
198  double tree_kin_Higgs(
199  Event const & ev
200  ) const;
201  double tree_kin_Higgs_first(
202  Event const & ev
203  ) const;
204  double tree_kin_Higgs_last(
205  Event const & ev
206  ) const;
207 
208  double tree_kin_Higgs_between(
209  Event const & ev
210  ) const;
211 
212  double tree_param_partons(
213  double alpha_s, double mur,
214  std::vector<Particle> const & partons
215  ) const;
216 
217  std::vector<int> in_extremal_jet_indices(
218  std::vector<fastjet::PseudoJet> const & partons
219  ) const;
220 
221  double MH2_backwardH(
222  ParticleID type_forward,
223  CLHEP::HepLorentzVector const & pa,
224  CLHEP::HepLorentzVector const & pb,
225  CLHEP::HepLorentzVector const & pH,
226  CLHEP::HepLorentzVector const & pn
227  ) const;
228 
229  double MH2_unob_forwardH(
230  CLHEP::HepLorentzVector const & pa,
231  CLHEP::HepLorentzVector const & pb,
232  CLHEP::HepLorentzVector const & pg,
233  CLHEP::HepLorentzVector const & p1,
234  CLHEP::HepLorentzVector const & pH
235  ) const;
236 
237  std::function<double (double)> alpha_s_;
238 
239  MatrixElementConfig param_;
240  };
241 
242 } // namespace HEJ
HEJ 2 configuration parameters.
Contains the Particle IDs of all relevant SM particles.
Containers for Parameter variations, e.g. different Weights.
An event with clustered jets.
Definition: Event.hh:51
Class to calculate the squares of matrix elements.
Definition: MatrixElement.hh:28
Weights tree_param(Event const &event) const
Scale-dependent part of tree-level matrix element squares.
Weights operator()(Event const &event) const
squares of regulated HEJ matrix elements
std::vector< Weights > virtual_corrections(Event const &event) const
Virtual corrections to matrix element squares.
std::vector< double > tree_kin(Event const &event) const
Kinematic part of tree-level matrix element squares.
Weights tree(Event const &event) const
Squares of HEJ tree-level matrix elements.
MatrixElement(std::function< double(double)> alpha_s, MatrixElementConfig conf)
MatrixElement Constructor.
Definition: MatrixElement.hh:19
ParticleID
The possible particle identities. We use PDG IDs as standard.
Definition: PDG_codes.hh:25
Main HEJ 2 Namespace.
Definition: mainpage.dox:1
Configuration options for the MatrixElement class.
Definition: Config.hh:173
Collection of parameters, e.g. Weights, assigned to a single event.
Definition: Parameters.hh:26
Class representing a particle.
Definition: Particle.hh:25