cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
FortranInterface.cpp
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
20#include "CepGen/Generator.h"
28#include "CepGen/Physics/PDG.h"
30
31static std::unordered_map<int, std::unique_ptr<cepgen::strfun::Parameterisation> > kBuiltStrFunsParameterisations;
32static std::unordered_map<int, std::unique_ptr<cepgen::KTFlux> > kBuiltKtFluxParameterisations;
33static std::unordered_map<std::string, std::unique_ptr<cepgen::Coupling> > kBuiltAlphaSParameterisations,
35
36#ifdef __cplusplus
37extern "C" {
38#endif
40void cepgen_structure_functions_(int& sfmode, double& xbj, double& q2, double& f2, double& fl) {
41 using namespace cepgen;
42 if (kBuiltStrFunsParameterisations.count(sfmode) == 0)
43 kBuiltStrFunsParameterisations[sfmode] = StructureFunctionsFactory::get().build(sfmode);
44 const auto& sf = kBuiltStrFunsParameterisations.at(sfmode);
45 f2 = sf->F2(xbj, q2);
46 fl = sf->FL(xbj, q2);
47}
48
56double cepgen_kt_flux_(int& fmode, double& x, double& kt2, int& sfmode, double& min, double& mout) {
57 using namespace cepgen;
58 if (kBuiltKtFluxParameterisations.count(fmode) == 0)
59 kBuiltKtFluxParameterisations[fmode] = KTFluxFactory::get().build(
60 fmode,
62 .set<double>("mass", min)
63 .set<ParametersList>("structureFunctions",
64 StructureFunctionsFactory::get().describeParameters(sfmode).parameters())
65 .set<ParametersList>(
66 "formFactors",
67 FormFactorsFactory::get()
68 .describeParameters(formfac::gFFStandardDipoleHandler) // use another argument for the modelling?
69 .parameters()));
70 return kBuiltKtFluxParameterisations.at(fmode)->fluxMX2(x, kt2, mout * mout);
71}
72
79double cepgen_kt_flux_hi_(int& fmode, double& x, double& kt2, int& a, int& z) {
80 using namespace cepgen;
81 if (kBuiltKtFluxParameterisations.count(fmode) == 0)
82 kBuiltKtFluxParameterisations[fmode] = KTFluxFactory::get().build(
83 fmode, ParametersList().setAs<pdgid_t, HeavyIon>("heavyIon", HeavyIon{(unsigned short)a, (Element)z}));
84 return kBuiltKtFluxParameterisations.at(fmode)->fluxMX2(x, kt2, 0.);
85}
86
88double cepgen_particle_mass_(int& pdg_id) {
89 try {
90 return cepgen::PDG::get().mass((cepgen::pdgid_t)pdg_id);
91 } catch (const cepgen::Exception& e) {
92 e.dump();
93 exit(0);
94 }
95}
96
98double cepgen_particle_charge_(int& pdg_id) {
99 try {
100 return cepgen::PDG::get().charge((cepgen::pdgid_t)pdg_id);
101 } catch (const cepgen::Exception& e) {
102 e.dump();
103 exit(0);
104 }
105}
106
108double cepgen_particle_colour_(int& pdg_id) {
109 try {
110 return cepgen::PDG::get().colours((cepgen::pdgid_t)pdg_id);
111 } catch (const cepgen::Exception& e) {
112 e.dump();
113 exit(0);
114 }
115}
116
118
119void cepgen_debug_(char* str, int size) { CG_DEBUG("fortran_process") << std::string(str, size); }
120
121void cepgen_warning_(char* str, int size) { CG_WARNING("fortran_process") << std::string(str, size); }
122
123void cepgen_error_(char* str, int size) { CG_ERROR("fortran_process") << std::string(str, size); }
124
125void cepgen_fatal_(char* str, int size) { throw CG_FATAL("fortran_process") << std::string(str, size); }
126
127double cepgen_alphas_(char* str, double& q, int size) {
128 const auto name = std::string(str, size);
129 if (name.empty())
130 throw CG_FATAL("cepgen_alphas") << "Invalid name retrieved for alphaS(q) parameterisation: '" << name << "'.";
131 if (kBuiltAlphaSParameterisations.count(name) == 0)
132 kBuiltAlphaSParameterisations[name] = cepgen::AlphaSFactory::get().build(name);
133 return kBuiltAlphaSParameterisations.at(name)->operator()(q);
134}
135
136double cepgen_alphaem_(char* str, double& q, int size) {
137 const auto name = std::string(str, size);
138 if (name.empty())
139 throw CG_FATAL("cepgen_alphaem") << "Invalid name retrieved for alphaEM(q) parameterisation: '" << name << "'.";
140 if (kBuiltAlphaEMParameterisations.count(name) == 0)
141 kBuiltAlphaEMParameterisations[name] = cepgen::AlphaEMFactory::get().build(name);
142 return kBuiltAlphaEMParameterisations.at(name)->operator()(q);
143}
144
145#ifdef __cplusplus
146}
147#endif
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_ERROR(mod)
Definition Exception.h:60
double cepgen_kt_flux_(int &fmode, double &x, double &kt2, int &sfmode, double &min, double &mout)
Compute a -dependent flux for single nucleons.
void cepgen_init_()
static std::unordered_map< int, std::unique_ptr< cepgen::strfun::Parameterisation > > kBuiltStrFunsParameterisations
double cepgen_particle_colour_(int &pdg_id)
Colour factor of a particle.
double cepgen_particle_mass_(int &pdg_id)
Mass of a particle, in GeV/c^2.
static std::unordered_map< int, std::unique_ptr< cepgen::KTFlux > > kBuiltKtFluxParameterisations
void cepgen_debug_(char *str, int size)
double cepgen_particle_charge_(int &pdg_id)
Charge of a particle, in e.
void cepgen_fatal_(char *str, int size)
static std::unordered_map< std::string, std::unique_ptr< cepgen::Coupling > > kBuiltAlphaEMParameterisations
static std::unordered_map< std::string, std::unique_ptr< cepgen::Coupling > > kBuiltAlphaSParameterisations
void cepgen_error_(char *str, int size)
double cepgen_alphaem_(char *str, double &q, int size)
double cepgen_kt_flux_hi_(int &fmode, double &x, double &kt2, int &a, int &z)
Compute a -dependent flux for heavy ions.
void cepgen_warning_(char *str, int size)
void cepgen_structure_functions_(int &sfmode, double &xbj, double &q2, double &f2, double &fl)
Expose structure functions calculators to Fortran.
double cepgen_alphas_(char *str, double &q, int size)
#define CG_WARNING(mod)
Definition Message.h:228
#define CG_DEBUG(mod)
Definition Message.h:220
void dump(std::ostream *os=nullptr) const noexcept override
Human-readable dump of the exception.
Definition Exception.cpp:41
Core generator object allowing for process definition, cross section computation, and event generatio...
Definition Generator.h:48
double colours(spdgid_t) const
Colour factor for this particle.
Definition PDG.cpp:88
double mass(spdgid_t) const
Particle mass (in GeV)
Definition PDG.cpp:90
static PDG & get()
Retrieve a unique instance of this particles info collection.
Definition PDG.cpp:41
double charge(spdgid_t) const
Electric charge (in ) for this particle.
Definition PDG.cpp:98
Common namespace for this Monte Carlo generator.
Element
Enumeration of chemical elements.
Definition HeavyIon.h:28
unsigned long long pdgid_t
Alias for the integer-like particle PDG id.
Heavy ion container (Z+A)
Definition HeavyIon.h:44