CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Py8InterfaceBase.cc
Go to the documentation of this file.
2 
5 
6 // EvtGen plugin
7 //
8 //#include "Pythia8Plugins/EvtGen.h"
9 
10 using namespace Pythia8;
11 
12 namespace gen {
13 
14 Py8InterfaceBase::Py8InterfaceBase( edm::ParameterSet const& ps ) :
15 useEvtGen(false), evtgenDecays(0)
16 {
17  fMasterGen.reset(new Pythia);
18  fDecayer.reset(new Pythia);
19 
20  fMasterGen->readString("Next:numberShowEvent = 0");
21  fDecayer->readString("Next:numberShowEvent = 0");
22 
23  fMasterGen->setRndmEnginePtr( &p8RndmEngine_ );
24  fDecayer->setRndmEnginePtr( &p8RndmEngine_ );
25 
26  fParameters = ps.getParameter<edm::ParameterSet>("PythiaParameters");
27 
28  pythiaPylistVerbosity = ps.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
29  pythiaHepMCVerbosity = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false);
30  pythiaHepMCVerbosityParticles = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosityParticles", false);
31  maxEventsToPrint = ps.getUntrackedParameter<int>("maxEventsToPrint", 0);
32 
33  if(pythiaHepMCVerbosityParticles)
34  ascii_io = new HepMC::IO_AsciiParticles("cout", std::ios::out);
35 
36  if ( ps.exists("useEvtGenPlugin") ) {
37 
38  useEvtGen = true;
39 
40  string evtgenpath(getenv("EVTGENDATA"));
41  evtgenDecFile = evtgenpath + string("/DECAY_2010.DEC");
42  evtgenPdlFile = evtgenpath + string("/evt.pdl");
43 
44  if ( ps.exists( "evtgenDecFile" ) )
45  evtgenDecFile = ps.getParameter<string>("evtgenDecFile");
46 
47  if ( ps.exists( "evtgenPdlFile" ) )
48  evtgenPdlFile = ps.getParameter<string>("evtgenPdlFile");
49 
50  if ( ps.exists( "evtgenUserFile" ) )
51  evtgenUserFiles = ps.getParameter< std::vector<std::string> >("evtgenUserFile");
52 
53  }
54 
55 }
56 
58 {
59 
61  line != fParameters.end(); ++line )
62  {
63  if (line->find("Random:") != std::string::npos)
64  throw cms::Exception("PythiaError") << "Attempted to set random number "
65  "using Pythia commands. Please use " "the RandomNumberGeneratorService."
66  << std::endl;
67 
68  if (!fMasterGen->readString(*line)) throw cms::Exception("PythiaError")
69  << "Pythia 8 did not accept \""
70  << *line << "\"." << std::endl;
71 
72  if (line->find("ParticleDecays:") != std::string::npos) {
73 
74  if (!fDecayer->readString(*line)) throw cms::Exception("PythiaError")
75  << "Pythia 8 Decayer did not accept \""
76  << *line << "\"." << std::endl;
77  }
78 
79  }
80 
81  return true;
82 
83 }
84 
85 bool Py8InterfaceBase::declareStableParticles( const std::vector<int>& pdgIds )
86 {
87 
88  for ( size_t i=0; i<pdgIds.size(); i++ )
89  {
90  // FIXME: need to double check if PID's are the same in Py6 & Py8,
91  // because the HepPDT translation tool is actually for **Py6**
92  //
93  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
94  //
95 // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
96  int PyID = pdgIds[i];
97  std::ostringstream pyCard ;
98  pyCard << PyID <<":mayDecay=false";
99 
100  if ( fMasterGen->particleData.isParticle( PyID ) ) {
101  fMasterGen->readString( pyCard.str() );
102  } else {
103 
104  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
105  << "recognize particle id = "
106  << PyID << std::endl;
107  }
108  // alternative:
109  // set the 2nd input argument warn=false
110  // - this way Py8 will NOT print warnings about unknown particle code(s)
111  // fMasterPtr->readString( pyCard.str(), false )
112  }
113 
114  return true;
115 
116 }
117 
118 bool Py8InterfaceBase:: declareSpecialSettings( const std::vector<std::string>& settings ){
119  for ( unsigned int iss=0; iss<settings.size(); iss++ ){
120  if ( settings[iss].find("QED-brem-off") != std::string::npos ){
121  fMasterGen->readString( "TimeShower:QEDshowerByL=off" );
122  }
123  else{
124  size_t fnd1 = settings[iss].find("Pythia8:");
125  if ( fnd1 != std::string::npos ){
126  std::string value = settings[iss].substr (fnd1+8);
127  fDecayer->readString(value);
128  }
129  }
130  }
131  return true;
132 }
133 
134 
136 {
137 
138  fMasterGen->stat();
139  return;
140 
141 }
142 
143 }
bool declareSpecialSettings(const std::vector< std::string > &)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
ParameterCollector fParameters
std::auto_ptr< Pythia8::Pythia > fMasterGen
HepMC::IO_AsciiParticles * ascii_io
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:7
std::vector< std::string > evtgenUserFiles
P8RndmEngine p8RndmEngine_
unsigned int pythiaPylistVerbosity
bool declareStableParticles(const std::vector< int > &)
unsigned int maxEventsToPrint
std::auto_ptr< Pythia8::Pythia > fDecayer
const_iterator end() const
const_iterator begin() const
volatile std::atomic< bool > shutdown_flag false