cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
CollinearFlux.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
19#include <Pythia8/PartonDistributions.h>
20
21#include <numeric>
22
26#include "CepGen/Physics/PDG.h"
27
28namespace cepgen {
29 namespace pythia8 {
30 class CollinearFlux final : public cepgen::CollinearFlux {
31 public:
32 explicit CollinearFlux(const ParametersList& params)
33 : cepgen::CollinearFlux(params), type_(steer<std::string>("type")), pdgid_(steer<pdgid_t>("partonPdgId")) {
34 if (type_ == "Lepton") {
35 auto lepton_params = steer<ParametersList>("leptonParams");
36 info_.reset(new Pythia8::Info);
37 if (const auto dil_sqrt_s = lepton_params.get<double>("sqrtS"); dil_sqrt_s > 0.)
38 info_->setECM(dil_sqrt_s);
39 else
40 CG_WARNING("pythia8:CollinearFlux") << "Beam-beam centre-of-mass energy is required (through the 'sqrtS' "
41 "parameter) for the 'Lepton' collinear flux mode.";
42 pdf_.reset(new Pythia8::Lepton(
43 lepton_params.get<pdgid_t>("beamPdgId"), lepton_params.get<double>("Q2max"), info_.get()));
44 } else if (type_ == "LHAGrid1")
45 pdf_.reset(new Pythia8::LHAGrid1);
46 else if (type_ == "MSTWpdf")
47 pdf_.reset(new Pythia8::MSTWpdf);
48 else if (type_ == "Proton2gammaDZ")
49 pdf_.reset(new Pythia8::Proton2gammaDZ);
50 else if (type_ == "ProtonPoint")
51 pdf_.reset(new Pythia8::ProtonPoint);
52 else
53 throw CG_FATAL("pythia8:CollinearFlux") << "Failed to initialise the Pythia 8 evaluator!\n"
54 << "Parameters: " << params_;
55
56 CG_INFO("pythia8:CollinearFlux") << "Pythia 8 '" << type_ << "' evaluator for collinear parton "
57 << "(" << (PDG::Id)pdgid_ << ") flux initialised.";
58 }
59
62 desc.setDescription("Pythia 8 coll.flux");
63 desc.add<std::string>("type", "Proton2gammaDZ").setDescription("type of PDF evaluator to use");
64 desc.add<pdgid_t>("partonPdgId", PDG::photon).setDescription("parton PDG identifier");
65 auto lepton_desc = ParametersDescription();
66 lepton_desc.add<pdgid_t>("beamPdgId", PDG::electron).setDescription("beam particle PDG identifier");
67 lepton_desc.add<double>("sqrtS", -1.);
68 lepton_desc.add<double>("Q2max", 50.);
69 desc.add<ParametersDescription>("leptonParams", lepton_desc);
70 return desc;
71 }
72
73 pdgid_t partonPdgId() const override { return pdgid_; }
74 bool fragmenting() const override { return true; }
75 double mass2() const override { return mp2_; }
76
77 double fluxQ2(double x, double q2) const override {
78 if (x == 0. || x < pdf_->getXmin())
79 return 0.;
80 return pdf_->xf((int)pdgid_, x, q2);
81 }
82
83 private:
84 std::unique_ptr<Pythia8::PDF> pdf_;
85 std::unique_ptr<Pythia8::Info> info_;
86 const std::string type_;
87 const pdgid_t pdgid_;
88 };
89 } // namespace pythia8
90} // namespace cepgen
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_WARNING(mod)
Definition Message.h:228
#define CG_INFO(mod)
Definition Message.h:216
#define REGISTER_COLLINEAR_FLUX(name, obj)
Add a generic collinear parton flux evaluator builder definition.
static ParametersDescription description()
A class-in-the-middle PDG identifier for printout operations.
Definition PDG.h:55
@ electron
Definition PDG.h:37
@ photon
Definition PDG.h:41
A description object for parameters collection.
const double mp2_
Definition PartonFlux.h:39
ParametersList params_
Module parameters.
Definition Steerable.h:50
T steer(const std::string &key) const
Retrieve a parameters as previously steered.
Definition Steerable.h:39
bool fragmenting() const override
Is initiator particle fragmenting after parton emission?
double fluxQ2(double x, double q2) const override
Compute the collinear flux for this x value and virtuality.
static ParametersDescription description()
CollinearFlux(const ParametersList &params)
pdgid_t partonPdgId() const override
Parton PDG identifier.
double mass2() const override
Initiator particle squared mass (in )
Common namespace for this Monte Carlo generator.
unsigned long long pdgid_t
Alias for the integer-like particle PDG id.