cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
N/A
Central exclusive processes event generator
Toggle main menu visibility
Main Page
Related Pages
Packages
Package List
Package Members
All
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
x
y
Functions
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
x
y
Variables
Typedefs
Enumerations
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Variables
a
b
c
d
e
f
g
h
i
k
l
m
n
p
q
r
s
t
u
v
w
x
y
z
Typedefs
Enumerations
Enumerator
b
c
d
e
g
h
i
l
m
n
p
r
t
u
w
x
y
z
Related Symbols
d
g
h
o
u
v
Files
File List
File Members
All
_
a
b
c
d
e
h
o
p
r
s
Functions
Variables
Macros
_
b
c
d
e
p
r
s
▼
CepGen
Reference manual
Bibliography
►
Packages
►
Classes
▼
Files
▼
File List
▼
include
▼
CepGen
►
Cards
►
Core
►
Event
►
EventFilter
►
FormFactors
►
Integration
►
Modules
►
PartonFluxes
►
Physics
▼
Process
►
Fortran
►
FactorisedProcess.h
►
FortranFactorisedProcess.h
►
PartonsPhaseSpaceGenerator.h
►
PhaseSpaceGenerator.h
►
Process.h
►
StructureFunctions
►
Utils
►
Generator.h
►
Version.h
►
CepGenBoost
►
CepGenHepMC2
►
CepGenHepMC3
►
CepGenHerwig6
►
CepGenMadGraph
►
CepGenPythia6
►
CepGenPythia8
►
CepGenPython
►
CepGenRoot
►
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
Loading...
Searching...
No Matches
Process.h
Go to the documentation of this file.
1
/*
2
* CepGen: a central exclusive processes event generator
3
* Copyright (C) 2013-2025 Laurent Forthomme
4
*
5
* This program is free software: you can redistribute it and/or modify
6
* it under the terms of the GNU General Public License as published by
7
* the Free Software Foundation, either version 3 of the License, or
8
* any later version.
9
*
10
* This program is distributed in the hope that it will be useful,
11
* but WITHOUT ANY WARRANTY; without even the implied warranty of
12
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13
* GNU General Public License for more details.
14
*
15
* You should have received a copy of the GNU General Public License
16
* along with this program. If not, see <http://www.gnu.org/licenses/>.
17
*/
18
19
#ifndef CepGen_Process_Process_h
20
#define CepGen_Process_Process_h
21
22
#include <cmath>
23
24
#include "
CepGen/Event/Event.h
"
25
#include "
CepGen/Physics/Coupling.h
"
26
#include "
CepGen/Physics/Kinematics.h
"
27
#include "
CepGen/Utils/ProcessVariablesAnalyser.h
"
28
30
namespace
cepgen::proc
{
34
class
Process
:
public
NamedModule
<Process> {
35
public
:
36
explicit
Process
(
const
ParametersList
&);
37
Process
(
const
Process
&);
38
~Process
()
override
=
default
;
39
Process
&
operator=
(
const
Process
&);
40
41
static
ParametersDescription
description
();
42
43
virtual
std::unique_ptr<Process>
clone
()
const
;
44
virtual
double
computeWeight
() = 0;
45
46
void
clear
();
47
void
clearEvent
();
48
void
initialise
();
49
50
inline
const
Kinematics
&
kinematics
()
const
{
return
kin_; }
51
inline
Kinematics
&
kinematics
() {
return
kin_; }
52
void
setKinematics
();
53
54
// debugging utilities
55
double
weight
(
const
std::vector<double>&);
56
void
dumpPoint
(std::ostream* =
nullptr
)
const
;
57
void
dumpVariables
(std::ostream* =
nullptr
)
const
;
58
59
inline
size_t
ndim
()
const
{
return
mapped_variables_.size(); }
60
62
inline
bool
hasEvent
()
const
{
return
static_cast<
bool
>
(event_); }
63
const
Event
&
event
()
const
;
64
Event
&
event
();
65
Event
*
eventPtr
();
66
67
const
Momentum
&
pA
()
const
;
68
inline
double
mA2
()
const
{
return
mA2_; }
69
inline
double
mA
()
const
{
return
std::sqrt(
mA2
()); }
70
const
Momentum
&
pB
()
const
;
71
inline
double
mB2
()
const
{
return
mB2_; }
72
inline
double
mB
()
const
{
return
std::sqrt(
mB2
()); }
73
const
Momentum
&
pX
()
const
;
74
inline
double
mX2
()
const
{
return
mX2_; }
75
inline
double
mX
()
const
{
return
std::sqrt(
mX2
()); }
76
const
Momentum
&
pY
()
const
;
77
inline
double
mY2
()
const
{
return
mY2_; }
78
inline
double
mY
()
const
{
return
std::sqrt(
mY2
()); }
79
80
inline
double
x1
()
const
{
return
x1_; }
81
inline
double
t1
()
const
{
return
t1_; }
82
const
Momentum
&
q1
()
const
;
83
inline
double
x2
()
const
{
return
x2_; }
84
inline
double
t2
()
const
{
return
t2_; }
85
const
Momentum
&
q2
()
const
;
86
87
Momentum
&
q1
();
88
Momentum
&
q2
();
89
90
inline
double
wCM
()
const
{
return
wcm_; }
91
92
Momentum
&
pA
();
93
Momentum
&
pB
();
94
Momentum
&
pX
();
95
Momentum
&
pY
();
96
97
inline
double
&
mX2
() {
return
mX2_; }
98
inline
double
&
mY2
() {
return
mY2_; }
99
inline
double
&
t1
() {
return
t1_; }
100
inline
double
&
x1
() {
return
x1_; }
101
inline
double
&
t2
() {
return
t2_; }
102
inline
double
&
x2
() {
return
x2_; }
103
104
Momentum
&
pc
(
size_t
);
105
const
Momentum
&
pc
(
size_t
)
const
;
106
107
inline
double
s
()
const
{
return
s_; }
108
inline
double
sqrtS
()
const
{
return
sqs_; }
109
inline
double
inverseS
()
const
{
return
inv_s_; }
110
inline
double
inverseSqrtS
()
const
{
return
inv_sqs_; }
111
112
protected
:
113
virtual
void
addEventContent
() = 0;
114
inline
virtual
void
prepareKinematics
() {}
115
virtual
void
fillKinematics
() = 0;
116
117
//--- Mandelstam variables
118
double
shat
()
const
;
119
120
double
mp_
;
121
double
mp2_
;
122
123
public
:
125
enum class
Mapping
{
126
linear
= 0,
127
exponential
,
128
square
,
129
power_law
130
};
125
enum class
Mapping
{
…
};
131
friend
std::ostream&
operator<<
(std::ostream&,
const
Mapping
&);
139
Process
&
defineVariable
(
double
& out,
140
const
Mapping
& type,
141
const
Limits
& lim,
142
const
std::string&
name
,
143
const
std::string&
description
=
""
);
145
double
variableValue
(
size_t
i,
double
x)
const
;
146
147
protected
:
151
double
generateVariables
()
const
;
152
154
void
setEventContent
(
const
std::unordered_map<Particle::Role, spdgids_t>&);
155
156
double
alphaEM
(
double
q)
const
;
157
double
alphaS
(
double
q)
const
;
158
159
inline
const
std::vector<double>&
lastCoordinates
()
const
{
return
point_coord_; }
160
161
private
:
162
double
s_{-1.};
163
double
sqs_{-1.};
164
double
inv_s_{-1.};
165
double
inv_sqs_{-1.};
166
double
wcm_{-1.};
167
double
mA2_{-1.};
168
double
mB2_{-1.};
169
double
mX2_{-1.};
170
double
mY2_{-1.};
171
double
t1_{-1.};
172
double
t2_{-1.};
173
double
x1_{0.};
174
double
x2_{0.};
175
std::unique_ptr<Coupling> alpha_em_;
176
std::unique_ptr<Coupling> alpha_s_;
178
struct
MappingVariable {
179
std::string name;
180
std::string description;
181
Limits limits;
182
double
& value;
183
Mapping
type;
184
size_t
index;
185
};
186
std::vector<MappingVariable> mapped_variables_;
187
std::vector<double> point_coord_;
188
double
base_jacobian_{1.};
189
Kinematics
kin_{ParametersList()};
190
std::unique_ptr<Event> event_;
191
friend
class
utils::ProcessVariablesAnalyser
;
192
};
193
using
ProcessPtr
= std::unique_ptr<Process>;
194
}
// namespace cepgen::proc
195
196
#endif
115
virtual
void
fillKinematics
() = 0; {
…
}
113
virtual
void
addEventContent
() = 0; {
…
}
110
inline
double
inverseSqrtS
()
const
{
return
inv_sqs_; } {
…
}
109
inline
double
inverseS
()
const
{
return
inv_s_; } {
…
}
108
inline
double
sqrtS
()
const
{
return
sqs_; } {
…
}
104
Momentum
&
pc
(
size_t
); {
…
}
102
inline
double
&
x2
() {
return
x2_; } {
…
}
101
inline
double
&
t2
() {
return
t2_; } {
…
}
100
inline
double
&
x1
() {
return
x1_; } {
…
}
99
inline
double
&
t1
() {
return
t1_; } {
…
}
98
inline
double
&
mY2
() {
return
mY2_; } {
…
}
92
Momentum
&
pA
(); {
…
}
85
const
Momentum
&
q2
()
const
; {
…
}
84
inline
double
t2
()
const
{
return
t2_; } {
…
}
82
const
Momentum
&
q1
()
const
; {
…
}
81
inline
double
t1
()
const
{
return
t1_; } {
…
}
80
inline
double
x1
()
const
{
return
x1_; } {
…
}
78
inline
double
mY
()
const
{
return
std::sqrt(
mY2
()); } {
…
}
76
const
Momentum
&
pY
()
const
; {
…
}
75
inline
double
mX
()
const
{
return
std::sqrt(
mX2
()); } {
…
}
73
const
Momentum
&
pX
()
const
; {
…
}
72
inline
double
mB
()
const
{
return
std::sqrt(
mB2
()); } {
…
}
70
const
Momentum
&
pB
()
const
; {
…
}
69
inline
double
mA
()
const
{
return
std::sqrt(
mA2
()); } {
…
}
52
void
setKinematics
(); {
…
}
51
inline
Kinematics
&
kinematics
() {
return
kin_; } {
…
}
34
class
Process
:
public
NamedModule
<Process> {
…
};
Coupling.h
Event.h
Kinematics.h
ProcessVariablesAnalyser.h
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:26
cepgen::Kinematics
List of kinematic constraints to apply on the process phase space.
Definition
Kinematics.h:27
cepgen::Limits
Validity interval for a variable.
Definition
Limits.h:28
cepgen::Momentum
Container for a particle's 4-momentum, along with useful methods to ease the development of any matri...
Definition
Momentum.h:30
cepgen::NamedModule
Base runtime module object.
Definition
NamedModule.h:28
cepgen::NamedModule< Process >::name
const std::string & name() const
Module unique indexing name.
Definition
NamedModule.h:42
cepgen::ParametersDescription
A description object for parameters collection.
Definition
ParametersDescription.h:26
cepgen::ParametersList
Definition
ParametersList.h:52
cepgen::proc::Process
Class template to define any process to compute using this MC integrator/events generator.
Definition
Process.h:34
cepgen::proc::Process::clearEvent
void clearEvent()
Restore the event object to its initial state.
cepgen::proc::Process::dumpPoint
void dumpPoint(std::ostream *=nullptr) const
Dump the coordinate of the phase-space point being evaluated.
cepgen::proc::Process::prepareKinematics
virtual void prepareKinematics()
Compute the incoming state kinematics.
Definition
Process.h:114
cepgen::proc::Process::pY
Momentum & pY()
Negative-z outgoing beam particle's 4-momentum.
cepgen::proc::Process::wCM
double wCM() const
Two-parton centre of mass energy.
Definition
Process.h:90
cepgen::proc::Process::mB2
double mB2() const
Negative-z incoming beam particle's squared mass.
Definition
Process.h:71
cepgen::proc::Process::kinematics
const Kinematics & kinematics() const
Constant reference to the process kinematics.
Definition
Process.h:50
cepgen::proc::Process::pA
const Momentum & pA() const
Positive-z incoming beam particle's 4-momentum.
cepgen::proc::Process::initialise
void initialise()
Initialise the process once the kinematics has been set.
cepgen::proc::Process::pB
Momentum & pB()
Negative-z incoming beam particle's 4-momentum.
cepgen::proc::Process::pX
const Momentum & pX() const
Positive-z outgoing beam particle's 4-momentum.
cepgen::proc::Process::t1
double t1() const
Positive-z incoming parton's squared mass.
Definition
Process.h:81
cepgen::proc::Process::x2
double x2() const
Negative-z incoming parton's fractional momentum.
Definition
Process.h:83
cepgen::proc::Process::inverseS
double inverseS() const
Inverse two-beam squared centre of mass energy.
Definition
Process.h:109
cepgen::proc::Process::mp_
double mp_
Proton mass, in GeV/c .
Definition
Process.h:120
cepgen::proc::Process::clone
virtual std::unique_ptr< Process > clone() const
Copy all process attributes into a new object.
cepgen::proc::Process::mX
double mX() const
Positive-z outgoing beam particle's mass.
Definition
Process.h:75
cepgen::proc::Process::lastCoordinates
const std::vector< double > & lastCoordinates() const
Last coordinates fed.
Definition
Process.h:159
cepgen::proc::Process::q2
const Momentum & q2() const
Negative-z incoming parton's 4-momentum.
cepgen::proc::Process::inverseSqrtS
double inverseSqrtS() const
Inverse two-beam centre of mass energy.
Definition
Process.h:110
cepgen::proc::Process::addEventContent
virtual void addEventContent()=0
Set the incoming and outgoing state to be expected in the process.
cepgen::proc::Process::fillKinematics
virtual void fillKinematics()=0
Fill the Event object with the particles' kinematics.
cepgen::proc::Process::event
Event & event()
Event object read/write accessor.
cepgen::proc::Process::operator<<
friend std::ostream & operator<<(std::ostream &, const Mapping &)
Human-friendly printout of the mapping type Register a variable to be handled and populated whenever ...
cepgen::proc::Process::defineVariable
Process & defineVariable(double &out, const Mapping &type, const Limits &lim, const std::string &name, const std::string &description="")
cepgen::proc::Process::t1
double & t1()
Positive-z incoming parton's squared mass.
Definition
Process.h:99
cepgen::proc::Process::dumpVariables
void dumpVariables(std::ostream *=nullptr) const
List all variables handled by this generic process.
cepgen::proc::Process::q2
Momentum & q2()
Negative-z incoming parton's 4-momentum.
cepgen::proc::Process::q1
Momentum & q1()
Positive-z incoming parton's 4-momentum.
cepgen::proc::Process::kinematics
Kinematics & kinematics()
Reference to the process kinematics.
Definition
Process.h:51
cepgen::proc::Process::Mapping
Mapping
Type of mapping to apply on the variable.
Definition
Process.h:125
cepgen::proc::Process::Mapping::square
@ square
a square mapping
cepgen::proc::Process::Mapping::linear
@ linear
a linear mapping
cepgen::proc::Process::Mapping::exponential
@ exponential
an exponential mapping
cepgen::proc::Process::Mapping::power_law
@ power_law
a power-law mapping inherited from LPAIR
cepgen::proc::Process::description
static ParametersDescription description()
cepgen::proc::Process::s
double s() const
Two-beam squared centre of mass energy.
Definition
Process.h:107
cepgen::proc::Process::x1
double x1() const
Positive-z incoming parton's fractional momentum.
Definition
Process.h:80
cepgen::proc::Process::t2
double & t2()
Negative-z incoming parton's squared mass.
Definition
Process.h:101
cepgen::proc::Process::shat
double shat() const
cepgen::proc::Process::mA
double mA() const
Positive-z incoming beam particle's mass.
Definition
Process.h:69
cepgen::proc::Process::mA2
double mA2() const
Positive-z incoming beam particle's squared mass.
Definition
Process.h:68
cepgen::proc::Process::mY2
double mY2() const
Negative-z outgoing beam particle's squared mass.
Definition
Process.h:77
cepgen::proc::Process::x1
double & x1()
Positive-z incoming parton's fractional momentum.
Definition
Process.h:100
cepgen::proc::Process::Process
Process(const ParametersList &)
cepgen::proc::Process::sqrtS
double sqrtS() const
Two-beam centre of mass energy.
Definition
Process.h:108
cepgen::proc::Process::alphaEM
double alphaEM(double q) const
Compute the electromagnetic running coupling algorithm at a given scale.
cepgen::proc::Process::pB
const Momentum & pB() const
Negative-z incoming beam particle's 4-momentum.
cepgen::proc::Process::setEventContent
void setEventContent(const std::unordered_map< Particle::Role, spdgids_t > &)
Set the incoming and outgoing states to be defined in this process (and prepare the Event object acco...
cepgen::proc::Process::mY
double mY() const
Negative-z outgoing beam particle's mass.
Definition
Process.h:78
cepgen::proc::Process::mB
double mB() const
Negative-z incoming beam particle's mass.
Definition
Process.h:72
cepgen::proc::Process::hasEvent
bool hasEvent() const
Does the process contain (and hold) an event?
Definition
Process.h:62
cepgen::proc::Process::mX2
double & mX2()
Positive-z outgoing beam particle's squared mass.
Definition
Process.h:97
cepgen::proc::Process::pc
Momentum & pc(size_t)
Central particle's 4-momentum.
cepgen::proc::Process::weight
double weight(const std::vector< double > &)
Compute the weight for a phase-space point.
cepgen::proc::Process::eventPtr
Event * eventPtr()
Event pointer read/write accessor.
cepgen::proc::Process::pY
const Momentum & pY() const
Negative-z outgoing beam particle's 4-momentum.
cepgen::proc::Process::pA
Momentum & pA()
Positive-z incoming beam particle's 4-momentum.
cepgen::proc::Process::pc
const Momentum & pc(size_t) const
Central particle's 4-momentum.
cepgen::proc::Process::Process
Process(const Process &)
Copy constructor for a user process.
cepgen::proc::Process::computeWeight
virtual double computeWeight()=0
Compute the phase space point weight.
cepgen::proc::Process::clear
void clear()
Reset process prior to the phase space and variables definition.
cepgen::proc::Process::setKinematics
void setKinematics()
cepgen::proc::Process::mX2
double mX2() const
Positive-z outgoing beam particle's squared mass.
Definition
Process.h:74
cepgen::proc::Process::event
const Event & event() const
Handled particles objects and their relationships.
cepgen::proc::Process::mY2
double & mY2()
Negative-z outgoing beam particle's squared mass.
Definition
Process.h:98
cepgen::proc::Process::t2
double t2() const
Negative-z incoming parton's squared mass.
Definition
Process.h:84
cepgen::proc::Process::ndim
size_t ndim() const
Number of dimensions to perform integration.
Definition
Process.h:59
cepgen::proc::Process::~Process
~Process() override=default
cepgen::proc::Process::alphaS
double alphaS(double q) const
Compute the strong coupling algorithm at a given scale.
cepgen::proc::Process::generateVariables
double generateVariables() const
Generate and initialise all variables handled by this process.
cepgen::proc::Process::q1
const Momentum & q1() const
Positive-z incoming parton's 4-momentum.
cepgen::proc::Process::x2
double & x2()
Negative-z incoming parton's fractional momentum.
Definition
Process.h:102
cepgen::proc::Process::mp2_
double mp2_
Squared proton mass, in GeV /c .
Definition
Process.h:121
cepgen::proc::Process::pX
Momentum & pX()
Positive-z outgoing beam particle's 4-momentum.
cepgen::proc::Process::operator=
Process & operator=(const Process &)
Assignment operator.
cepgen::proc::Process::variableValue
double variableValue(size_t i, double x) const
Retrieve the physical value for one variable.
cepgen::utils::ProcessVariablesAnalyser
Definition
ProcessVariablesAnalyser.h:35
cepgen::mode::Kinematics
Kinematics
Type of scattering.
Definition
Modes.h:27
cepgen::proc
Location for all physics processes to be generated.
Definition
GeneratorWorker.h:31
cepgen::proc::ProcessPtr
std::unique_ptr< Process > ProcessPtr
Helper for a Process unique pointer.
Definition
Process.h:193
include
CepGen
Process
Process.h
Generated on Tue Apr 22 2025 for CepGen by
1.10.0