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::Py8MassGun gen::Py8PtAndDxyGun gen::Py8PtGun gen::Py8PtotGun

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 makeTmpSLHA (const std::string &)
 
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 () noexcept(false)
 

Protected Attributes

HepMC::IO_AsciiParticles * ascii_io
 
std::shared_ptr< Pythia8::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 27 of file Py8InterfaceBase.h.

Constructor & Destructor Documentation

◆ Py8InterfaceBase()

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

Definition at line 19 of file Py8InterfaceBase.cc.

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 
29  ascii_io = new HepMC::IO_AsciiParticles("cout", std::ios::out);
30 
31  if (ps.exists("useEvtGenPlugin")) {
32  useEvtGen = true;
33  string evtgenpath(std::getenv("EVTGENDATA"));
34  evtgenDecFile = evtgenpath + string("/DECAY_2010.DEC");
35  evtgenPdlFile = evtgenpath + string("/evt.pdl");
36 
37  if (ps.exists("evtgenDecFile")) {
38  edm::FileInPath decay_table(ps.getParameter<std::string>("evtgenDecFile"));
39  evtgenDecFile = decay_table.fullPath();
40  }
41 
42  if (ps.exists("evtgenPdlFile")) {
43  edm::FileInPath pdt(ps.getParameter<std::string>("evtgenPdlFile"));
44  evtgenPdlFile = pdt.fullPath();
45  }
46 
47  if (ps.exists("evtgenUserFile")) {
48  std::vector<std::string> user_decays = ps.getParameter<std::vector<std::string> >("evtgenUserFile");
49  for (unsigned int i = 0; i < user_decays.size(); i++) {
50  edm::FileInPath user_decay(user_decays.at(i));
51  evtgenUserFiles.push_back(user_decay.fullPath());
52  }
53  //evtgenUserFiles = ps.getParameter< std::vector<std::string> >("evtgenUserFile");
54  }
55 
56  if (ps.exists("evtgenUserFileEmbedded")) {
57  std::vector<std::string> user_decay_lines =
58  ps.getParameter<std::vector<std::string> >("evtgenUserFileEmbedded");
59  char tempslhaname[] = "pythia8evtgenXXXXXX";
60  int fd = mkstemp(tempslhaname);
61 
62  for (unsigned int i = 0; i < user_decay_lines.size(); i++) {
63  user_decay_lines.at(i) += "\n";
64  write(fd, user_decay_lines.at(i).c_str(), user_decay_lines.at(i).size());
65  }
66  close(fd);
67  evtgenUserFiles.push_back(std::string(tempslhaname));
68  }
69  }
70  }

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

◆ ~Py8InterfaceBase()

gen::Py8InterfaceBase::~Py8InterfaceBase ( )
inlineoverride

Definition at line 30 of file Py8InterfaceBase.h.

30 {}

Member Function Documentation

◆ classname()

virtual const char* gen::Py8InterfaceBase::classname ( ) const
pure virtual

◆ decay()

bool gen::Py8InterfaceBase::decay ( )
inline

Definition at line 33 of file Py8InterfaceBase.h.

33 { return true; } // NOT used - let's call it "design imperfection"

◆ declareSpecialSettings()

bool gen::Py8InterfaceBase::declareSpecialSettings ( const std::vector< std::string > &  settings)

Definition at line 212 of file Py8InterfaceBase.cc.

212  {
213  for (unsigned int iss = 0; iss < settings.size(); iss++) {
214  if (settings[iss].find("QED-brem-off") != std::string::npos) {
215  fMasterGen->readString("TimeShower:QEDshowerByL=off");
216  } else {
217  size_t fnd1 = settings[iss].find("Pythia8:");
218  if (fnd1 != std::string::npos) {
219  std::string value = settings[iss].substr(fnd1 + 8);
220  fDecayer->readString(value);
221  }
222  }
223  }
224  return true;
225  }

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

◆ declareStableParticles()

bool gen::Py8InterfaceBase::declareStableParticles ( const std::vector< int > &  pdgIds)

Definition at line 185 of file Py8InterfaceBase.cc.

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

References fMasterGen, mps_fire::i, and CosmicGenFilterHelix_cfi::pdgIds.

◆ finalizeEvent()

virtual void gen::Py8InterfaceBase::finalizeEvent ( )
pure virtual

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

◆ generatePartonsAndHadronize()

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

◆ initializeForInternalPartons()

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

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

◆ makeTmpSLHA()

