cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
PythiaEventInterface.h
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2016-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#ifndef CepGenAddOns_EventInterfaces_PythiaEventInterface_h
20#define CepGenAddOns_EventInterfaces_PythiaEventInterface_h
21
22#include <Pythia8/Pythia.h>
23
24#include <unordered_map>
25
26namespace cepgen {
27 class RunParameters;
28 class Event;
29 class Particle;
30} // namespace cepgen
31
32namespace Pythia8 {
34 class CepGenEvent : public LHAup {
35 public:
37 enum struct Type {
41 };
42 explicit CepGenEvent();
48 void feedEvent(const cepgen::Event& ev, const Type& type);
53 void setCrossSection(int id, double cross_section, double cross_section_err);
60 void setProcess(int id, double cross_section, double q2_scale, double alpha_qed, double alpha_qcd);
61
63 void addComments(const std::string& comments);
64
68 unsigned short cepgenId(unsigned short py_id) const;
72 unsigned short pythiaId(unsigned short cg_id) const;
74 void addCepGenParticle(const cepgen::Particle& part,
75 int status = INVALID_ID,
76 const std::pair<int, int>& mothers = {0, 0},
77 const std::pair<int, int>& colours = {0, 0});
81 void addCorresp(unsigned short py_id, unsigned short cg_id);
83 void dumpCorresp() const;
84
85 static constexpr unsigned short INVALID_ID = 999;
86 static constexpr unsigned short MIN_COLOUR_INDEX = 501;
87
88 inline bool setInit() override { return true; }
89#if defined(PYTHIA_VERSION_INTEGER) && PYTHIA_VERSION_INTEGER >= 8200
90 bool setEvent(int) override { return true; }
91#else
92 bool setEvent(int, double) override { return true; }
93#endif
94
95 private:
96 std::pair<int, int> findMothers(const cepgen::Event& ev, const cepgen::Particle& p) const;
97 const double mp_, mp2_;
98 bool inel1_{false}, inel2_{false};
99 std::unordered_map<unsigned short, unsigned short> py_cg_corresp_;
100 const cepgen::RunParameters* params_{nullptr}; // borrowed
101 };
102} // namespace Pythia8
103#endif
Interfacing between CepGen and Pythia8 event definitions.
void addCorresp(unsigned short py_id, unsigned short cg_id)
Register a new Pythia8 / CepGen particle mapping.
bool setEvent(int, double) override
void feedEvent(const cepgen::Event &ev, const Type &type)
Feed a new CepGen event to this conversion object.
Type
List of particles to be included to the event content.
@ centralAndPartons
only include initiators and central system
@ centralAndFullBeamRemnants
include dissociated beam remnants and central system
@ centralAndBeamRemnants
include undissociated beam remnants and central system
void setCrossSection(int id, double cross_section, double cross_section_err)
Set the cross section for a given process.
static constexpr unsigned short INVALID_ID
Invalid id association.
unsigned short cepgenId(unsigned short py_id) const
Retrieve the CepGen particle index given its Pythia8 event id.
void initialise(const cepgen::RunParameters &)
Initialise this conversion object with CepGen parameters.
void addComments(const std::string &comments)
Feed comments to the LHEF block.
void addCepGenParticle(const cepgen::Particle &part, int status=INVALID_ID, const std::pair< int, int > &mothers={0, 0}, const std::pair< int, int > &colours={0, 0})
Add a CepGen particle to the event content.
static constexpr unsigned short MIN_COLOUR_INDEX
Minimal colour indexing number.
unsigned short pythiaId(unsigned short cg_id) const
Retrieve the Pythia8 particle index given its CepGen event id.
void dumpCorresp() const
Print all Pythia8 / CepGen Particles correspondences.
void setProcess(int id, double cross_section, double q2_scale, double alpha_qed, double alpha_qcd)
Specify new process attributes.
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Kinematic information for one particle.
Definition Particle.h:33
List of parameters used to start and run the simulation job.
Common namespace for this Monte Carlo generator.