cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
PythonWrapper.cpp
Go to the documentation of this file.
1
/*
2
* CepGen: a central exclusive processes event generator
3
* Copyright (C) 2023 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
#include "
CepGen/FormFactors/Parameterisation.h
"
20
#include "
CepGen/Generator.h
"
21
#include "
CepGen/Modules/FormFactorsFactory.h
"
22
#include "
CepGen/Modules/PartonFluxFactory.h
"
23
#include "
CepGen/Modules/StructureFunctionsFactory.h
"
24
#include "
CepGen/Physics/PDG.h
"
25
#include "
CepGen/StructureFunctions/Parameterisation.h
"
26
#include "
CepGen/StructureFunctions/SigmaRatio.h
"
27
#include "
CepGenAddOns/BoostWrapper/PythonObjectsWrappers.h
"
28
#include "
CepGenAddOns/BoostWrapper/PythonUtils.h
"
29
30
namespace
{
32
BOOST_PYTHON_MODULE(pycepgen) {
33
cepgen::initialise
();
34
35
py::class_<cepgen::Steerable>(
"_Steerable"
,
"base steerable object"
, py::no_init)
36
.add_property(
37
"parameters"
,
38
+[](
const
cepgen::Steerable
& st) {
return
plist_to_py_dict
(st.
parameters
()); },
39
&
cepgen::strfun::Parameterisation::setParameters
,
40
"Operational parameters"
)
41
.add_property(
42
"name"
,
43
+[](
const
cepgen::Steerable
& st) {
return
st.
parameters
().
getString
(
cepgen::MODULE_NAME
); },
44
"Module name"
);
45
46
py::class_<cepgen::strfun::Parameterisation, py::bases<cepgen::Steerable>, boost::noncopyable>(
47
"_StructureFunctions"
,
"nucleon structure functions modelling"
, py::no_init)
48
.add_static_property(
"name"
,
49
py::make_function(&
cepgen::strfun::Parameterisation::name
,
50
py::return_value_policy<py::copy_const_reference>()))
51
.def(
"F2"
, &
cepgen::strfun::Parameterisation::F2
)
52
.def(
"FL"
, &
cepgen::strfun::Parameterisation::FL
)
53
.def(
"F1"
, &
cepgen::strfun::Parameterisation::F1
);
54
55
EXPOSE_FACTORY
(cepgen::StructureFunctionsFactory,
56
std::string,
57
"StructureFunctionsFactory"
,
58
"a structure functions evaluator objects factory"
)
59
.def(
"build"
,
adapt_unique
(+[](
int
mod) {
return
cepgen::StructureFunctionsFactory::get().build(mod); }));
60
61
py::class_<cepgen::sigrat::Parameterisation, py::bases<cepgen::Steerable>, boost::noncopyable>(
62
"_SigmaRatio"
,
"L/T cross section ratio modelling"
, py::no_init)
63
.def(
64
"__call__"
, +[](
const
cepgen::sigrat::Parameterisation
& par,
double
xbj,
double
q2) {
65
double
unc{0.};
66
const
auto
sig_rat = par(xbj, q2, unc);
67
return
std_vector_to_py_tuple
(std::vector<double>{sig_rat, unc});
68
});
69
70
EXPOSE_FACTORY
(cepgen::SigmaRatiosFactory,
71
int
,
72
"SigmaRatiosFactory"
,
73
"a longitudinal-to-transverse cross section ratio evaluator objects factory"
);
74
75
py::class_<cepgen::formfac::Parameterisation, py::bases<cepgen::Steerable>, boost::noncopyable>(
76
"_FormFactors"
,
"nucleon electromagnetic form factors modelling"
, py::no_init)
77
.def(
"__call__"
,
78
py::make_function(&cepgen::formfac::Parameterisation::operator(), py::return_internal_reference()));
79
80
py::class_<cepgen::formfac::FormFactors>(
"FormFactors"
,
"nucleon electromagnetic form factors values"
)
81
.add_property(
"FE"
, &
cepgen::formfac::FormFactors::FE
,
"Electric form factor"
)
82
.add_property(
"FM"
, &
cepgen::formfac::FormFactors::FM
,
"Magnetic form factor"
)
83
.add_property(
"GE"
, &
cepgen::formfac::FormFactors::GE
,
"Sachs electric form factor"
)
84
.add_property(
"GM"
, &
cepgen::formfac::FormFactors::GM
,
"Sachs magnetic form factor"
);
85
86
EXPOSE_FACTORY
(cepgen::FormFactorsFactory,
87
std::string,
88
"FormFactorsFactory"
,
89
"an electromagnetic form factors evaluator objects factory"
);
90
91
py::class_<PartonFluxWrap, py::bases<cepgen::Steerable>, boost::noncopyable>(
92
"_PartonFlux"
,
"generic parton flux evaluator"
, py::no_init)
93
.add_property(
"partonPdgId"
, &
cepgen::PartonFlux::partonPdgId
)
94
.add_property(
"fragmenting"
, &
cepgen::PartonFlux::fragmenting
)
95
.add_property(
"ktFactorised"
, &
cepgen::PartonFlux::ktFactorised
)
96
.def(
97
"__call__"
,
98
+[](
const
cepgen::PartonFlux
& flux) {
99
if
(flux.
ktFactorised
())
100
return
adapt_reference
(&
dynamic_cast<
const
cepgen::KTFlux
&
>
(flux));
101
return
adapt_reference
(&
dynamic_cast<
const
cepgen::CollinearFlux
&
>
(flux));
102
},
103
"Expose the flux evaluator object from its type"
);
104
105
py::class_<CollinearFluxWrap, py::bases<PartonFluxWrap>, boost::noncopyable>(
106
"_CollinearFlux"
,
"fractional momentum/virtuality-dependent parton flux evaluator"
, py::no_init)
107
.def(
"fluxMX2"
, py::pure_virtual(&
cepgen::CollinearFlux::fluxMX2
))
108
.def(
"fluxQ2"
, py::pure_virtual(&
cepgen::CollinearFlux::fluxQ2
));
109
110
EXPOSE_FACTORY
(cepgen::CollinearFluxFactory,
111
std::string,
112
"CollinearFluxFactory"
,
113
"a collinear parton fluxes evaluator objects factory"
);
114
115
py::class_<KTFluxWrap, py::bases<PartonFluxWrap>, boost::noncopyable>(
116
"_KTFlux"
,
117
"fractional momentum/(transverse-longitudinal) virtuality-dependent parton flux evaluator"
,
118
py::no_init)
119
.def(
"fluxMX2"
, py::pure_virtual(&
cepgen::KTFlux::fluxMX2
))
120
.def(
"fluxQ2"
, py::pure_virtual(&
cepgen::KTFlux::fluxQ2
));
121
122
EXPOSE_FACTORY
(
123
cepgen::KTFluxFactory, std::string,
"KTFluxFactory"
,
"a kt-factorised parton fluxes evaluator objects factory"
);
124
125
py::class_<cepgen::PDG, boost::noncopyable>(
"PDG"
,
"collection of particle definitions and properties"
, py::no_init)
126
.def(
127
"colours"
, +[](
cepgen::pdgid_t
pdgid) {
return
cepgen::PDG::get
().
colours
(pdgid); })
128
.def(
129
"mass"
, +[](
cepgen::pdgid_t
pdgid) {
return
cepgen::PDG::get
().
mass
(pdgid); })
130
.def(
131
"width"
, +[](
cepgen::pdgid_t
pdgid) {
return
cepgen::PDG::get
().
width
(pdgid); })
132
.def(
133
"charge"
, +[](
cepgen::pdgid_t
pdgid) {
return
cepgen::PDG::get
().
charge
(pdgid); });
134
}
135
}
// namespace
FormFactorsFactory.h
Parameterisation.h
Generator.h
PDG.h
PartonFluxFactory.h
PythonObjectsWrappers.h
plist_to_py_dict
py::dict plist_to_py_dict(const cepgen::ParametersList &plist)
Definition
PythonUtils.cpp:71
PythonUtils.h
EXPOSE_FACTORY
#define EXPOSE_FACTORY(obj, key, name, description)
Definition
PythonUtils.h:25
adapt_reference
py::object adapt_reference(T *ptr)
Definition
PythonUtils.h:94
std_vector_to_py_tuple
py::tuple std_vector_to_py_tuple(const std::vector< T > &vec)
Definition
PythonUtils.h:72
adapt_unique
py::object adapt_unique(std::unique_ptr< T >(*fn)(Args...))
Definition
PythonUtils.h:80
SigmaRatio.h
StructureFunctionsFactory.h
Parameterisation.h
cepgen::CollinearFlux
Definition
CollinearFlux.h:25
cepgen::CollinearFlux::fluxQ2
virtual double fluxQ2(double x, double q2) const
Compute the collinear flux for this x value and virtuality.
Definition
CollinearFlux.cpp:31
cepgen::CollinearFlux::fluxMX2
virtual double fluxMX2(double x, double mf2=0.) const
Compute the collinear flux for this x value and remnant mass.
Definition
CollinearFlux.cpp:35
cepgen::KTFlux
Definition
KTFlux.h:25
cepgen::KTFlux::fluxMX2
virtual double fluxMX2(double x, double kt2, double mf2) const
Compute the kt-dependent flux for this x value and remnant mass.
Definition
KTFlux.cpp:36
cepgen::KTFlux::fluxQ2
virtual double fluxQ2(double x, double kt2, double q2) const
Compute the kt-dependent flux for this x value and virtuality.
Definition
KTFlux.cpp:32
cepgen::NamedModule< Parameterisation >::name
const std::string & name() const
Module unique indexing name.
Definition
NamedModule.h:42
cepgen::PDG::colours
double colours(spdgid_t) const
Colour factor for this particle.
Definition
PDG.cpp:88
cepgen::PDG::width
double width(spdgid_t) const
Resonance width (in GeV)
Definition
PDG.cpp:96
cepgen::PDG::mass
double mass(spdgid_t) const
Particle mass (in GeV)
Definition
PDG.cpp:90
cepgen::PDG::get
static PDG & get()
Retrieve a unique instance of this particles info collection.
Definition
PDG.cpp:41
cepgen::PDG::charge
double charge(spdgid_t) const
Electric charge (in ) for this particle.
Definition
PDG.cpp:98
cepgen::ParametersList::getString
std::string getString(const std::string &key, bool wrap=false) const
Get a string-converted version of a value.
Definition
ParametersList.cpp:313
cepgen::PartonFlux
Definition
PartonFlux.h:26
cepgen::PartonFlux::partonPdgId
virtual pdgid_t partonPdgId() const =0
Parton PDG identifier.
cepgen::PartonFlux::fragmenting
virtual bool fragmenting() const =0
Is initiator particle fragmenting after parton emission?
cepgen::PartonFlux::ktFactorised
virtual bool ktFactorised() const
Is the flux parton kT-dependent?
Definition
PartonFlux.h:32
cepgen::Steerable
Base runtime module object.
Definition
Steerable.h:26
cepgen::Steerable::parameters
virtual const ParametersList & parameters() const
Module parameters.
Definition
Steerable.h:33
cepgen::Steerable::setParameters
virtual void setParameters(const ParametersList &)
Set module parameters.
Definition
Steerable.cpp:28
cepgen::sigrat::Parameterisation
A generic modelling of the ratio.
Definition
SigmaRatio.h:28
cepgen::strfun::Parameterisation::F2
double F2(double xbj, double q2)
Transverse structure function.
Definition
Parameterisation.cpp:57
cepgen::strfun::Parameterisation::F1
double F1(double xbj, double q2)
structure function
Definition
Parameterisation.cpp:73
cepgen::strfun::Parameterisation::FL
double FL(double xbj, double q2)
Longitudinal structure function.
Definition
Parameterisation.cpp:59
cepgen::formfac::FormFactors::GE
double GE
Sachs electric form factor.
Definition
FormFactors.h:29
cepgen::formfac::FormFactors::GM
double GM
Sachs magnetic form factor.
Definition
FormFactors.h:30
cepgen::formfac::FormFactors::FE
double FE
Electric form factor.
Definition
FormFactors.h:27
cepgen::formfac::FormFactors::FM
double FM
Magnetic form factor.
Definition
FormFactors.h:28
cepgen::initialise
void initialise(bool safe_mode)
Definition
GlobalFunctions.cpp:91
cepgen::pdgid_t
unsigned long long pdgid_t
Alias for the integer-like particle PDG id.
Definition
ParticleProperties.h:26
cepgen::MODULE_NAME
const char *const MODULE_NAME
Indexing key for the module name Parameters container.
Definition
ParametersList.h:50
CepGenAddOns
BoostWrapper
PythonWrapper.cpp
Generated on Mon Jul 29 2024 for CepGen by
1.9.7