cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
InelasticKTFluxes.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2023-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
23#include "CepGen/Physics/PDG.h"
26#include "CepGen/Utils/Math.h"
27
28namespace cepgen {
30 public:
31 explicit InelasticNucleonKTFlux(const ParametersList& params)
32 : KTFlux(params), sf_(StructureFunctionsFactory::get().build(steer<ParametersList>("structureFunctions"))) {
33 if (!sf_)
34 throw CG_FATAL("InelasticNucleonKTFlux") << "Inelastic kT flux requires a modelling of structure functions!";
35 CG_DEBUG("InelasticNucleonKTFlux") << "Inelastic KT-dependent flux initialised with '"
36 << steer<ParametersList>("structureFunctions")
37 << "' structure functions modelling.";
38 }
39
41 auto desc = KTFlux::description();
42 desc.setDescription("Nucl. inel. photon emission");
43 desc.add<ParametersDescription>("structureFunctions",
44 StructureFunctionsFactory::get().describeParameters("LUXLike"));
45 return desc;
46 }
47
48 double mass2() const override { return mp2_; }
49 bool fragmenting() const override final { return true; }
50 pdgid_t partonPdgId() const override { return PDG::photon; }
51 double fluxMX2(double x, double kt2, double mx2) const override {
52 if (!x_range_.contains(x, true))
53 return 0.;
54 if (!utils::positive(mx2)) {
55 CG_WARNING("InelasticNucleonKTFlux") << "Invalid diffractive mass squared mX^2 specified: " << mx2 << ".";
56 return 0.;
57 }
58 const auto q2 = utils::kt::q2(x, kt2, mass2(), mx2), q2min = q2 - kt2 / (1. - x);
59 const auto xbj = utils::xBj(q2, mass2(), mx2), qnorm = 1. - q2min / q2;
60 return prefactor_ * sf_->F2(xbj, q2) * (xbj / q2) * qnorm * qnorm * (1. - x) / q2;
61 }
62
63 protected:
64 std::unique_ptr<strfun::Parameterisation> sf_;
65 };
66
71 desc.setDescription("Nucl. inel. photon emission (Budnev flux)");
72 return desc;
73 }
74 double fluxMX2(double x, double kt2, double mx2) const override {
75 if (!x_range_.contains(x, true))
76 return 0.;
77 if (!utils::positive(mx2)) {
78 CG_WARNING("InelasticNucleonKTFlux") << "Invalid diffractive mass squared mX^2 specified: " << mx2 << ".";
79 return 0.;
80 }
81 const auto q2 = utils::kt::q2(x, kt2, mass2(), mx2), q2min = q2 - kt2 / (1. - x);
82 const auto xbj = utils::xBj(q2, mass2(), mx2), qnorm = 1. - q2min / q2;
83 const double f_D = sf_->F2(xbj, q2) * (xbj / q2) * (1. - x) * qnorm;
84 const double f_C = sf_->F1(xbj, q2) * 2. / q2;
85 return prefactor_ * (f_D + 0.5 * x * x * f_C) * (1. - x) / q2;
86 }
87 };
88} // namespace cepgen
89
90REGISTER_KT_FLUX("Inelastic", 1, InelasticNucleonKTFlux);
91REGISTER_KT_FLUX("BudnevInelastic", 11, BudnevInelasticNucleonKTFlux);
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_WARNING(mod)
Definition Message.h:228
#define CG_DEBUG(mod)
Definition Message.h:220
#define REGISTER_KT_FLUX(name, id, obj)
Add a generic KT-factorised flux evaluator builder definition.
bool fragmenting() const override final
Is initiator particle fragmenting after parton emission?
static ParametersDescription description()
InelasticNucleonKTFlux(const ParametersList &params)
std::unique_ptr< strfun::Parameterisation > sf_
double fluxMX2(double x, double kt2, double mx2) const override
Compute the kt-dependent flux for this x value and remnant mass.
pdgid_t partonPdgId() const override
Parton PDG identifier.
double mass2() const override
Initiator particle squared mass (in )
static ParametersDescription description()
Definition KTFlux.cpp:26
bool contains(double val, bool exclude_boundaries=false) const
Check if value is inside limits' boundaries.
Definition Limits.cpp:77
@ photon
Definition PDG.h:41
A description object for parameters collection.
const double prefactor_
Definition PartonFlux.h:38
const double mp2_
Definition PartonFlux.h:39
const Limits x_range_
Definition PartonFlux.h:40
T steer(const std::string &key) const
Retrieve a parameters as previously steered.
Definition Steerable.h:39
double q2(double x, double kt2, double mi2, double mx2)
Compute the virtuality from longitudinal loss/transverse virtuality/diffractive mass.
Definition Utils.cpp:60
double xBj(double q2, double mp2, double mx2)
Compute Bjorken x from virtuality/diffractive mass.
Definition Utils.cpp:41
bool positive(const T &val)
Check if a number is positive and finite.
Definition Math.cpp:26
Common namespace for this Monte Carlo generator.
unsigned long long pdgid_t
Alias for the integer-like particle PDG id.
static ParametersDescription description()
double fluxMX2(double x, double kt2, double mx2) const override
Compute the kt-dependent flux for this x value and remnant mass.