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
22
#include "
CepGen/Event/Particle.h
"
23
24
namespace
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
35
static
ParametersDescription
description
();
36
bool
contain
(
const
Particles
&,
const
Event
* evt =
nullptr
)
const
;
37
38
Limits
pt_single
;
39
Limits
eta_single
;
40
Limits
rapidity_single
;
41
Limits
energy_single
;
42
Limits
mass_single
;
43
Limits
pt_sum
;
44
Limits
eta_sum
;
45
Limits
energy_sum
;
46
Limits
mass_sum
;
47
Limits
pt_diff
;
48
Limits
phi_diff
;
49
Limits
rapidity_diff
;
50
};
51
53
struct
Initial
final :
public
SteeredObject
<Initial> {
54
explicit
Initial
(
const
ParametersList
&);
55
56
static
ParametersDescription
description
();
57
bool
contain
(
const
Particles
&,
const
Event
* evt =
nullptr
)
const
;
58
59
void
setParameters
(
const
ParametersList
&)
override
;
60
61
std::vector<Limits>
q2
;
62
Limits
qt
;
63
Limits
phi
;
64
};
65
67
struct
Remnants
final :
public
SteeredObject
<Remnants> {
68
explicit
Remnants
(
const
ParametersList
&);
69
70
static
ParametersDescription
description
();
71
bool
contain
(
const
Particles
&,
const
Event
* evt =
nullptr
)
const
;
72
73
void
setParameters
(
const
ParametersList
&)
override
;
74
75
Limits
mx
;
76
Limits
yj
;
77
Limits
xi
;
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
88
static
ParametersDescription
description
();
89
90
void
setParameters
(
const
ParametersList
&)
override
;
91
92
friend
std::ostream&
operator<<
(std::ostream&,
const
CutsList
&);
93
94
cuts::Initial
initial
;
95
cuts::Central
central
;
96
PerIdCuts
central_particles
;
97
cuts::Remnants
remnants
;
98
};
99
}
// namespace cepgen
100
101
#endif
Particle.h
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:28
cepgen::Limits
Validity interval for a variable.
Definition
Limits.h:28
cepgen::ParametersDescription
A description object for parameters collection.
Definition
ParametersDescription.h:26
cepgen::ParametersList
Definition
ParametersList.h:52
cepgen::SteeredObject
Base user-steerable object.
Definition
SteeredObject.h:41
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
cepgen::PerIdCuts
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
cepgen::Particles
std::vector< Particle > Particles
List of Particle objects.
Definition
Particle.h:169
cepgen::CutsList
A collection of cuts to apply on the physical phase space.
Definition
Cuts.h:85
cepgen::CutsList::setParameters
void setParameters(const ParametersList &) override
Set module parameters.
Definition
Cuts.cpp:187
cepgen::CutsList::initial
cuts::Initial initial
Cuts on the initial particles kinematics.
Definition
Cuts.h:94
cepgen::CutsList::remnants
cuts::Remnants remnants
Cuts on the beam remnants system.
Definition
Cuts.h:97
cepgen::CutsList::central
cuts::Central central
Cuts on the central system produced.
Definition
Cuts.h:95
cepgen::CutsList::description
static ParametersDescription description()
Definition
Cuts.cpp:224
cepgen::CutsList::operator<<
friend std::ostream & operator<<(std::ostream &, const CutsList &)
Definition
Cuts.cpp:232
cepgen::CutsList::central_particles
PerIdCuts central_particles
Cuts on the central individual particles.
Definition
Cuts.h:96
cepgen::cuts::Central
Centrally produced particles phase space cuts.
Definition
Cuts.h:31
cepgen::cuts::Central::energy_sum
Limits energy_sum
multiparticle system energy
Definition
Cuts.h:45
cepgen::cuts::Central::eta_single
Limits eta_single
single particle pseudo-rapidity
Definition
Cuts.h:39
cepgen::cuts::Central::energy_single
Limits energy_single
single particle energy
Definition
Cuts.h:41
cepgen::cuts::Central::mass_single
Limits mass_single
single particle mass
Definition
Cuts.h:42
cepgen::cuts::Central::pt_diff
Limits pt_diff
transverse momentum balance between the central particles
Definition
Cuts.h:47
cepgen::cuts::Central::rapidity_single
Limits rapidity_single
single particle rapidity
Definition
Cuts.h:40
cepgen::cuts::Central::pt_sum
Limits pt_sum
multiparticle system transverse momentum
Definition
Cuts.h:43
cepgen::cuts::Central::phi_diff
Limits phi_diff
azimuthal angles difference between the central particles
Definition
Cuts.h:48
cepgen::cuts::Central::mass_sum
Limits mass_sum
multiparticle system invariant mass
Definition
Cuts.h:46
cepgen::cuts::Central::contain
bool contain(const Particles &, const Event *evt=nullptr) const
Definition
Cuts.cpp:55
cepgen::cuts::Central::eta_sum
Limits eta_sum
multiparticle system pseudo-rapidity
Definition
Cuts.h:44
cepgen::cuts::Central::Central
Central()
Definition
Cuts.cpp:32
cepgen::cuts::Central::pt_single
Limits pt_single
single particle transverse momentum
Definition
Cuts.h:38
cepgen::cuts::Central::description
static ParametersDescription description()
Definition
Cuts.cpp:77
cepgen::cuts::Central::rapidity_diff
Limits rapidity_diff
rapidity balance between the central particles
Definition
Cuts.h:49
cepgen::cuts::Initial
Initial parton-like particles phase space cuts.
Definition
Cuts.h:53
cepgen::cuts::Initial::phi
Limits phi
parton azimuthal angle
Definition
Cuts.h:63
cepgen::cuts::Initial::setParameters
void setParameters(const ParametersList &) override
Set module parameters.
Definition
Cuts.cpp:100
cepgen::cuts::Initial::contain
bool contain(const Particles &, const Event *evt=nullptr) const
Definition
Cuts.cpp:114
cepgen::cuts::Initial::q2
std::vector< Limits > q2
parton virtualities
Definition
Cuts.h:61
cepgen::cuts::Initial::description
static ParametersDescription description()
Definition
Cuts.cpp:130
cepgen::cuts::Initial::qt
Limits qt
parton transverse virtuality
Definition
Cuts.h:62
cepgen::cuts::Remnants
Outgoing beam remnant-like particles phase space cuts.
Definition
Cuts.h:67
cepgen::cuts::Remnants::mx
Limits mx
diffractive mass
Definition
Cuts.h:75
cepgen::cuts::Remnants::setParameters
void setParameters(const ParametersList &) override
Set module parameters.
Definition
Cuts.cpp:144
cepgen::cuts::Remnants::contain
bool contain(const Particles &, const Event *evt=nullptr) const
Definition
Cuts.cpp:155
cepgen::cuts::Remnants::xi
Limits xi
longitudinal momentum fraction
Definition
Cuts.h:77
cepgen::cuts::Remnants::description
static ParametersDescription description()
Definition
Cuts.cpp:168
cepgen::cuts::Remnants::yj
Limits yj
diffractive jet rapidity
Definition
Cuts.h:76
cepgen::cuts::Remnants::MX_MIN
static constexpr double MX_MIN
Minimal diffractive mass for dissociative proton treatment.
Definition
Cuts.h:79
CepGen
Physics
Cuts.h
Generated on Mon Jul 29 2024 for CepGen by
1.9.7