ATLAS Offline Software
Loading...
Searching...
No Matches
Pythia8_i.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 Author: James Monk
4*/
5
6#ifndef GENERATOR_PYTHIA8_H
7#define GENERATOR_PYTHIA8_H
8
10
11// calls to fortran routines
12#include "CLHEP/Random/RandFlat.h"
13#include "GaudiKernel/ToolHandle.h"
16
17//#include "Pythia8/../Pythia8Plugins/HepMC2.h"
18#ifdef HEPMC3
19#include "Pythia8Plugins/HepMC3.h"
20namespace HepMC {
21 typedef HepMC3::Pythia8ToHepMC3 Pythia8ToHepMC;
22}
23#else
24#include "Pythia8Plugins/HepMC2.h"
25#endif
26
27#include <stdexcept>
28
29
33
34namespace Pythia8{
35 class Sigma2Process;
36}
37
38
39class customRndm : public Pythia8::RndmEngine {
40public:
41
42 // Constructor.
44
45 // Routine for generating a random number.
46 inline double flat(){
47 m_RNC++;
48 return CLHEP::RandFlat::shoot(m_engine);
49 };
50
51 // Initialisation Routine
52 inline void init(CLHEP::HepRandomEngine* engine) {m_engine=engine; m_RNC=0;}
53 inline CLHEP::HepRandomEngine* getEngine() { return m_engine; }
54 inline unsigned long long int getRNCalls() {return m_RNC;}
55private:
56 unsigned long long int m_RNC{};
57 CLHEP::HepRandomEngine* m_engine{};
58};
59
60
61//namespace Generator{
62class Pythia8_i: public GenModule{
63
64public:
65 Pythia8_i(const std::string &name, ISvcLocator *pSvcLocator);
66
67 ~Pythia8_i();
68
69 class CommandException : public std::runtime_error{
70 public:
71
72 CommandException(const std::string &cmd): std::runtime_error("Cannot interpret command: " + cmd){
73 }
74 };
75
76 virtual StatusCode genInitialize();
77 virtual StatusCode callGenerator();
78 virtual StatusCode fillEvt(HepMC::GenEvent *evt);
79 virtual StatusCode fillWeights(HepMC::GenEvent *evt);
80 virtual StatusCode genFinalize();
81
82 double pythiaVersion() const;
83
84 static const std::string& pythia_stream();
85 static std::string xmlpath();
86
87protected:
88
89 bool useRndmGenSvc() const { return m_useRndmGenSvc; }
90 bool useReseed() const {return m_useReseed; }
91
92 std::unique_ptr<Pythia8::Pythia> m_pythia{};
93 HepMC::Pythia8ToHepMC m_pythiaToHepMC;
94 UnsignedIntegerProperty m_maxFailures{this, "MaxFailures", 10};
95
96 BooleanProperty m_useRndmGenSvc{this, "useRndmGenSvc", true, "the max number of consecutive failures"};
97 std::shared_ptr<customRndm> m_atlasRndmEngine{};
98
99 BooleanProperty m_useReseed{this,"useReseed", false};
100
101 IntegerProperty m_dsid{this, "Dsid", 999999, "Dataset ID number"};
102 StringArrayProperty m_userHooks{this, "UserHooks", {} };
103//for Py8B
104 DoubleProperty m_pt0timesMPI{this,"pT0timesMPI", 1.0};
105 DoubleProperty m_numberAlphaS{this,"numberAlphaS", 3.0};
106 BooleanProperty m_sameAlphaSAsMPI{this,"useSameAlphaSasMPI", false};
107
108private:
109
110 static std::string findValue(const std::string &command, const std::string &key);
111
113
114 double m_version{-1.};
115
116 StringArrayProperty m_commands{this, "Commands", {} };
117 std::vector<std::string> m_userParams;
118 std::vector<std::string> m_userModes;
119
120 enum PDGID {PROTON=2212, ANTIPROTON=-2212, LEAD=1000822080, OXYGEN=1000080160, HELIUM= 1000020040, NEUTRON=2112, ANTINEUTRON=-2112, MUON=13, ANTIMUON=-13, ELECTRON=11, POSITRON=-11, INVALID=0};
121
122 DoubleProperty m_collisionEnergy{this, "CollisionEnergy", 14000.0};
123
124
125 StringProperty m_beam1{this, "Beam1", "PROTON"};
126 StringProperty m_beam2{this, "Beam2", "PROTON"};
128
129 StringProperty m_lheFile{this, "LHEFile", ""};
130
131 BooleanProperty m_doCKKWLAcceptance{this, "CKKWLAcceptance", false};
132 BooleanProperty m_doFxFxXS{this, "FxFxXS", false};
133 BooleanProperty m_computeEfficiency{this, "computeEfficiency", false};
134 double m_nAccepted{0.};
135 double m_nMerged{0.};
136 double m_sigmaTotal{0.};
137 double m_conversion{1.};
138
139 unsigned int m_failureCount{0};
140
141 std::map<std::string, PDGID> m_particleIDs;
142
143 std::vector<long int> m_seeds{};
144
145 StringProperty m_userProcess{this, "UserProcess", ""};
146
147 // ptr to possible user process
148 std::shared_ptr<Pythia8::Sigma2Process> m_procPtr{};
149
150 std::vector<UserHooksPtrType> m_userHooksPtrs{};
151
152 StringProperty m_userResonances{this, "UserResonances", ""};
153
154 std::vector<std::shared_ptr<Pythia8::ResonanceWidths> > m_userResonancePtrs;
155
156 BooleanProperty m_useLHAPDF{this, "UseLHAPDF", true};
157
158 StringProperty m_particleDataFile{this, "ParticleData", ""};
159 StringProperty m_outputParticleDataFile{this, "OutputParticleData", "ParticleData.local.xml"};
160
161 double m_mergingWeight{1.0};
162 double m_enhanceWeight{1.0};
163 std::vector<std::string> m_weightIDs{};
164 std::vector<std::string> m_weightNames{};
165 bool m_doLHE3Weights{false};
166 std::vector<std::string> m_weightCommands{};
167 std::vector<std::string> m_showerWeightNames{"nominal"};
168 StringArrayProperty m_showerWeightNamesProp{this, "ShowerWeightNames", {} };
169
170 PublicToolHandle<IPythia8Custom> m_athenaTool{this, "CustomInterface", ""};
171
172 BooleanProperty m_saveLHE{this, "SaveLHERecord", false};
173
174 static int s_allowedTunes(double version);
175
176 Pythia8::SuppressSmallPT *m_SuppressSmallPT{};
177
178};
179
180#endif
GenModule(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition GenModule.cxx:14
CommandException(const std::string &cmd)
Definition Pythia8_i.h:72
bool m_override_transform_beamenergy
Definition Pythia8_i.h:127
virtual StatusCode genInitialize()
For initializing the generator, if required.
Definition Pythia8_i.cxx:98
IntegerProperty m_dsid
Definition Pythia8_i.h:101
BooleanProperty m_sameAlphaSAsMPI
Definition Pythia8_i.h:106
std::vector< std::string > m_weightCommands
Definition Pythia8_i.h:166
virtual StatusCode fillWeights(HepMC::GenEvent *evt)
virtual StatusCode fillEvt(HepMC::GenEvent *evt)
For filling the HepMC event object.
StringProperty m_userProcess
Definition Pythia8_i.h:145
std::vector< std::string > m_userParams
Definition Pythia8_i.h:117
BooleanProperty m_useReseed
Definition Pythia8_i.h:99
static std::string xmlpath()
static int s_allowedTunes(double version)
DoubleProperty m_numberAlphaS
Definition Pythia8_i.h:105
bool useRndmGenSvc() const
Definition Pythia8_i.h:89
Pythia8::SuppressSmallPT * m_SuppressSmallPT
Definition Pythia8_i.h:176
double m_nAccepted
Definition Pythia8_i.h:134
unsigned int m_failureCount
Definition Pythia8_i.h:139
bool m_doLHE3Weights
Definition Pythia8_i.h:165
DoubleProperty m_collisionEnergy
Definition Pythia8_i.h:122
StringProperty m_lheFile
Definition Pythia8_i.h:129
StringProperty m_beam1
Definition Pythia8_i.h:125
BooleanProperty m_computeEfficiency
Definition Pythia8_i.h:133
std::vector< std::shared_ptr< Pythia8::ResonanceWidths > > m_userResonancePtrs
Definition Pythia8_i.h:154
std::shared_ptr< Pythia8::Sigma2Process > m_procPtr
Definition Pythia8_i.h:148
double m_sigmaTotal
Definition Pythia8_i.h:136
BooleanProperty m_saveLHE
Definition Pythia8_i.h:172
HepMC::Pythia8ToHepMC m_pythiaToHepMC
Definition Pythia8_i.h:93
double m_nMerged
Definition Pythia8_i.h:135
PublicToolHandle< IPythia8Custom > m_athenaTool
Definition Pythia8_i.h:170
StringProperty m_userResonances
Definition Pythia8_i.h:152
virtual StatusCode callGenerator()
For calling the generator on each iteration of the event loop.
BooleanProperty m_doFxFxXS
Definition Pythia8_i.h:132
std::vector< long int > m_seeds
Definition Pythia8_i.h:143
std::vector< std::string > m_showerWeightNames
Definition Pythia8_i.h:167
std::vector< std::string > m_weightIDs
Definition Pythia8_i.h:163
DoubleProperty m_pt0timesMPI
Definition Pythia8_i.h:104
std::vector< std::string > m_weightNames
Definition Pythia8_i.h:164
int m_internal_event_number
Definition Pythia8_i.h:112
static const std::string & pythia_stream()
std::shared_ptr< customRndm > m_atlasRndmEngine
Definition Pythia8_i.h:97
std::vector< std::string > m_userModes
Definition Pythia8_i.h:118
double pythiaVersion() const
StringArrayProperty m_userHooks
Definition Pythia8_i.h:102
StringProperty m_beam2
Definition Pythia8_i.h:126
BooleanProperty m_useLHAPDF
Definition Pythia8_i.h:156
static std::string findValue(const std::string &command, const std::string &key)
virtual StatusCode genFinalize()
For finalising the generator, if required.
std::unique_ptr< Pythia8::Pythia > m_pythia
Definition Pythia8_i.h:92
double m_enhanceWeight
Definition Pythia8_i.h:162
double m_version
Definition Pythia8_i.h:114
double m_conversion
Definition Pythia8_i.h:137
std::map< std::string, PDGID > m_particleIDs
Definition Pythia8_i.h:141
BooleanProperty m_doCKKWLAcceptance
Definition Pythia8_i.h:131
double m_mergingWeight
Definition Pythia8_i.h:161
UnsignedIntegerProperty m_maxFailures
Definition Pythia8_i.h:94
BooleanProperty m_useRndmGenSvc
Definition Pythia8_i.h:96
StringProperty m_particleDataFile
Definition Pythia8_i.h:158
std::vector< UserHooksPtrType > m_userHooksPtrs
Definition Pythia8_i.h:150
StringArrayProperty m_showerWeightNamesProp
Definition Pythia8_i.h:168
StringProperty m_outputParticleDataFile
Definition Pythia8_i.h:159
StringArrayProperty m_commands
Definition Pythia8_i.h:116
bool useReseed() const
Definition Pythia8_i.h:90
Pythia8_i(const std::string &name, ISvcLocator *pSvcLocator)
Definition Pythia8_i.cxx:68
CLHEP::HepRandomEngine * m_engine
Definition Pythia8_i.h:57
double flat()
Definition Pythia8_i.h:46
CLHEP::HepRandomEngine * getEngine()
Definition Pythia8_i.h:53
void init(CLHEP::HepRandomEngine *engine)
Definition Pythia8_i.h:52
unsigned long long int getRNCalls()
Definition Pythia8_i.h:54
unsigned long long int m_RNC
Definition Pythia8_i.h:56
STL class.
Author: James Monk (jmonk@cern.ch)
STL namespace.