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  maxEventsToPrint = ps.getUntrackedParameter<int>("maxEventsToPrint", 0);
26 
27 }
28 
29 bool Py8InterfaceBase::readSettings( int )
30 {
31 
32  for ( ParameterCollector::const_iterator line = fParameters.begin();
33  line != fParameters.end(); ++line )
34  {
35  if (line->find("Random:") != std::string::npos)
36  throw cms::Exception("PythiaError") << "Attempted to set random number "
37  "using Pythia commands. Please use " "the RandomNumberGeneratorService."
38  << std::endl;
39 
40  if (!fMasterGen->readString(*line)) throw cms::Exception("PythiaError")
41  << "Pythia 8 did not accept \""
42  << *line << "\"." << std::endl;
43  }
44 
45  return true;
46 
47 }
48 
49 bool Py8InterfaceBase::declareStableParticles( const std::vector<int>& pdgIds )
50 {
51 
52  for ( size_t i=0; i<pdgIds.size(); i++ )
53  {
54  // FIXME: need to double check if PID's are the same in Py6 & Py8,
55  // because the HepPDT translation tool is actually for **Py6**
56  //
57  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
58  //
59 // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
60  int PyID = pdgIds[i];
61  std::ostringstream pyCard ;
62  pyCard << PyID <<":mayDecay=false";
63 
64  if ( fMasterGen->particleData.isParticle( PyID ) ) {
65  fMasterGen->readString( pyCard.str() );
66  } else {
67 
68  edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
69  << "recognize particle id = "
70  << PyID << std::endl;
71  }
72  // alternative:
73  // set the 2nd input argument warn=false
74  // - this way Py8 will NOT print warnings about unknown particle code(s)
75  // fMasterPtr->readString( pyCard.str(), false )
76  }
77 
78  return true;
79 
80 }
81 
82 bool Py8InterfaceBase:: declareSpecialSettings( const std::vector<std::string>& settings )
83 {
84 
85  for ( unsigned int iss=0; iss<settings.size(); iss++ )
86  {
87  if ( settings[iss].find("QED-brem-off") == std::string::npos ) continue;
88  fMasterGen->readString( "TimeShower:QEDshowerByL=off" );
89  }
90 
91  return true;
92 }
93 
94 
96 {
97 
98  fMasterGen->statistics();
99  return;
100 
101 }
102 
103 }
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 pythiaPylistVerbosity
tuple maxEventsToPrint
tuple pythiaHepMCVerbosity