CMS 3D CMS Logo

Py8InterfaceBase.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 
6 
9 #include <filesystem>
10 
11 // EvtGen plugin
12 //
13 //#include "Pythia8Plugins/EvtGen.h"
14 
15 using namespace Pythia8;
16 
17 namespace gen {
18 
19  Py8InterfaceBase::Py8InterfaceBase(edm::ParameterSet const& ps)
20  : BaseHadronizer(ps), useEvtGen(false), evtgenDecays(nullptr) {
21  fParameters = ps;
22 
23  pythiaPylistVerbosity = ps.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
24  pythiaHepMCVerbosity = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false);
25  pythiaHepMCVerbosityParticles = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosityParticles", false);
26  maxEventsToPrint = ps.getUntrackedParameter<int>("maxEventsToPrint", 0);
27  p8RndmEngine_ = std::make_shared<P8RndmEngine>();
28 
30  ascii_io = new HepMC::IO_AsciiParticles("cout", std::ios::out);
31 
32  if (ps.exists("useEvtGenPlugin")) {
33  useEvtGen = true;
34  auto env = std::getenv("EVTGENDATA");
35  if (not env) {
36  throw cms::Exception("EvtGenMissingEnv") << "The environment variable EVTGENDATA must be defined";
37  }
38  string evtgenpath(env);
39  evtgenDecFile = evtgenpath + string("/DECAY_2010.DEC");
40  evtgenPdlFile = evtgenpath + string("/evt.pdl");
41 
42  if (ps.exists("evtgenDecFile")) {
43  edm::FileInPath decay_table(ps.getParameter<std::string>("evtgenDecFile"));
44  evtgenDecFile = decay_table.fullPath();
45  }
46 
47  if (ps.exists("evtgenPdlFile")) {
48  edm::FileInPath pdt(ps.getParameter<std::string>("evtgenPdlFile"));
49  evtgenPdlFile = pdt.fullPath();
50  }
51 
52  if (ps.exists("evtgenUserFile")) {
53  std::vector<std::string> user_decays = ps.getParameter<std::vector<std::string> >("evtgenUserFile");
54  for (unsigned int i = 0; i < user_decays.size(); i++) {
55  edm::FileInPath user_decay(user_decays.at(i));
56  evtgenUserFiles.push_back(user_decay.fullPath());
57  }
58  //evtgenUserFiles = ps.getParameter< std::vector<std::string> >("evtgenUserFile");
59  }
60 
61  if (ps.exists("evtgenUserFileEmbedded")) {
62  std::vector<std::string> user_decay_lines =
63  ps.getParameter<std::vector<std::string> >("evtgenUserFileEmbedded");
64  char tempslhaname[] = "pythia8evtgenXXXXXX";
65  int fd = mkstemp(tempslhaname);
66 
67  for (unsigned int i = 0; i < user_decay_lines.size(); i++) {
68  user_decay_lines.at(i) += "\n";
69  write(fd, user_decay_lines.at(i).c_str(), user_decay_lines.at(i).size());
70  }
71  close(fd);
72  evtgenUserFiles.push_back(std::string(tempslhaname));
73  }
74  }
75  }
76 
78  //Pythia 8's default value for first argument to constructor
79  const string xmlDir = "../share/Pythia8/xmldoc";
80  bool printBanner = true;
81  if (fParameters.exists("printBanner")) {
82  printBanner = fParameters.getUntrackedParameter<bool>("printBanner");
83  }
84  if (!fMasterGen.get())
85  fMasterGen = std::make_unique<Pythia>(xmlDir, printBanner);
86  fDecayer = std::make_unique<Pythia>(xmlDir, printBanner);
87 
88  //add settings for resonance decay filter
89  fMasterGen->settings.addFlag("BiasedTauDecayer:filter", false);
90  fMasterGen->settings.addFlag("BiasedTauDecayer:eDecays", true);
91  fMasterGen->settings.addFlag("BiasedTauDecayer:muDecays", true);
92 
93  //add settings for resonance decay filter
94  fMasterGen->settings.addFlag("ResonanceDecayFilter:filter", false);
95  fMasterGen->settings.addFlag("ResonanceDecayFilter:exclusive", false);
96  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuAsEquivalent", false);
97  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuTauAsEquivalent", false);
98  fMasterGen->settings.addFlag("ResonanceDecayFilter:allNuAsEquivalent", false);
99  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscAsEquivalent", false);
100  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscbAsEquivalent", false);
101  fMasterGen->settings.addFlag("ResonanceDecayFilter:wzAsEquivalent", false);
102  fMasterGen->settings.addMVec("ResonanceDecayFilter:mothers", std::vector<int>(), false, false, 0, 0);
103  fMasterGen->settings.addMVec("ResonanceDecayFilter:daughters", std::vector<int>(), false, false, 0, 0);
104 
105  //add settings for PT filter
106  fMasterGen->settings.addFlag("PTFilter:filter", false);
107  fMasterGen->settings.addMode("PTFilter:quarkToFilter", 5, true, true, 3, 6);
108  fMasterGen->settings.addParm("PTFilter:scaleToFilter", 0.4, true, true, 0.0, 10.);
109  fMasterGen->settings.addParm("PTFilter:quarkRapidity", 10.0, true, true, 0.0, 10.);
110  fMasterGen->settings.addParm("PTFilter:quarkPt", -.1, true, true, -.1, 100.);
111 
112  //add settings for RecoilToTop tool
113  fMasterGen->settings.addFlag("TopRecoilHook:doTopRecoilIn", false);
114  fMasterGen->settings.addFlag("TopRecoilHook:useOldDipoleIn", false);
115  fMasterGen->settings.addFlag("TopRecoilHook:doListIn", false);
116 
117  //add settings for powheg resonance scale calculation
118  fMasterGen->settings.addFlag("POWHEGres:calcScales", false);
119  fMasterGen->settings.addFlag("POWHEG:bb4l", false);
120  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:veto", false);
121  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoQED", false);
122  fMasterGen->settings.addFlag("POWHEG:bb4l:DEBUG", false);
123  fMasterGen->settings.addFlag("POWHEG:bb4l:ScaleResonance:veto", false);
124  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoDipoleFrame", false);
125  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:pTpythiaVeto", false);
126  fMasterGen->settings.addParm("POWHEG:bb4l:pTminVeto", 10.0, true, true, 0.0, 10.);
127  fMasterGen->settings.addFlag("POWHEG:bb4l:vetoAllRadtypes", false);
128 
129  fMasterGen->setRndmEnginePtr(p8RndmEngine_);
130  fDecayer->setRndmEnginePtr(p8RndmEngine_);
131 
132  fMasterGen->readString("Next:numberShowEvent = 0");
133  fDecayer->readString("Next:numberShowEvent = 0");
134 
135  edm::ParameterSet currentParameters;
136  if (randomIndex() >= 0) {
137  std::vector<edm::ParameterSet> randomizedParameters =
138  fParameters.getParameter<std::vector<edm::ParameterSet> >("RandomizedParameters");
139  currentParameters = randomizedParameters[randomIndex()];
140  } else {
141  currentParameters = fParameters;
142  }
143 
144  ParameterCollector pCollector = currentParameters.getParameter<edm::ParameterSet>("PythiaParameters");
145 
146  for (ParameterCollector::const_iterator line = pCollector.begin(); line != pCollector.end(); ++line) {
147  if (line->find("Random:") != std::string::npos)
148  throw cms::Exception("PythiaError") << "Attempted to set random number "
149  "using Pythia commands. Please use "
150  "the RandomNumberGeneratorService."
151  << std::endl;
152 
153  if (!fMasterGen->readString(*line))
154  throw cms::Exception("PythiaError") << "Pythia 8 did not accept \"" << *line << "\"." << std::endl;
155 
156  if (line->find("ParticleDecays:") != std::string::npos) {
157  if (!fDecayer->readString(*line))
158  throw cms::Exception("PythiaError") << "Pythia 8 Decayer did not accept \"" << *line << "\"." << std::endl;
159  }
160  }
161 
162  slhafile_.clear();
163 
164  if (currentParameters.exists("SLHAFileForPythia8")) {
165  std::string slhafilenameshort = currentParameters.getParameter<std::string>("SLHAFileForPythia8");
166  edm::FileInPath f1(slhafilenameshort);
167 
168  fMasterGen->settings.mode("SLHA:readFrom", 2);
169  fMasterGen->settings.word("SLHA:file", f1.fullPath());
170  } else if (currentParameters.exists("SLHATableForPythia8")) {
171  std::string slhatable = currentParameters.getParameter<std::string>("SLHATableForPythia8");
172 
173  makeTmpSLHA(slhatable);
174  } else if (currentParameters.exists("SLHATreeForPythia8")) {
175  auto slhaReaderParams = currentParameters.getParameter<edm::ParameterSet>("SLHATreeForPythia8");
176  std::unique_ptr<SLHAReaderBase> reader =
177  SLHAReaderFactory::get()->create(slhaReaderParams.getParameter<std::string>("name"), slhaReaderParams);
178 
179  makeTmpSLHA(reader->getSLHA(currentParameters.getParameter<std::string>("ConfigDescription")));
180  }
181 
182  return true;
183  }
184 
186  char tempslhaname[] = "pythia8SLHAtableXXXXXX";
187  int fd = mkstemp(tempslhaname);
188  write(fd, slhatable.c_str(), slhatable.size());
189  close(fd);
190 
191  slhafile_ = tempslhaname;
192 
193  fMasterGen->settings.mode("SLHA:readFrom", 2);
194  fMasterGen->settings.word("SLHA:file", slhafile_);
195  }
196 
197  bool Py8InterfaceBase::declareStableParticles(const std::vector<int>& pdgIds) {
198  for (size_t i = 0; i < pdgIds.size(); i++) {
199  // FIXME: need to double check if PID's are the same in Py6 & Py8,
200  // because the HepPDT translation tool is actually for **Py6**
201  //
202  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
203  //
204  // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
205  int PyID = pdgIds[i];
206  std::ostringstream pyCard;
207  pyCard << PyID << ":mayDecay=false";
208 
209  if (fMasterGen->particleData.isParticle(PyID)) {
210  fMasterGen->readString(pyCard.str());
211  } else {
212  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
213  << "recognize particle id = " << PyID << std::endl;
214  }
215  // alternative:
216  // set the 2nd input argument warn=false
217  // - this way Py8 will NOT print warnings about unknown particle code(s)
218  // fMasterPtr->readString( pyCard.str(), false )
219  }
220 
221  return true;
222  }
223 
224  bool Py8InterfaceBase::declareSpecialSettings(const std::vector<std::string>& settings) {
225  for (unsigned int iss = 0; iss < settings.size(); iss++) {
226  if (settings[iss].find("QED-brem-off") != std::string::npos) {
227  fMasterGen->readString("TimeShower:QEDshowerByL=off");
228  } else {
229  size_t fnd1 = settings[iss].find("Pythia8:");
230  if (fnd1 != std::string::npos) {
231  std::string value = settings[iss].substr(fnd1 + 8);
232  fDecayer->readString(value);
233  }
234  }
235  }
236  return true;
237  }
238 
240  fMasterGen->stat();
241  return;
242  }
243 
244 } // namespace gen
bool declareSpecialSettings(const std::vector< std::string > &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::shared_ptr< P8RndmEngine > p8RndmEngine_
std::unique_ptr< Pythia8::Pythia > fDecayer
void makeTmpSLHA(const std::string &)
HepMC::IO_AsciiParticles * ascii_io
const_iterator begin() const
bool exists(std::string const &parameterName) const
checks if a parameter exists
const_iterator end() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< std::string > evtgenUserFiles
T getUntrackedParameter(std::string const &, T const &) const
int randomIndex() const
unsigned int pythiaPylistVerbosity
Definition: value.py:1
bool declareStableParticles(const std::vector< int > &)
unsigned int maxEventsToPrint
std::unique_ptr< Pythia8::Pythia > fMasterGen
const std::string & fullPath() const
Definition: FileInPath.cc:144
#define get
Log< level::Warning, false > LogWarning
edm::ParameterSet fParameters
fd
Definition: ztee.py:136