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  }
51 
52 }
53 
55 {
56 
58  line != fParameters.end(); ++line )
59  {
60  if (line->find("Random:") != std::string::npos)
61  throw cms::Exception("PythiaError") << "Attempted to set random number "
62  "using Pythia commands. Please use " "the RandomNumberGeneratorService."
63  << std::endl;
64 
65  if (!fMasterGen->readString(*line)) throw cms::Exception("PythiaError")
66  << "Pythia 8 did not accept \""
67  << *line << "\"." << std::endl;
68 
69  if (line->find("ParticleDecays:") != std::string::npos) {
70 
71  if (!fDecayer->readString(*line)) throw cms::Exception("PythiaError")
72  << "Pythia 8 Decayer did not accept \""
73  << *line << "\"." << std::endl;
74  }
75 
76  }
77 
78  return true;
79 
80 }
81 
82 bool Py8InterfaceBase::declareStableParticles( const std::vector<int>& pdgIds )
83 {
84 
85  for ( size_t i=0; i<pdgIds.size(); i++ )
86  {
87  // FIXME: need to double check if PID's are the same in Py6 & Py8,
88  // because the HepPDT translation tool is actually for **Py6**
89  //
90  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
91  //
92 // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
93  int PyID = pdgIds[i];
94  std::ostringstream pyCard ;
95  pyCard << PyID <<":mayDecay=false";
96 
97  if ( fMasterGen->particleData.isParticle( PyID ) ) {
98  fMasterGen->readString( pyCard.str() );
99  } else {
100 
101  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
102  << "recognize particle id = "
103  << PyID << std::endl;
104  }
105  // alternative:
106  // set the 2nd input argument warn=false
107  // - this way Py8 will NOT print warnings about unknown particle code(s)
108  // fMasterPtr->readString( pyCard.str(), false )
109  }
110 
111  return true;
112 
113 }
114 
115 bool Py8InterfaceBase:: declareSpecialSettings( const std::vector<std::string>& settings ){
116  for ( unsigned int iss=0; iss<settings.size(); iss++ ){
117  if ( settings[iss].find("QED-brem-off") != std::string::npos ){
118  fMasterGen->readString( "TimeShower:QEDshowerByL=off" );
119  }
120  else{
121  size_t fnd1 = settings[iss].find("Pythia8:");
122  if ( fnd1 != std::string::npos ){
123  std::string value = settings[iss].substr (fnd1+8);
124  fDecayer->readString(value);
125  }
126  }
127  }
128  return true;
129 }
130 
131 
133 {
134 
135  fMasterGen->stat();
136  return;
137 
138 }
139 
140 }
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
P8RndmEngine p8RndmEngine_
unsigned int pythiaPylistVerbosity
bool declareStableParticles(const std::vector< int > &)
tuple out
Definition: dbtoconf.py:99
unsigned int maxEventsToPrint
std::auto_ptr< Pythia8::Pythia > fDecayer
const_iterator end() const
const_iterator begin() const
volatile std::atomic< bool > shutdown_flag false