CMS 3D CMS Logo

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