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
►
Event.h
►
Particle.h
►
EventFilter
►
FormFactors
►
Integration
►
Modules
►
PartonFluxes
►
Physics
►
Process
►
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
Event.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_Event_Event_h
20
#define CepGen_Event_Event_h
21
22
#include "
CepGen/Event/Particle.h
"
23
24
namespace
cepgen
{
26
class
Event
{
27
public
:
28
explicit
Event
(
bool
compressed
=
false
);
29
Event
(
const
Event
&);
30
31
Event
&
operator=
(
const
Event
&);
32
bool
operator==
(
const
Event
&)
const
;
33
36
static
Event
minimal
(
size_t
num_out_particles = 1);
37
38
void
clear
();
39
void
freeze
();
40
void
restore
();
41
bool
compressed
()
const
;
42
Event
compress
()
const
;
43
void
updateRoles
();
44
46
friend
std::ostream&
operator<<
(std::ostream&,
const
Event
&);
48
void
dump
()
const
;
50
double
cmEnergy
()
const
;
51
52
//----- particles adders
53
57
ParticleRef
addParticle
(
Particle
& particle,
bool
replace =
false
);
61
ParticleRef
addParticle
(
Particle::Role
role,
bool
replace =
false
);
62
63
//----- particles retrievers
64
65
size_t
size
()
const
;
66
bool
empty
()
const
;
67
Particles
particles
()
const
;
68
Particles
stableParticles
()
const
;
69
Particles
stableParticlesWithRole
(
Particle::Role
)
const
;
70
ParticlesMap
&
map
() {
return
particles_; }
71
74
ParticlesRefs
operator[]
(
Particle::Role
role);
76
const
Particles
&
operator()
(
Particle::Role
role)
const
;
78
ParticlesIds
ids
(
Particle::Role
role)
const
;
80
inline
bool
hasRole
(
Particle::Role
role)
const
{
return
particles_.count(role) != 0; }
83
Particle
&
oneWithRole
(
Particle::Role
role);
85
const
Particle
&
oneWithRole
(
Particle::Role
role)
const
;
88
Particle
&
operator[]
(
int
id
);
91
const
Particle
&
operator()
(
int
id
)
const
;
94
ParticlesRefs
operator[]
(
const
ParticlesIds
&
ids
);
97
Particles
operator()
(
const
ParticlesIds
&
ids
)
const
;
98
100
Momentum
missingMomentum
()
const
;
101
102
//----- general particles information retriever
103
106
Particles
mothers
(
const
Particle
& particle)
const
;
109
ParticlesRefs
mothers
(
const
Particle
& particle);
111
void
clearMothers
(
Particle
& particle);
114
Particles
children
(
const
Particle
& particle)
const
;
117
ParticlesRefs
children
(
const
Particle
& particle);
119
Particles
stableChildren
(
const
Particle
& particle,
bool
recursive =
false
)
const
;
121
ParticlesRefs
stableChildren
(
const
Particle
& particle,
bool
recursive =
false
);
123
void
clearChildren
(
Particle
& particle);
125
ParticleRoles
roles
()
const
;
126
128
struct
EventMetadata
: std::unordered_map<std::string, float> {
129
EventMetadata
();
131
float
operator()
(
const
std::string& key)
const
{
return
count(key) > 0 ? at(key) : -1.; }
132
};
128
struct
EventMetadata
: std::unordered_map<std::string, float> {
…
};
134
EventMetadata
metadata
;
135
136
private
:
137
static
constexpr
double
MIN_PRECISION = 1.e-10;
138
void
checkKinematics()
const
;
139
ParticlesMap
particles_;
141
struct
NumParticles {
142
size_t
cs{0};
143
size_t
op1{0};
144
size_t
op2{0};
145
} event_content_{};
146
bool
compressed_{
false
};
147
};
148
}
// namespace cepgen
149
150
#endif
74
ParticlesRefs
operator[]
(
Particle::Role
role); {
…
}
26
class
Event
{
…
};
Particle.h
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:26
cepgen::Event::stableParticles
Particles stableParticles() const
Vector of all stable particles in the event.
cepgen::Event::clearMothers
void clearMothers(Particle &particle)
Remove all mothers from a given particle (also affects the mothers' filiation)
cepgen::Event::compress
Event compress() const
Compress the event record.
cepgen::Event::map
ParticlesMap & map()
Internal particles map retrieval operator.
Definition
Event.h:70
cepgen::Event::Event
Event(const Event &)
Copy constructor.
cepgen::Event::children
Particles children(const Particle &particle) const
List of all the daughters from a particle.
cepgen::Event::operator()
const Particles & operator()(Particle::Role role) const
Get a list of constant Particle objects corresponding to a certain role in the process kinematics.
cepgen::Event::stableParticlesWithRole
Particles stableParticlesWithRole(Particle::Role) const
Vector of all stable particles with a given role.
cepgen::Event::stableChildren
ParticlesRefs stableChildren(const Particle &particle, bool recursive=false)
List all the stable daughters of a particle in this event.
cepgen::Event::updateRoles
void updateRoles()
Update the table of particles roles.
cepgen::Event::mothers
ParticlesRefs mothers(const Particle &particle)
List of all parent Particle object for this given particle.
cepgen::Event::oneWithRole
Particle & oneWithRole(Particle::Role role)
First Particle object with a given role in the event.
cepgen::Event::operator<<
friend std::ostream & operator<<(std::ostream &, const Event &)
Human-readable version of the event content.
cepgen::Event::size
size_t size() const
Number of particles in the event.
cepgen::Event::operator()
const Particle & operator()(int id) const
Constant reference to the Particle object corresponding to a unique identifier in the event.
cepgen::Event::cmEnergy
double cmEnergy() const
Incoming beams centre-of-mass energy, in GeV.
cepgen::Event::roles
ParticleRoles roles() const
List of roles defined for the given event (really process-dependant for the central system)
cepgen::Event::ids
ParticlesIds ids(Particle::Role role) const
Get a list of particle identifiers in Event corresponding to a certain role in the process kinematics...
cepgen::Event::operator()
Particles operator()(const ParticlesIds &ids) const
Particle objects corresponding to the unique identifiers in the event.
cepgen::Event::empty
bool empty() const
Is the particles map empty?
cepgen::Event::missingMomentum
Momentum missingMomentum() const
Compute the missing momentum for central particles in this event.
cepgen::Event::children
ParticlesRefs children(const Particle &particle)
List of all the daughters from a particle.
cepgen::Event::particles
Particles particles() const
Vector of all particles in the event.
cepgen::Event::addParticle
ParticleRef addParticle(Particle::Role role, bool replace=false)
Create a new particle in the event, with no kinematic information but the role it has to play in the ...
cepgen::Event::operator==
bool operator==(const Event &) const
Equality operator.
cepgen::Event::metadata
EventMetadata metadata
List of auxiliary information.
Definition
Event.h:134
cepgen::Event::mothers
Particles mothers(const Particle &particle) const
List of all parent Particle object for this given particle.
cepgen::Event::Event
Event(bool compressed=false)
Build an empty event.
cepgen::Event::oneWithRole
const Particle & oneWithRole(Particle::Role role) const
First constant Particle object with a given role in the event.
cepgen::Event::compressed
bool compressed() const
Is the event already without intermediate-channel information?
cepgen::Event::clearChildren
void clearChildren(Particle &particle)
Remove all daughters from a given particle (also affects the daughters' parentage)
cepgen::Event::addParticle
ParticleRef addParticle(Particle &particle, bool replace=false)
Set the information on one particle in the process.
cepgen::Event::minimal
static Event minimal(size_t num_out_particles=1)
Build a trivial event with the minimal information.
cepgen::Event::operator[]
Particle & operator[](int id)
Reference to the Particle object corresponding to a unique identifier in the event.
cepgen::Event::hasRole
bool hasRole(Particle::Role role) const
Check whether a particle role is represented in this event.
Definition
Event.h:80
cepgen::Event::freeze
void freeze()
Store a snapshot of the primordial event block.
cepgen::Event::clear
void clear()
Empty the whole event content.
cepgen::Event::operator[]
ParticlesRefs operator[](const ParticlesIds &ids)
References to the Particle objects corresponding to the unique identifiers in the event.
cepgen::Event::stableChildren
Particles stableChildren(const Particle &particle, bool recursive=false) const
List all the stable daughters of a particle in this event.
cepgen::Event::operator[]
ParticlesRefs operator[](Particle::Role role)
List of references to Particle objects corresponding to a certain role in the process kinematics.
cepgen::Event::operator=
Event & operator=(const Event &)
Assignment operator.
cepgen::Event::dump
void dump() const
Dump all the known information on every Particle object contained in this Event container in the outp...
cepgen::Event::restore
void restore()
Restore the event to its "empty" state.
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::Particle
Kinematic information for one particle.
Definition
Particle.h:32
cepgen::Particle::Role
Role
Definition
Particle.h:49
cepgen::ParticlesMap
Map between a particle's role and its associated Particle object.
Definition
Particle.h:173
cepgen
Common namespace for this Monte Carlo generator.
Definition
Handler.h:26
cepgen::ParticlesIds
std::set< int > ParticlesIds
A set of integer-type particle identifiers.
Definition
Particle.h:29
cepgen::ParticleRoles
std::vector< Particle::Role > ParticleRoles
List of particles' roles.
Definition
Particle.h:170
cepgen::ParticlesRefs
std::vector< ParticleRef > ParticlesRefs
List of references to Particle objects.
Definition
Particle.h:169
cepgen::Particles
std::vector< Particle > Particles
List of Particle objects.
Definition
Particle.h:168
cepgen::ParticleRef
std::reference_wrapper< Particle > ParticleRef
Reference to a Particle object.
Definition
Particle.h:167
cepgen::Event::EventMetadata
Collection of key -> value pairs storing event metadata.
Definition
Event.h:128
cepgen::Event::EventMetadata::operator()
float operator()(const std::string &key) const
Retrieve the metadata value associated with a key.
Definition
Event.h:131
cepgen::Event::EventMetadata::EventMetadata
EventMetadata()
include
CepGen
Event
Event.h
Generated on Tue Apr 22 2025 for CepGen by
1.10.0