CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes | Private Attributes
gen::Py8InterfaceBase Class Referenceabstract

#include <Py8InterfaceBase.h>

Inheritance diagram for gen::Py8InterfaceBase:
gen::BaseHadronizer gen::Py8GunBase Pythia8Hadronizer gen::Py8EGun gen::Py8JetGun gen::Py8PtGun

Public Member Functions

virtual const char * classname () const =0
 
bool decay ()
 
bool declareSpecialSettings (const std::vector< std::string > &)
 
bool declareStableParticles (const std::vector< int > &)
 
virtual void finalizeEvent ()=0
 
virtual bool generatePartonsAndHadronize ()=0
 
virtual bool initializeForInternalPartons ()=0
 
void p8SetRandomEngine (CLHEP::HepRandomEngine *v)
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
P8RndmEnginerandomEngine ()
 
bool readSettings (int)
 
virtual void statistics ()
 
 ~Py8InterfaceBase () override
 
- Public Member Functions inherited from gen::BaseHadronizer
 BaseHadronizer (edm::ParameterSet const &ps)
 
void cleanLHE ()
 
void generateLHE (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine, unsigned int ncpu)
 
edm::EventgetEDMEvent () const
 
std::unique_ptr< HepMC::GenEventgetGenEvent ()
 
std::unique_ptr< GenEventInfoProductgetGenEventInfo ()
 
virtual std::unique_ptr< GenLumiInfoHeadergetGenLumiInfoHeader () const
 
GenRunInfoProductgetGenRunInfo ()
 
std::unique_ptr< lhef::LHEEventgetLHEEvent ()
 
const std::shared_ptr< lhef::LHERunInfo > & getLHERunInfo () const
 
const std::string & gridpackPath () const
 
int randomIndex () const
 
const std::string & randomInitConfigDescription () const
 
void randomizeIndex (edm::LuminosityBlock const &lumi, CLHEP::HepRandomEngine *rengine)
 
void resetEvent (std::unique_ptr< HepMC::GenEvent > event)
 
void resetEventInfo (std::unique_ptr< GenEventInfoProduct > eventInfo)
 
virtual bool select (HepMC::GenEvent *) const
 
void setEDMEvent (edm::Event &event)
 
void setLHEEvent (std::unique_ptr< lhef::LHEEvent > event)
 
void setLHERunInfo (std::unique_ptr< lhef::LHERunInfo > runInfo)
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
virtual ~BaseHadronizer ()(false)
 

Protected Attributes

HepMC::IO_AsciiParticles * ascii_io
 
std::shared_ptr< EvtGenDecays > evtgenDecays
 
std::string evtgenDecFile
 
std::string evtgenPdlFile
 
std::vector< std::string > evtgenUserFiles
 
std::unique_ptr< Pythia8::Pythia > fDecayer
 
std::unique_ptr< Pythia8::Pythia > fMasterGen
 
edm::ParameterSet fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
bool pythiaHepMCVerbosityParticles
 
unsigned int pythiaPylistVerbosity
 
std::string slhafile_
 
HepMC::Pythia8ToHepMC toHepMC
 
bool useEvtGen
 
- Protected Attributes inherited from gen::BaseHadronizer
std::string lheFile_
 
int randomIndex_
 

Private Attributes

P8RndmEngine p8RndmEngine_
 

Additional Inherited Members

- Protected Member Functions inherited from gen::BaseHadronizer
std::unique_ptr< HepMC::GenEvent > & event ()
 
std::unique_ptr< GenEventInfoProduct > & eventInfo ()
 
lhef::LHEEventlheEvent ()
 
lhef::LHERunInfolheRunInfo ()
 
GenRunInfoProductrunInfo ()
 

Detailed Description

Definition at line 25 of file Py8InterfaceBase.h.

Constructor & Destructor Documentation

gen::Py8InterfaceBase::Py8InterfaceBase ( edm::ParameterSet const &  ps)

Definition at line 17 of file Py8InterfaceBase.cc.

