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
▼
EventFilter
►
EventBrowser.h
►
EventExporter.h
►
EventHandler.h
►
EventHarvester.h
►
EventImporter.h
►
EventModifier.h
►
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
EventBrowser.h
Go to the documentation of this file.
1
/*
2
* CepGen: a central exclusive processes event generator
3
* Copyright (C) 2019-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_EventFilter_EventBrowser_h
20
#define CepGen_EventFilter_EventBrowser_h
21
22
#include <regex>
23
24
#include "
CepGen/Event/Particle.h
"
25
26
namespace
cepgen
{
27
class
Event;
28
}
29
namespace
cepgen::utils
{
33
class
EventBrowser
{
34
public
:
35
EventBrowser
() =
default
;
36
double
get
(
const
Event
& event,
const
std::string& variable_name)
const
;
37
38
private
:
39
double
variable(
const
Event
&,
const
Particle
&,
const
std::string&)
const
;
40
double
variable(
const
Event
&,
41
const
Particle
&,
42
const
Particle
&,
43
const
std::string&)
const
;
44
static
double
variable(
const
Event
&,
const
std::string&);
45
46
static
const
std::regex rgx_select_id_;
47
static
const
std::regex rgx_select_id2_;
48
static
const
std::regex rgx_select_role_;
49
static
const
std::regex rgx_select_role2_;
50
static
constexpr
double
INVALID_OUTPUT = -999.;
51
52
//--- auxiliary helper maps
53
const
std::unordered_map<std::string, Particle::Role> role_str_ = {{
"ib1"
,
Particle::Role::IncomingBeam1
},
54
{
"ib2"
,
Particle::Role::IncomingBeam2
},
55
{
"ob1"
,
Particle::Role::OutgoingBeam1
},
56
{
"ob2"
,
Particle::Role::OutgoingBeam2
},
57
{
"pa1"
,
Particle::Role::Parton1
},
58
{
"pa2"
,
Particle::Role::Parton2
},
59
{
"cs"
,
Particle::Role::CentralSystem
},
60
{
"int"
,
Particle::Role::Intermediate
}};
61
typedef
double (Momentum::*pMethod)() const;
63
const
std::unordered_map<std::string, pMethod> m_mom_str_ = {
64
{
"px"
, &
Momentum::px
}, {
"py"
, &
Momentum::py
}, {
"pz"
, &
Momentum::pz
},
65
{
"pt"
, &
Momentum::pt
}, {
"pt2"
, &
Momentum::pt2
}, {
"eta"
, &
Momentum::eta
},
66
{
"phi"
, &
Momentum::phi
}, {
"m"
, &
Momentum::mass
}, {
"m2"
, &
Momentum::mass2
},
67
{
"mt"
, &
Momentum::massT
}, {
"mt2"
, &
Momentum::massT2
}, {
"e"
, &
Momentum::energy
},
68
{
"e2"
, &
Momentum::energy2
}, {
"et"
, &
Momentum::energyT
}, {
"et2"
, &
Momentum::energyT2
},
69
{
"p"
, &
Momentum::p
}, {
"p2"
, &
Momentum::p2
}, {
"th"
, &
Momentum::theta
},
70
{
"y"
, &
Momentum::rapidity
}, {
"beta"
, &
Momentum::beta
}, {
"gamma"
, &
Momentum::gamma
},
71
{
"gamma2"
, &
Momentum::gamma2
}};
72
typedef
double (Momentum::*pMethodOth)(const Momentum&) const;
73
const
std::unordered_map<std::string, pMethodOth> m_two_mom_str_ = {{
"deta"
, &
Momentum::deltaEta
},
74
{
"dphi"
, &
Momentum::deltaPhi
},
75
{
"dpt"
, &
Momentum::deltaPt
},
76
{
"dr"
, &
Momentum::deltaR
}};
77
};
33
class
EventBrowser
{
…
};
78
}
// namespace cepgen::utils
79
80
#endif
Particle.h
cepgen::Event
Container for the information on the in- and outgoing particles' kinematics.
Definition
Event.h:26
cepgen::Momentum::pt2
double pt2() const
Squared transverse momentum (in GeV )
cepgen::Momentum::pz
double pz() const
Longitudinal momentum (in GeV)
Definition
Momentum.h:99
cepgen::Momentum::deltaR
double deltaR(const Momentum &) const
Angular distance between two momenta.
cepgen::Momentum::deltaEta
double deltaEta(const Momentum &) const
Pseudo-rapidity distance between two momenta.
cepgen::Momentum::energyT2
double energyT2() const
Squared transverse energy component (in GeV )
cepgen::Momentum::rapidity
double rapidity() const
Rapidity.
cepgen::Momentum::deltaPt
double deltaPt(const Momentum &) const
Transverse momentum distance between two momenta.
cepgen::Momentum::gamma2
double gamma2() const
Squared gamma scalar value.
cepgen::Momentum::deltaPhi
double deltaPhi(const Momentum &) const
Azimuthal angle opening between two momenta.
cepgen::Momentum::eta
double eta() const
Pseudo-rapidity.
cepgen::Momentum::px
double px() const
Momentum along the -axis (in GeV)
Definition
Momentum.h:97
cepgen::Momentum::gamma
double gamma() const
Gamma scalar value.
cepgen::Momentum::p2
double p2() const
Squared 3-momentum norm (in GeV )
Definition
Momentum.h:105
cepgen::Momentum::py
double py() const
Momentum along the -axis (in GeV)
Definition
Momentum.h:98
cepgen::Momentum::beta
double beta() const
Beta scalar value.
cepgen::Momentum::massT2
double massT2() const
Squared transverse mass (in GeV )
cepgen::Momentum::energy2
double energy2() const
Squared energy (in GeV )
Definition
Momentum.h:108
cepgen::Momentum::massT
double massT() const
Transverse mass (in GeV)
cepgen::Momentum::pt
double pt() const
Transverse momentum (in GeV)
cepgen::Momentum::phi
double phi() const
Azimuthal angle (angle in the transverse plane)
cepgen::Momentum::p
double p() const
3-momentum norm (in GeV)
Definition
Momentum.h:104
cepgen::Momentum::mass
double mass() const
cepgen::Momentum::energyT
double energyT() const
Transverse energy component (in GeV)
cepgen::Momentum::energy
double energy() const
Energy (in GeV)
Definition
Momentum.h:107
cepgen::Momentum::theta
double theta() const
Polar angle (angle with respect to the longitudinal direction)
cepgen::Momentum::mass2
double mass2() const
Squared mass (in GeV ) as computed from its energy and momentum.
cepgen::Particle
Kinematic information for one particle.
Definition
Particle.h:32
cepgen::Particle::Role::OutgoingBeam2
@ OutgoingBeam2
outgoing beam state/particle
cepgen::Particle::Role::CentralSystem
@ CentralSystem
Central particles system.
cepgen::Particle::Role::Parton2
@ Parton2
beam incoming parton
cepgen::Particle::Role::OutgoingBeam1
@ OutgoingBeam1
outgoing beam state/particle
cepgen::Particle::Role::IncomingBeam2
@ IncomingBeam2
incoming beam particle
cepgen::Particle::Role::Intermediate
@ Intermediate
Intermediate two-parton system.
cepgen::Particle::Role::Parton1
@ Parton1
beam incoming parton
cepgen::Particle::Role::IncomingBeam1
@ IncomingBeam1
incoming beam particle
cepgen::utils::EventBrowser
A user-friendly browser for the Event content.
Definition
EventBrowser.h:33
cepgen::utils::EventBrowser::EventBrowser
EventBrowser()=default
cepgen::utils::EventBrowser::get
double get(const Event &event, const std::string &variable_name) const
Get/compute a variable value.
cepgen::utils
Collection of utilities.
Definition
RunParameters.h:33
cepgen
Common namespace for this Monte Carlo generator.
Definition
Handler.h:26
include
CepGen
EventFilter
EventBrowser.h
Generated on Tue Apr 22 2025 for CepGen by
1.10.0