cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
1.2.3
A generic central exclusive processes event generator
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-2024 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 <memory>
23
24
#include "
CepGen/Event/Particle.h
"
25
26
namespace
cepgen
{
28
class
Event
{
29
public
:
30
explicit
Event
(
bool
compressed
=
false
);
31
Event
(
const
Event
&);
32
33
Event
&
operator=
(
const
Event
&);
34
37
static
Event
minimal
(
size_t
num_out_particles = 1);
38
39
void
clear
();
40
void
freeze
();
41
void
restore
();
42
bool
compressed
()
const
;
43
Event
compress
()
const
;
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
& part,
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
ParticlesMap
&
map
() {
return
particles_; }
70
73
ParticlesRefs
operator[]
(
Particle::Role
role);
75
const
Particles
&
operator()
(
Particle::Role
role)
const
;
77
ParticlesIds
ids
(
Particle::Role
role)
const
;
79
inline
bool
hasRole
(
Particle::Role
role)
const
{
return
particles_.count(role) != 0; }
82
Particle
&
oneWithRole
(
Particle::Role
role);
84
const
Particle
&
oneWithRole
(
Particle::Role
role)
const
;
87
Particle
&
operator[]
(
int
id
);
90
const
Particle
&
operator()
(
int
id
)
const
;
93
ParticlesRefs
operator[]
(
const
ParticlesIds
&
ids
);
96
Particles
operator()
(
const
ParticlesIds
&
ids
)
const
;
97
99
Momentum
missingMomentum
()
const
;
100
101
//----- general particles information retriever
102
105
Particles
mothers
(
const
Particle
& part)
const
;
108
ParticlesRefs
mothers
(
const
Particle
& part);
111
Particles
daughters
(
const
Particle
& part)
const
;
114
ParticlesRefs
daughters
(
const
Particle
& part);
116
Particles
stableDaughters
(
const
Particle
& part,
bool
recursive =
false
)
const
;
118
ParticlesRefs
stableDaughters
(
const
Particle
& part,
bool
recursive =
false
);
120
ParticleRoles
roles
()
const
;
121
123
struct
EventMetadata
: std::unordered_map<std::string, float> {
124
EventMetadata
();
126
float
operator()
(
const
std::string& key)
const
{
return
count(key) > 0 ? at(key) : -1.; }
127
};
129
EventMetadata
metadata
;
130
131
private
:
132
static
constexpr
double
MIN_PRECISION = 1.e-10;
133
void
checkKinematics()
const
;
134
ParticlesMap
particles_;
136
struct
NumParticles {
137
size_t
cs{0};
138
size_t
op1{0};
139
size_t
op2{0};
140
} evtcontent_{};
141
bool
compressed_{
false
};
142
};
143
}
// namespace cepgen
144
145
#endif
Particle.h
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:28
cepgen::Event::stableParticles
Particles stableParticles() const
Vector of all stable particles in the event.
Definition
Event.cpp:320
cepgen::Event::compress
Event compress() const
Compress the event record.
Definition
Event.cpp:105
cepgen::Event::map
ParticlesMap & map()
Internal particles map retrieval operator.
Definition
Event.h:69
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.
Definition
Event.cpp:162
cepgen::Event::oneWithRole
Particle & oneWithRole(Particle::Role role)
First Particle object with a given role in the event.
Definition
Event.cpp:181
cepgen::Event::operator<<
friend std::ostream & operator<<(std::ostream &, const Event &)
Human-readable version of the event content.
Definition
Event.cpp:372
cepgen::Event::size
size_t size() const
Number of particles in the event.
Definition
Event.cpp:303
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)
Definition
Event.cpp:274
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...
Definition
Event.cpp:169
cepgen::Event::empty
bool empty() const
Is the particles map empty?
Definition
Event.cpp:309
cepgen::Event::missingMomentum
Momentum missingMomentum() const
Compute the missing momentum for central particles in this event.
Definition
Event.cpp:331
cepgen::Event::daughters
Particles daughters(const Particle &part) const
List of all the daughters from a particle.
Definition
Event.cpp:242
cepgen::Event::addParticle
ParticleRef addParticle(Particle &part, bool replace=false)
Set the information on one particle in the process.
Definition
Event.cpp:281
cepgen::Event::particles
Particles particles() const
Vector of all particles in the event.
Definition
Event.cpp:311
cepgen::Event::metadata
EventMetadata metadata
List of auxiliary information.
Definition
Event.h:129
cepgen::Event::compressed
bool compressed() const
Is the event already without intermediate-channel information?
Definition
Event.cpp:103
cepgen::Event::stableDaughters
Particles stableDaughters(const Particle &part, bool recursive=false) const
List all the stable daughters of a particle in this event.
Definition
Event.cpp:246
cepgen::Event::hasRole
bool hasRole(Particle::Role role) const
Check whether a particle role is represented in this event.
Definition
Event.h:79
cepgen::Event::freeze
void freeze()
Store a snapshot of the primordial event block.
Definition
Event.cpp:85
cepgen::Event::mothers
Particles mothers(const Particle &part) const
List of all parent Particle object for this given particle.
Definition
Event.cpp:238
cepgen::Event::clear
void clear()
Empty the whole event content.
Definition
Event.cpp:80
cepgen::Event::operator[]
ParticlesRefs operator[](Particle::Role role)
List of references to Particle objects corresponding to a certain role in the process kinematics.
Definition
Event.cpp:154
cepgen::Event::operator=
Event & operator=(const Event &)
Assignment operator.
Definition
Event.cpp:33
cepgen::Event::dump
void dump() const
Dump all the known information on every Particle object contained in this Event container in the outp...
Definition
Event.cpp:370
cepgen::Event::minimal
static Event minimal(size_t num_out_particles=1)
Build a trivial event with the minimal information.
Definition
Event.cpp:42
cepgen::Event::restore
void restore()
Restore the event to its "empty" state.
Definition
Event.cpp:94
cepgen::Momentum
Container for a particle's 4-momentum, along with useful methods to ease the development of any matri...
Definition
Momentum.h:33
cepgen::Particle
Kinematic information for one particle.
Definition
Particle.h:33
cepgen::Particle::Role
Role
Definition
Particle.h:50
cepgen::ParticlesMap
Map between a particle's role and its associated Particle object.
Definition
Particle.h:177
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
cepgen::ParticleRoles
std::vector< Particle::Role > ParticleRoles
List of particles' roles.
Definition
Particle.h:174
cepgen::ParticleRef
std::reference_wrapper< Particle > ParticleRef
Reference to a Particle object.
Definition
Particle.h:171
cepgen::ParticlesIds
std::set< int > ParticlesIds
A set of integer-type particle identifiers.
Definition
Particle.h:30
cepgen::ParticlesRefs
std::vector< ParticleRef > ParticlesRefs
List of references to Particle objects.
Definition
Particle.h:173
cepgen::Particles
std::vector< Particle > Particles
List of Particle objects.
Definition
Particle.h:172
cepgen::Event::EventMetadata
Collection of key -> value pairs storing event metadata.
Definition
Event.h:123
cepgen::Event::EventMetadata::operator()
float operator()(const std::string &key) const
Retrieve the metadata value associated to a key.
Definition
Event.h:126
cepgen::Event::EventMetadata::EventMetadata
EventMetadata()
Definition
Event.cpp:462
CepGen
Event
Event.h
Generated on Thu Apr 25 2024 for CepGen by
1.9.7