References ascii_io, evtgenDecFile, evtgenPdlFile, evtgenUserFiles, edm::ParameterSet::exists(), cmsRelvalreport::exit, fParameters, edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), mps_fire::i, maxEventsToPrint, MillePedeFileConverter_cfg::out, pythiaHepMCVerbosity, pythiaHepMCVerbosityParticles, pythiaPylistVerbosity, AlCaHLTBitMon_QueryRunRegistry::string, and useEvtGen.

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  }
std::shared_ptr< EvtGenDecays > evtgenDecays
HepMC::IO_AsciiParticles * ascii_io
BaseHadronizer(edm::ParameterSet const &ps)
std::vector< std::string > evtgenUserFiles
unsigned int pythiaPylistVerbosity
unsigned int maxEventsToPrint
edm::ParameterSet fParameters
gen::Py8InterfaceBase::~Py8InterfaceBase ( )
inlineoverride

Definition at line 30 of file Py8InterfaceBase.h.

30 {}

Member Function Documentation

virtual const char* gen::Py8InterfaceBase::classname ( ) const
pure virtual
bool gen::Py8InterfaceBase::decay ( )
inline

Definition at line 33 of file Py8InterfaceBase.h.

References data-class-funcs::classname, and myMessageLogger_cff::statistics.

33 { return true; } // NOT used - let's call it "design imperfection"
bool gen::Py8InterfaceBase::declareSpecialSettings ( const std::vector< std::string > &  settings)

Definition at line 207 of file Py8InterfaceBase.cc.

References fDecayer, spr::find(), fMasterGen, and AlCaHLTBitMon_QueryRunRegistry::string.

207  {
208  for (unsigned int iss = 0; iss < settings.size(); iss++) {
209  if (settings[iss].find("QED-brem-off") != std::string::npos) {
210  fMasterGen->readString("TimeShower:QEDshowerByL=off");
211  } else {
212  size_t fnd1 = settings[iss].find("Pythia8:");
213  if (fnd1 != std::string::npos) {
214  std::string value = settings[iss].substr(fnd1 + 8);
215  fDecayer->readString(value);
216  }
217  }
218  }
219  return true;
220  }
std::unique_ptr< Pythia8::Pythia > fDecayer
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
Definition: value.py:1
std::unique_ptr< Pythia8::Pythia > fMasterGen
bool gen::Py8InterfaceBase::declareStableParticles ( const std::vector< int > &  pdgIds)

Definition at line 180 of file Py8InterfaceBase.cc.

References fMasterGen, and mps_fire::i.

180  {
181  for (size_t i = 0; i < pdgIds.size(); i++) {
182  // FIXME: need to double check if PID's are the same in Py6 & Py8,
183  // because the HepPDT translation tool is actually for **Py6**
184  //
185  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
186  //
187  // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
188  int PyID = pdgIds[i];
189  std::ostringstream pyCard;
190  pyCard << PyID << ":mayDecay=false";
191 
192  if (fMasterGen->particleData.isParticle(PyID)) {
193  fMasterGen->readString(pyCard.str());
194  } else {
195  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
196  << "recognize particle id = " << PyID << std::endl;
197  }
198  // alternative:
199  // set the 2nd input argument warn=false
200  // - this way Py8 will NOT print warnings about unknown particle code(s)
201  // fMasterPtr->readString( pyCard.str(), false )
202  }
203 
204  return true;
205  }
std::unique_ptr< Pythia8::Pythia > fMasterGen
virtual void gen::Py8InterfaceBase::finalizeEvent ( )
pure virtual

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

virtual bool gen::Py8InterfaceBase::generatePartonsAndHadronize ( )
pure virtual
virtual bool gen::Py8InterfaceBase::initializeForInternalPartons ( )
pure virtual

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

void gen::Py8InterfaceBase::p8SetRandomEngine ( CLHEP::HepRandomEngine *  v)
inline

Definition at line 42 of file Py8InterfaceBase.h.

