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 using namespace Pythia8;
7 
8 namespace gen {
9 
10 Py8InterfaceBase::Py8InterfaceBase( edm::ParameterSet const& ps )
11 {
12  fMasterGen.reset(new Pythia);
13  fDecayer.reset(new Pythia);
14 
15  fMasterGen->readString("Next:numberShowEvent = 0");
16  fDecayer->readString("Next:numberShowEvent = 0");
17 
18  fMasterGen->setRndmEnginePtr( &p8RndmEngine_ );
19  fDecayer->setRndmEnginePtr( &p8RndmEngine_ );
20 
21  fParameters = ps.getParameter<edm::ParameterSet>("PythiaParameters");
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 
28  if(pythiaHepMCVerbosityParticles)
29  ascii_io = new HepMC::IO_AsciiParticles("cout", std::ios::out);
30 }
31 
32 bool Py8InterfaceBase::readSettings( int )
33 {
34 
35  for ( ParameterCollector::const_iterator line = fParameters.begin();
36  line != fParameters.end(); ++line )
37  {
38  if (line->find("Random:") != std::string::npos)
39  throw cms::Exception("PythiaError") << "Attempted to set random number "
40  "using Pythia commands. Please use " "the RandomNumberGeneratorService."
41  << std::endl;
42 
43  if (!fMasterGen->readString(*line)) throw cms::Exception("PythiaError")
44  << "Pythia 8 did not accept \""
45  << *line << "\"." << std::endl;
46 
47  if (line->find("ParticleDecays:") != std::string::npos) {
48 
49  if (!fDecayer->readString(*line)) throw cms::Exception("PythiaError")
50  << "Pythia 8 Decayer did not accept \""
51  << *line << "\"." << std::endl;
52  }
53 
54  }
55 
56  return true;
57 
58 }
59 
60 bool Py8InterfaceBase::declareStableParticles( const std::vector<int>& pdgIds )
61 {
62 
63  for ( size_t i=0; i<pdgIds.size(); i++ )
64  {
65  // FIXME: need to double check if PID's are the same in Py6 & Py8,
66  // because the HepPDT translation tool is actually for **Py6**
67  //
68  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
69  //
70 // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
71  int PyID = pdgIds[i];
72  std::ostringstream pyCard ;
73  pyCard << PyID <<":mayDecay=false";
74 
75  if ( fMasterGen->particleData.isParticle( PyID ) ) {
76  fMasterGen->readString( pyCard.str() );
77  } else {
78 
79  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
80  << "recognize particle id = "
81  << PyID << std::endl;
82  }
83  // alternative:
84  // set the 2nd input argument warn=false
85  // - this way Py8 will NOT print warnings about unknown particle code(s)
86  // fMasterPtr->readString( pyCard.str(), false )
87  }
88 
89  return true;
90 
91 }
92 
93 bool Py8InterfaceBase:: declareSpecialSettings( const std::vector<std::string>& settings ){
94  for ( unsigned int iss=0; iss<settings.size(); iss++ ){
95  if ( settings[iss].find("QED-brem-off") != std::string::npos ){
96  fMasterGen->readString( "TimeShower:QEDshowerByL=off" );
97  }
98  else{
99  size_t fnd1 = settings[iss].find("Pythia8:");
100  if ( fnd1 != std::string::npos ){
101  std::string value = settings[iss].substr (fnd1+8);
102  fDecayer->readString(value);
103  }
104  }
105  }
106  return true;
107 }
108 
109 
111 {
112 
113  fMasterGen->stat();
114  return;
115 
116 }
117 
118 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
tuple out
Definition: dbtoconf.py:99