cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
Cuts.h
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2017-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_Physics_Cuts_h
20#define CepGen_Physics_Cuts_h
21
23
24namespace cepgen {
25 class Event;
26 class ParametersList;
27
29 namespace cuts {
31 struct Central final : public SteeredObject<Central> {
32 Central();
33 explicit Central(const ParametersList&);
34
36 bool contain(const Particles&, const Event* evt = nullptr) const;
37
50 };
51
53 struct Initial final : public SteeredObject<Initial> {
54 explicit Initial(const ParametersList&);
55
57 bool contain(const Particles&, const Event* evt = nullptr) const;
58
59 void setParameters(const ParametersList&) override;
60
61 std::vector<Limits> q2;
64 };
65
67 struct Remnants final : public SteeredObject<Remnants> {
68 explicit Remnants(const ParametersList&);
69
71 bool contain(const Particles&, const Event* evt = nullptr) const;
72
73 void setParameters(const ParametersList&) override;
74
78
79 static constexpr double MX_MIN = 1.07;
80 };
81 } // namespace cuts
83 typedef std::unordered_map<pdgid_t, cuts::Central> PerIdCuts;
85 struct CutsList final : SteeredObject<CutsList> {
86 explicit CutsList(const ParametersList&);
87
89
90 void setParameters(const ParametersList&) override;
91
92 friend std::ostream& operator<<(std::ostream&, const CutsList&);
93
98 };
99} // namespace cepgen
100
101#endif
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Validity interval for a variable.
Definition Limits.h:28
A description object for parameters collection.
Base user-steerable object.
Common namespace for this Monte Carlo generator.
std::unordered_map< pdgid_t, cuts::Central > PerIdCuts
Collection of cuts to be applied on all particle with a given PDG id.
Definition Cuts.h:83
std::vector< Particle > Particles
List of Particle objects.
Definition Particle.h:169
A collection of cuts to apply on the physical phase space.
Definition Cuts.h:85
void setParameters(const ParametersList &) override
Set module parameters.
Definition Cuts.cpp:187
cuts::Initial initial
Cuts on the initial particles kinematics.
Definition Cuts.h:94
cuts::Remnants remnants
Cuts on the beam remnants system.
Definition Cuts.h:97
cuts::Central central
Cuts on the central system produced.
Definition Cuts.h:95
static ParametersDescription description()
Definition Cuts.cpp:224
friend std::ostream & operator<<(std::ostream &, const CutsList &)
Definition Cuts.cpp:232
PerIdCuts central_particles
Cuts on the central individual particles.
Definition Cuts.h:96
Centrally produced particles phase space cuts.
Definition Cuts.h:31
Limits energy_sum
multiparticle system energy
Definition Cuts.h:45
Limits eta_single
single particle pseudo-rapidity
Definition Cuts.h:39
Limits energy_single
single particle energy
Definition Cuts.h:41
Limits mass_single
single particle mass
Definition Cuts.h:42
Limits pt_diff
transverse momentum balance between the central particles
Definition Cuts.h:47
Limits rapidity_single
single particle rapidity
Definition Cuts.h:40
Limits pt_sum
multiparticle system transverse momentum
Definition Cuts.h:43
Limits phi_diff
azimuthal angles difference between the central particles
Definition Cuts.h:48
Limits mass_sum
multiparticle system invariant mass
Definition Cuts.h:46
bool contain(const Particles &, const Event *evt=nullptr) const
Definition Cuts.cpp:55
Limits eta_sum
multiparticle system pseudo-rapidity
Definition Cuts.h:44
Limits pt_single
single particle transverse momentum
Definition Cuts.h:38
static ParametersDescription description()
Definition Cuts.cpp:77
Limits rapidity_diff
rapidity balance between the central particles
Definition Cuts.h:49
Initial parton-like particles phase space cuts.
Definition Cuts.h:53
Limits phi
parton azimuthal angle
Definition Cuts.h:63
void setParameters(const ParametersList &) override
Set module parameters.
Definition Cuts.cpp:100
bool contain(const Particles &, const Event *evt=nullptr) const
Definition Cuts.cpp:114
std::vector< Limits > q2
parton virtualities
Definition Cuts.h:61
static ParametersDescription description()
Definition Cuts.cpp:130
Limits qt
parton transverse virtuality
Definition Cuts.h:62
Outgoing beam remnant-like particles phase space cuts.
Definition Cuts.h:67
Limits mx
diffractive mass
Definition Cuts.h:75
void setParameters(const ParametersList &) override
Set module parameters.
Definition Cuts.cpp:144
bool contain(const Particles &, const Event *evt=nullptr) const
Definition Cuts.cpp:155
Limits xi
longitudinal momentum fraction
Definition Cuts.h:77
static ParametersDescription description()
Definition Cuts.cpp:168
Limits yj
diffractive jet rapidity
Definition Cuts.h:76
static constexpr double MX_MIN
Minimal diffractive mass for dissociative proton treatment.
Definition Cuts.h:79