CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes
gen::Py8InterfaceBase Class Referenceabstract

#include <Py8InterfaceBase.h>

Inheritance diagram for gen::Py8InterfaceBase:
gen::Py8GunBase Pythia8Hadronizer gen::Py8EGun gen::Py8JetGun gen::Py8PtGun

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
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
bool readSettings (int)
 
virtual void statistics ()
 
 ~Py8InterfaceBase ()
 

Protected Attributes

std::auto_ptr< Pythia8::Pythia > fDecayer
 
std::auto_ptr< Pythia8::Pythia > fMasterGen
 
ParameterCollector fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
unsigned int pythiaPylistVerbosity
 
HepMC::I_Pythia8 toHepMC
 

Detailed Description

Definition at line 14 of file Py8InterfaceBase.h.

Constructor & Destructor Documentation

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

Definition at line 11 of file Py8InterfaceBase.cc.

References gen::getEngineReference(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), step1_ZMM_7Tev::maxEventsToPrint, step1_ZMM_7Tev::pythiaHepMCVerbosity, step1_ZMM_7Tev::pythiaPylistVerbosity, and randomEngine.

12 {
13 
15 
16  fMasterGen.reset(new Pythia);
17  fDecayer.reset(new Pythia);
18 
19  fMasterGen->readString("Next:numberShowEvent = 0");
20  fDecayer->readString("Next:numberShowEvent = 0");
21 
22  // RandomP8* RP8 = new RandomP8();
23  fMasterGen->setRndmEnginePtr( new RandomP8() );
24  fDecayer->setRndmEnginePtr( new RandomP8() );
25 
26  fParameters = ps.getParameter<edm::ParameterSet>("PythiaParameters");
27 
28  pythiaPylistVerbosity = ps.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
29  pythiaHepMCVerbosity = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false);
30  maxEventsToPrint = ps.getUntrackedParameter<int>("maxEventsToPrint", 0);
31 
32 }
ParameterCollector fParameters
std::auto_ptr< Pythia8::Pythia > fMasterGen
CLHEP::HepRandomEngine & getEngineReference()
unsigned int pythiaPylistVerbosity
CLHEP::HepRandomEngine * randomEngine
Definition: RandomP8.cc:5
unsigned int maxEventsToPrint
std::auto_ptr< Pythia8::Pythia > fDecayer
gen::Py8InterfaceBase::~Py8InterfaceBase ( )
inline

Definition at line 19 of file Py8InterfaceBase.h.

19 {}

Member Function Documentation

virtual const char* gen::Py8InterfaceBase::classname ( ) const
pure virtual
bool gen::Py8InterfaceBase::decay ( )
inline

Definition at line 22 of file Py8InterfaceBase.h.

22 { return true; } // NOT used - let's call it "design imperfection"
bool gen::Py8InterfaceBase::declareSpecialSettings ( const std::vector< std::string > &  settings)

Definition at line 79 of file Py8InterfaceBase.cc.

References spr::find().

80 {
81 
82  for ( unsigned int iss=0; iss<settings.size(); iss++ )
83  {
84  if ( settings[iss].find("QED-brem-off") == std::string::npos ) continue;
85  fMasterGen->readString( "TimeShower:QEDshowerByL=off" );
86  }
87 
88  return true;
89 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
bool gen::Py8InterfaceBase::declareStableParticles ( const std::vector< int > &  pdgIds)

Definition at line 54 of file Py8InterfaceBase.cc.

References i.

55 {
56 
57  for ( size_t i=0; i<pdgIds.size(); i++ )
58  {
59  // FIXME: need to double check if PID's are the same in Py6 & Py8,
60  // because the HepPDT translation tool is actually for **Py6**
61  //
62  // well, actually it looks like Py8 operates in PDT id's rather than Py6's
63  //
64  // int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
65  int PyID = pdgIds[i];
66  std::ostringstream pyCard ;
67  pyCard << PyID <<":mayDecay=false";
68  fMasterGen->readString( pyCard.str() );
69  // alternative:
70  // set the 2nd input argument warn=false
71  // - this way Py8 will NOT print warnings about unknown particle code(s)
72  // fMasterPtr->readString( pyCard.str(), false )
73  }
74 
75  return true;
76 
77 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia8::Pythia > fMasterGen
virtual void gen::Py8InterfaceBase::finalizeEvent ( )
pure virtual

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

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

Implemented in Pythia8Hadronizer, and gen::Py8GunBase.

bool gen::Py8InterfaceBase::readSettings ( int  )

Definition at line 34 of file Py8InterfaceBase.cc.

References geometryCSVtoXML::line.

35 {
36 
38  line != fParameters.end(); ++line )
39  {
40  if (line->find("Random:") != std::string::npos)
41  throw cms::Exception("PythiaError") << "Attempted to set random number "
42  "using Pythia commands. Please use " "the RandomNumberGeneratorService."
43  << std::endl;
44 
45  if (!fMasterGen->readString(*line)) throw cms::Exception("PythiaError")
46  << "Pythia 8 did not accept \""
47  << *line << "\"." << std::endl;
48  }
49 
50  return true;
51 
52 }
ParameterCollector fParameters
std::auto_ptr< Pythia8::Pythia > fMasterGen
const_iterator end() const
const_iterator begin() const
void gen::Py8InterfaceBase::statistics ( )
virtual

Reimplemented in Pythia8Hadronizer, and gen::Py8GunBase.

Definition at line 91 of file Py8InterfaceBase.cc.

92 {
93 
94  fMasterGen->statistics();
95  return;
96 
97 }
std::auto_ptr< Pythia8::Pythia > fMasterGen

Member Data Documentation

std::auto_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fDecayer
protected
std::auto_ptr<Pythia8::Pythia> gen::Py8InterfaceBase::fMasterGen
protected
ParameterCollector gen::Py8InterfaceBase::fParameters
protected

Definition at line 36 of file Py8InterfaceBase.h.

Referenced by Pythia8Hadronizer::Pythia8Hadronizer().

unsigned int gen::Py8InterfaceBase::maxEventsToPrint
protected
bool gen::Py8InterfaceBase::pythiaHepMCVerbosity
protected
unsigned int gen::Py8InterfaceBase::pythiaPylistVerbosity
protected
HepMC::I_Pythia8 gen::Py8InterfaceBase::toHepMC
protected