void gen::Py8InterfaceBase::makeTmpSLHA ( const std::string &  slhatable)

Definition at line 173 of file Py8InterfaceBase.cc.

173  {
174  char tempslhaname[] = "pythia8SLHAtableXXXXXX";
175  int fd = mkstemp(tempslhaname);
176  write(fd, slhatable.c_str(), slhatable.size());
177  close(fd);
178 
179  slhafile_ = tempslhaname;
180 
181  fMasterGen->settings.mode("SLHA:readFrom", 2);
182  fMasterGen->settings.word("SLHA:file", slhafile_);
183  }

References ztee::fd, fMasterGen, slhafile_, and writeEcalDQMStatus::write.

Referenced by readSettings().

◆ p8SetRandomEngine()

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

◆ randomEngine()

P8RndmEngine& gen::Py8InterfaceBase::randomEngine ( )
inline

◆ readSettings()

bool gen::Py8InterfaceBase::readSettings ( int  )

Definition at line 72 of file Py8InterfaceBase.cc.

72  {
73  if (!fMasterGen.get())
74  fMasterGen = std::make_unique<Pythia>();
75  fDecayer = std::make_unique<Pythia>();
76 
77  //add settings for resonance decay filter
78  fMasterGen->settings.addFlag("BiasedTauDecayer:filter", false);
79  fMasterGen->settings.addFlag("BiasedTauDecayer:eDecays", true);
80  fMasterGen->settings.addFlag("BiasedTauDecayer:muDecays", true);
81 
82  //add settings for resonance decay filter
83  fMasterGen->settings.addFlag("ResonanceDecayFilter:filter", false);
84  fMasterGen->settings.addFlag("ResonanceDecayFilter:exclusive", false);
85  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuAsEquivalent", false);
86  fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuTauAsEquivalent", false);
87  fMasterGen->settings.addFlag("ResonanceDecayFilter:allNuAsEquivalent", false);
88  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscAsEquivalent", false);
89  fMasterGen->settings.addFlag("ResonanceDecayFilter:udscbAsEquivalent", false);
90  fMasterGen->settings.addFlag("ResonanceDecayFilter:wzAsEquivalent", false);
91  fMasterGen->settings.addMVec("ResonanceDecayFilter:mothers", std::vector<int>(), false, false, 0, 0);
92  fMasterGen->settings.addMVec("ResonanceDecayFilter:daughters", std::vector<int>(), false, false, 0, 0);
93 
94  //add settings for PT filter
95  fMasterGen->settings.addFlag("PTFilter:filter", false);
96  fMasterGen->settings.addMode("PTFilter:quarkToFilter", 5, true, true, 3, 6);
97  fMasterGen->settings.addParm("PTFilter:scaleToFilter", 0.4, true, true, 0.0, 10.);
98  fMasterGen->settings.addParm("PTFilter:quarkRapidity", 10.0, true, true, 0.0, 10.);
99  fMasterGen->settings.addParm("PTFilter:quarkPt", -.1, true, true, -.1, 100.);
100 
101  //add settings for powheg resonance scale calculation
102  fMasterGen->settings.addFlag("POWHEGres:calcScales", false);
103  fMasterGen->settings.addFlag("POWHEG:bb4l", false);
104  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:onlyDistance1", false);
105  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:veto", false);
106  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:dryRun", false);
107  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoAtPL", false);
108  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoQED", false);
109  fMasterGen->settings.addFlag("POWHEG:bb4l:PartonLevel:veto", false);
110  fMasterGen->settings.addFlag("POWHEG:bb4l:PartonLevel:excludeFSRConflicting", false);
111  fMasterGen->settings.addFlag("POWHEG:bb4l:DEBUG", false);
112  fMasterGen->settings.addFlag("POWHEG:bb4l:ScaleResonance:veto", false);
113  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoDipoleFrame", false);
114  fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:pTpythiaVeto", false);
115  fMasterGen->settings.addParm("POWHEG:bb4l:pTminVeto", 10.0, true, true, 0.0, 10.);
116 
117  fMasterGen->setRndmEnginePtr(&p8RndmEngine_);
118  fDecayer->setRndmEnginePtr(&p8RndmEngine_);
119 
120  fMasterGen->readString("Next:numberShowEvent = 0");
121  fDecayer->readString("Next:numberShowEvent = 0");
122 
123  edm::ParameterSet currentParameters;
124  if (randomIndex() >= 0) {
125  std::vector<edm::ParameterSet> randomizedParameters =
126  fParameters.getParameter<std::vector<edm::ParameterSet> >("RandomizedParameters");
127  currentParameters = randomizedParameters[randomIndex()];
128  } else {
129  currentParameters = fParameters;
130  }
131 
132  ParameterCollector pCollector = currentParameters.getParameter<edm::ParameterSet>("PythiaParameters");
133 
134  for (ParameterCollector::const_iterator line = pCollector.begin(); line != pCollector.end(); ++line) {
135  if (line->find("Random:") != std::string::npos)
136  throw cms::Exception("PythiaError") << "Attempted to set random number "
137  "using Pythia commands. Please use "
138  "the RandomNumberGeneratorService."
139  << std::endl;
140 
141  if (!fMasterGen->readString(*line))
142  throw cms::Exception("PythiaError") << "Pythia 8 did not accept \"" << *line << "\"." << std::endl;
143 
144  if (line->find("ParticleDecays:") != std::string::npos) {
145  if (!fDecayer->readString(*line))
146  throw cms::Exception("PythiaError") << "Pythia 8 Decayer did not accept \"" << *line << "\"." << std::endl;
147  }
148  }
149 
150  slhafile_.clear();
151 
152  if (currentParameters.exists("SLHAFileForPythia8")) {
153  std::string slhafilenameshort = currentParameters.getParameter<std::string>("SLHAFileForPythia8");
154  edm::FileInPath f1(slhafilenameshort);
155 
156  fMasterGen->settings.mode("SLHA:readFrom", 2);
157  fMasterGen->settings.word("SLHA:file", f1.fullPath());
158  } else if (currentParameters.exists("SLHATableForPythia8")) {
159  std::string slhatable = currentParameters.getParameter<std::string>("SLHATableForPythia8");
160 
161  makeTmpSLHA(slhatable);
162  } else if (currentParameters.exists("SLHATreeForPythia8")) {
163  auto slhaReaderParams = currentParameters.getParameter<edm::ParameterSet>("SLHATreeForPythia8");
164  std::unique_ptr<SLHAReaderBase> reader =
165  SLHAReaderFactory::get()->create(slhaReaderParams.getParameter<std::string>("name"), slhaReaderParams);
166 
167  makeTmpSLHA(reader->getSLHA(currentParameters.getParameter<std::string>("ConfigDescription")));
168  }
169 
170  return true;
171  }