double v[5][pyjets_maxn]
void setRandomEngine(CLHEP::HepRandomEngine *v)
Definition: P8RndmEngine.h:35
P8RndmEngine p8RndmEngine_
P8RndmEngine& gen::Py8InterfaceBase::randomEngine ( )
inline
bool gen::Py8InterfaceBase::readSettings ( int  )

Definition at line 78 of file Py8InterfaceBase.cc.

References gen::ParameterCollector::begin(), gen::ParameterCollector::end(), edm::ParameterSet::exists(), connectstrParser::f1, fDecayer, fMasterGen, fParameters, edm::ParameterSet::getParameter(), mps_splice::line, p8RndmEngine_, gen::BaseHadronizer::randomIndex(), slhafile_, AlCaHLTBitMon_QueryRunRegistry::string, and TriggerAnalyzer::write().

78  {
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:eMuDecays", true);
86 
87  //add settings for resonance decay filter
88  fMasterGen->settings.addFlag("ResonanceDecayFilter:filter", false);
89  fMasterGen->settings.addFlag("ResonanceDecayFilter:exclusive", false);
90  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuAsEquivalent", false);
91  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuTauAsEquivalent", false);
92  fMasterGen->settings.addFlag("ResonanceDecayFilter:allNuAsEquivalent", false);
93  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscAsEquivalent", false);
94  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscbAsEquivalent", false);
95  fMasterGen->settings.addFlag("ResonanceDecayFilter:wzAsEquivalent", false);
96  fMasterGen->settings.addMVec("ResonanceDecayFilter:mothers", std::vector<int>(), false, false, 0, 0);
97  fMasterGen->settings.addMVec("ResonanceDecayFilter:daughters", std::vector<int>(), false, false, 0, 0);
98 
99  //add settings for PT filter
100  fMasterGen->settings.addFlag("PTFilter:filter", false);
101  fMasterGen->settings.addMode("PTFilter:quarkToFilter", 5, true, true, 3, 6);
102  fMasterGen->settings.addParm("PTFilter:scaleToFilter", 0.4, true, true, 0.0, 10.);
103  fMasterGen->settings.addParm("PTFilter:quarkRapidity", 10.0, true, true, 0.0, 10.);
104  fMasterGen->settings.addParm("PTFilter:quarkPt", -.1, true, true, -.1, 100.);
105 
106  //add settings for powheg resonance scale calculation
107  fMasterGen->settings.addFlag("POWHEGres:calcScales", false);
108  fMasterGen->settings.addFlag("POWHEG:bb4l", false);
109  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:onlyDistance1", false);
110  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:veto", false);
111  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:dryRun", false);
112  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoAtPL", false);
113  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoQED", false);
114  fMasterGen->settings.addFlag("POWHEG:bb4l:PartonLevel:veto", false);
115  fMasterGen->settings.addFlag("POWHEG:bb4l:PartonLevel:excludeFSRConflicting", false);
116  fMasterGen->settings.addFlag("POWHEG:bb4l:DEBUG", false);
117  fMasterGen->settings.addFlag("POWHEG:bb4l:ScaleResonance:veto", false);
118  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoDipoleFrame", false);
119  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:pTpythiaVeto", false);
120  fMasterGen->settings.addParm("POWHEG:bb4l:pTminVeto", 10.0, true, true, 0.0, 10.);
121 
122  fMasterGen->setRndmEnginePtr(&p8RndmEngine_);
123  fDecayer->setRndmEnginePtr(&p8RndmEngine_);
124 
125  fMasterGen->readString("Next:numberShowEvent = 0");
126  fDecayer->readString("Next:numberShowEvent = 0");
127 
128  edm::ParameterSet currentParameters;
129  if (randomIndex() >= 0) {
130  std::vector<edm::ParameterSet> randomizedParameters =
131  fParameters.getParameter<std::vector<edm::ParameterSet> >("RandomizedParameters");
132  currentParameters = randomizedParameters[randomIndex()];
133  } else {
134  currentParameters = fParameters;
135  }
136 
137  ParameterCollector pCollector = currentParameters.getParameter<edm::ParameterSet>("PythiaParameters");
138 
139  for (ParameterCollector::const_iterator line = pCollector.begin(); line != pCollector.end(); ++line) {
140  if (line->find("Random:") != std::string::npos)
141  throw cms::Exception("PythiaError") << "Attempted to set random number "
142  "using Pythia commands. Please use "
143  "the RandomNumberGeneratorService."
144  << std::endl;
145 
146  if (!fMasterGen->readString(*line))
147  throw cms::Exception("PythiaError") << "Pythia 8 did not accept \"" << *line << "\"." << std::endl;
148 
149  if (line->find("ParticleDecays:") != std::string::npos) {
150  if (!fDecayer->readString(*line))
151  throw cms::Exception("PythiaError") << "Pythia 8 Decayer did not accept \"" << *line << "\"." << std::endl;
152  }
153  }
154 
155  slhafile_.clear();
156 
157  if (currentParameters.exists("SLHAFileForPythia8")) {
158  std::string slhafilenameshort = currentParameters.getParameter<std::string>("SLHAFileForPythia8");
159  edm::FileInPath f1(slhafilenameshort);
160 
161  fMasterGen->settings.mode("SLHA:readFrom", 2);
162  fMasterGen->settings.word("SLHA:file", f1.fullPath());
163  } else if (currentParameters.exists("SLHATableForPythia8")) {
164  std::string slhatable = currentParameters.getParameter<std::string>("SLHATableForPythia8");
165 
166  char tempslhaname[] = "pythia8SLHAtableXXXXXX";
167  int fd = mkstemp(tempslhaname);
168  write(fd, slhatable.c_str(), slhatable.size());
169  close(fd);
170 
171  slhafile_ = tempslhaname;
172 
173  fMasterGen->settings.mode("SLHA:readFrom", 2);
174  fMasterGen->settings.word("SLHA:file", slhafile_);
175  }
176 
177  return true;
178  }
T getParameter(std::string const &) const
std::unique_ptr< Pythia8::Pythia > fDecayer
bool exists(std::string const &parameterName) const
checks if a parameter exists
P8RndmEngine p8RndmEngine_
int randomIndex() const
std::unique_ptr< Pythia8::Pythia > fMasterGen
def write(self, setup)
edm::ParameterSet fParameters
void gen::Py8InterfaceBase::statistics ( )
virtual

