cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
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
25
26namespace cepgen {
28 class Event {
29 public:
30 explicit Event(bool compressed = false);
31 Event(const Event&);
32
33 Event& operator=(const Event&);
34 bool operator==(const Event&) const;
35
38 static Event minimal(size_t num_out_particles = 1);
39
40 void clear();
41 void freeze();
42 void restore();
43 bool compressed() const;
44 Event compress() const;
45 void updateRoles();
46
48 friend std::ostream& operator<<(std::ostream&, const Event&);
50 void dump() const;
52 double cmEnergy() const;
53
54 //----- particles adders
55
59 ParticleRef addParticle(Particle& part, bool replace = false);
63 ParticleRef addParticle(Particle::Role role, bool replace = false);
64
65 //----- particles retrievers
66
67 size_t size() const;
68 bool empty() const;
69 Particles particles() const;
72 ParticlesMap& map() { return particles_; }
73
78 const Particles& operator()(Particle::Role role) const;
82 inline bool hasRole(Particle::Role role) const { return particles_.count(role) != 0; }
87 const Particle& oneWithRole(Particle::Role role) const;
90 Particle& operator[](int id);
93 const Particle& operator()(int id) const;
100
103
104 //----- general particles information retriever
105
108 Particles mothers(const Particle& part) const;
111 ParticlesRefs mothers(const Particle& part);
113 void clearMothers(Particle& part);
116 Particles daughters(const Particle& part) const;
119 ParticlesRefs daughters(const Particle& part);
121 Particles stableDaughters(const Particle& part, bool recursive = false) const;
123 ParticlesRefs stableDaughters(const Particle& part, bool recursive = false);
125 void clearDaughters(Particle& part);
127 ParticleRoles roles() const;
128
130 struct EventMetadata : std::unordered_map<std::string, float> {
133 float operator()(const std::string& key) const { return count(key) > 0 ? at(key) : -1.; }
134 };
137
138 private:
139 static constexpr double MIN_PRECISION = 1.e-10;
140 void checkKinematics() const;
141 ParticlesMap particles_;
143 struct NumParticles {
144 size_t cs{0};
145 size_t op1{0};
146 size_t op2{0};
147 } evtcontent_{};
148 bool compressed_{false};
149 };
150} // namespace cepgen
151
152#endif
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Particles stableParticles() const
Vector of all stable particles in the event.
Definition Event.cpp:362
Event compress() const
Compress the event record.
Definition Event.cpp:115
ParticlesMap & map()
Internal particles map retrieval operator.
Definition Event.h:72
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:171
Particles stableParticlesWithRole(Particle::Role) const
Vector of all stable particles of a given role.
Definition Event.cpp:373
void updateRoles()
Update the table of particles roles.
Definition Event.cpp:310
Particle & oneWithRole(Particle::Role role)
First Particle object with a given role in the event.
Definition Event.cpp:190
friend std::ostream & operator<<(std::ostream &, const Event &)
Human-readable version of the event content.
Definition Event.cpp:431
size_t size() const
Number of particles in the event.
Definition Event.cpp:345
double cmEnergy() const
Incoming beams centre-of-mass energy, in GeV.
Definition Event.cpp:425
ParticleRoles roles() const
List of roles defined for the given event (really process-dependant for the central system)
Definition Event.cpp:303
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:178
void clearMothers(Particle &part)
Remove all mothers from a given particle (also affects the mothers' filiation)
Definition Event.cpp:251
bool empty() const
Is the particles map empty?
Definition Event.cpp:351
void clearDaughters(Particle &part)
Remove all daughters from a given particle (also affects the daughters' parentage)
Definition Event.cpp:297
Momentum missingMomentum() const
Compute the missing momentum for central particles in this event.
Definition Event.cpp:384
Particles daughters(const Particle &part) const
List of all the daughters from a particle.
Definition Event.cpp:257
ParticleRef addParticle(Particle &part, bool replace=false)
Set the information on one particle in the process.
Definition Event.cpp:323
Particles particles() const
Vector of all particles in the event.
Definition Event.cpp:353
bool operator==(const Event &) const
Equality operator.
Definition Event.cpp:42
EventMetadata metadata
List of auxiliary information.
Definition Event.h:136
bool compressed() const
Is the event already without intermediate-channel information?
Definition Event.cpp:113
Particles stableDaughters(const Particle &part, bool recursive=false) const
List all the stable daughters of a particle in this event.
Definition Event.cpp:261
bool hasRole(Particle::Role role) const
Check whether a particle role is represented in this event.
Definition Event.h:82
void freeze()
Store a snapshot of the primordial event block.
Definition Event.cpp:95
Particles mothers(const Particle &part) const
List of all parent Particle object for this given particle.
Definition Event.cpp:247
void clear()
Empty the whole event content.
Definition Event.cpp:90
ParticlesRefs operator[](Particle::Role role)
List of references to Particle objects corresponding to a certain role in the process kinematics.
Definition Event.cpp:163
Event & operator=(const Event &)
Assignment operator.
Definition Event.cpp:33
void dump() const
Dump all the known information on every Particle object contained in this Event container in the outp...
Definition Event.cpp:423
static Event minimal(size_t num_out_particles=1)
Build a trivial event with the minimal information.
Definition Event.cpp:52
void restore()
Restore the event to its "empty" state.
Definition Event.cpp:104
Container for a particle's 4-momentum, along with useful methods to ease the development of any matri...
Definition Momentum.h:33
Kinematic information for one particle.
Definition Particle.h:33
Map between a particle's role and its associated Particle object.
Definition Particle.h:174
Common namespace for this Monte Carlo generator.
std::vector< Particle::Role > ParticleRoles
List of particles' roles.
Definition Particle.h:171
std::reference_wrapper< Particle > ParticleRef
Reference to a Particle object.
Definition Particle.h:168
std::set< int > ParticlesIds
A set of integer-type particle identifiers.
Definition Particle.h:30
std::vector< ParticleRef > ParticlesRefs
List of references to Particle objects.
Definition Particle.h:170
std::vector< Particle > Particles
List of Particle objects.
Definition Particle.h:169
Collection of key -> value pairs storing event metadata.
Definition Event.h:130
float operator()(const std::string &key) const
Retrieve the metadata value associated to a key.
Definition Event.h:133