References gen::ParameterCollector::begin(), gen::ParameterCollector::end(), edm::ParameterSet::exists(), DeadROC_duringRun::f1, fDecayer, fMasterGen, fParameters, get, edm::ParameterSet::getParameter(), mps_splice::line, makeTmpSLHA(), p8RndmEngine_, gen::BaseHadronizer::randomIndex(), DQM::reader, slhafile_, and AlCaHLTBitMon_QueryRunRegistry::string.

◆ statistics()

void gen::Py8InterfaceBase::statistics ( )
virtual

Reimplemented in Pythia8Hadronizer, and gen::Py8GunBase.

Definition at line 227 of file Py8InterfaceBase.cc.

227  {
228  fMasterGen->stat();
229  return;
230  }

References fMasterGen.

Member Data Documentation

◆ ascii_io

HepMC::IO_AsciiParticles* gen::Py8InterfaceBase::ascii_io
protected

◆ evtgenDecays

std::shared_ptr<Pythia8::EvtGenDecays> gen::Py8InterfaceBase::evtgenDecays
protected

◆ evtgenDecFile

std::string gen::Py8InterfaceBase::evtgenDecFile
protected

◆ evtgenPdlFile

std::string gen::Py8InterfaceBase::evtgenPdlFile
protected

◆ evtgenUserFiles

std::vector<std::string> gen::Py8InterfaceBase::evtgenUserFiles
protected

◆ fDecayer

std::unique_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fDecayer
protected

◆ fMasterGen

std::unique_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fMasterGen
protected

◆ fParameters

edm::ParameterSet gen::Py8InterfaceBase::fParameters
protected

Definition at line 51 of file Py8InterfaceBase.h.

Referenced by Py8InterfaceBase(), and readSettings().

◆ maxEventsToPrint

unsigned int gen::Py8InterfaceBase::maxEventsToPrint
protected

◆ p8RndmEngine_

P8RndmEngine gen::Py8InterfaceBase::p8RndmEngine_
private

Definition at line 70 of file Py8InterfaceBase.h.

Referenced by p8SetRandomEngine(), randomEngine(), and readSettings().

◆ pythiaHepMCVerbosity