Reimplemented in Pythia8Hadronizer, and gen::Py8GunBase.

Definition at line 222 of file Py8InterfaceBase.cc.

References fMasterGen.

222  {
223  fMasterGen->stat();
224  return;
225  }
std::unique_ptr< Pythia8::Pythia > fMasterGen

Member Data Documentation

HepMC::IO_AsciiParticles* gen::Py8InterfaceBase::ascii_io
protected
std::shared_ptr<EvtGenDecays> gen::Py8InterfaceBase::evtgenDecays
protected
std::string gen::Py8InterfaceBase::evtgenDecFile
protected
std::string gen::Py8InterfaceBase::evtgenPdlFile
protected
std::vector<std::string> gen::Py8InterfaceBase::evtgenUserFiles
protected
std::unique_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fDecayer
protected
std::unique_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fMasterGen
protected
edm::ParameterSet gen::Py8InterfaceBase::fParameters
protected

Definition at line 51 of file Py8InterfaceBase.h.

Referenced by Py8InterfaceBase(), and readSettings().

unsigned int gen::Py8InterfaceBase::maxEventsToPrint
protected
P8RndmEngine gen::Py8InterfaceBase::p8RndmEngine_
private

Definition at line 71 of file Py8InterfaceBase.h.

Referenced by readSettings().

bool gen::Py8InterfaceBase::pythiaHepMCVerbosity
protected
bool gen::Py8InterfaceBase::pythiaHepMCVerbosityParticles
protected
unsigned int gen::Py8InterfaceBase::pythiaPylistVerbosity
protected
std::string gen::Py8InterfaceBase::slhafile_
protected
HepMC::Pythia8ToHepMC gen::Py8InterfaceBase::toHepMC
protected
bool gen::Py8InterfaceBase::useEvtGen
protected