bool gen::Py8InterfaceBase::pythiaHepMCVerbosity
protected

◆ pythiaHepMCVerbosityParticles

bool gen::Py8InterfaceBase::pythiaHepMCVerbosityParticles
protected

◆ pythiaPylistVerbosity

unsigned int gen::Py8InterfaceBase::pythiaPylistVerbosity
protected

◆ slhafile_

std::string gen::Py8InterfaceBase::slhafile_
protected

◆ toHepMC

HepMC::Pythia8ToHepMC gen::Py8InterfaceBase::toHepMC
protected

◆ useEvtGen

bool gen::Py8InterfaceBase::useEvtGen
protected
gen::Py8InterfaceBase::fMasterGen
std::unique_ptr< Pythia8::Pythia > fMasterGen
Definition: Py8InterfaceBase.h:47
gen::Py8InterfaceBase::evtgenDecays
std::shared_ptr< Pythia8::EvtGenDecays > evtgenDecays
Definition: Py8InterfaceBase.h:62
mps_fire.i
i
Definition: mps_fire.py:428
gen::Py8InterfaceBase::makeTmpSLHA
void makeTmpSLHA(const std::string &)
Definition: Py8InterfaceBase.cc:173
gen::Py8InterfaceBase::ascii_io
HepMC::IO_AsciiParticles * ascii_io
Definition: Py8InterfaceBase.h:57
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
gen::Py8InterfaceBase::p8RndmEngine_
P8RndmEngine p8RndmEngine_
Definition: Py8InterfaceBase.h:70
edm::FileInPath
Definition: FileInPath.h:61
ztee.fd
fd
Definition: ztee.py:136
DQM.reader
reader
Definition: DQM.py:105
gen::Py8InterfaceBase::useEvtGen
bool useEvtGen
Definition: Py8InterfaceBase.h:61
gen::Py8InterfaceBase::pythiaHepMCVerbosityParticles
bool pythiaHepMCVerbosityParticles
Definition: Py8InterfaceBase.h:55
gen::Py8InterfaceBase::evtgenUserFiles
std::vector< std::string > evtgenUserFiles
Definition: Py8InterfaceBase.h:65
gen::Py8InterfaceBase::pythiaPylistVerbosity
unsigned int pythiaPylistVerbosity
Definition: Py8InterfaceBase.h:53
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
edm::ParameterSet
Definition: ParameterSet.h:47
gen::Py8InterfaceBase::evtgenPdlFile
std::string evtgenPdlFile
Definition: Py8InterfaceBase.h:64
gen::Py8InterfaceBase::maxEventsToPrint
unsigned int maxEventsToPrint
Definition: Py8InterfaceBase.h:56
gen::BaseHadronizer::randomIndex
int randomIndex() const
Definition: BaseHadronizer.h:76
gen::v
double v[5][pyjets_maxn]
Definition: Cascade2Hadronizer.cc:76
CosmicGenFilterHelix_cfi.pdgIds
pdgIds
Definition: CosmicGenFilterHelix_cfi.py:5
value
Definition: value.py:1
get
#define get
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
writeEcalDQMStatus.write
write
Definition: writeEcalDQMStatus.py:48
gen::Py8InterfaceBase::pythiaHepMCVerbosity
bool pythiaHepMCVerbosity
Definition: Py8InterfaceBase.h:54
gen::BaseHadronizer::BaseHadronizer
BaseHadronizer(edm::ParameterSet const &ps)
Definition: BaseHadronizer.cc:12
gen::P8RndmEngine::setRandomEngine
void setRandomEngine(CLHEP::HepRandomEngine *v)
Definition: P8RndmEngine.h:34
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
gen::Py8InterfaceBase::fParameters
edm::ParameterSet fParameters
Definition: Py8InterfaceBase.h:51
cms::Exception
Definition: Exception.h:70
DeadROC_duringRun.f1
f1
Definition: DeadROC_duringRun.py:219
mps_splice.line
line
Definition: mps_splice.py:76
gen::Py8InterfaceBase::fDecayer
std::unique_ptr< Pythia8::Pythia > fDecayer
Definition: Py8InterfaceBase.h:48
gen::Py8InterfaceBase::evtgenDecFile
std::string evtgenDecFile
Definition: Py8InterfaceBase.h:63
gen::Py8InterfaceBase::slhafile_
std::string slhafile_
Definition: Py8InterfaceBase.h:67
gen::ParameterCollector::const_iterator
friend class const_iterator
Definition: ParameterCollector.h:83