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 | Private 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
 
void p8SetRandomEngine (CLHEP::HepRandomEngine *v)
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
P8RndmEnginerandomEngine ()
 
bool readSettings (int)
 
virtual void statistics ()
 
 ~Py8InterfaceBase ()
 

Protected Attributes

HepMC::IO_AsciiParticles * ascii_io
 
std::auto_ptr< Pythia8::Pythia > fDecayer
 
std::auto_ptr< Pythia8::Pythia > fMasterGen
 
ParameterCollector fParameters
 
unsigned int maxEventsToPrint
 
bool pythiaHepMCVerbosity
 
bool pythiaHepMCVerbosityParticles
 
unsigned int pythiaPylistVerbosity
 
HepMC::Pythia8ToHepMC toHepMC
 

Private Attributes

P8RndmEngine p8RndmEngine_
 

Detailed Description

Definition at line 21 of file Py8InterfaceBase.h.

Constructor & Destructor Documentation

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

Definition at line 10 of file Py8InterfaceBase.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), Hadronizer_TuneCUETP8M1_13TeV_MLM_5f_max4j_LHE_pythia8_EvtGen_cff::maxEventsToPrint, dbtoconf::out, Hadronizer_TuneCUETP8M1_13TeV_MLM_5f_max4j_LHE_pythia8_EvtGen_cff::pythiaHepMCVerbosity, and Hadronizer_TuneCUETP8M1_13TeV_MLM_5f_max4j_LHE_pythia8_EvtGen_cff::pythiaPylistVerbosity.

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 }
ParameterCollector fParameters
std::auto_ptr< Pythia8::Pythia > fMasterGen
HepMC::IO_AsciiParticles * ascii_io
P8RndmEngine p8RndmEngine_
unsigned int pythiaPylistVerbosity
tuple out
Definition: dbtoconf.py:99
unsigned int maxEventsToPrint
std::auto_ptr< Pythia8::Pythia > fDecayer
gen::Py8InterfaceBase::~Py8InterfaceBase ( )
inline

Definition at line 26 of file Py8InterfaceBase.h.

26 {}

Member Function Documentation

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

Definition at line 29 of file Py8InterfaceBase.h.

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

Definition at line 93 of file Py8InterfaceBase.cc.

References spr::find(), AlCaHLTBitMon_QueryRunRegistry::string, and relativeConstraints::value.

93  {
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 }
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
std::auto_ptr< Pythia8::Pythia > fDecayer
bool gen::Py8InterfaceBase::declareStableParticles ( const std::vector< int > &  pdgIds)

Definition at line 60 of file Py8InterfaceBase.cc.

References i.

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 }
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.

void gen::Py8InterfaceBase::p8SetRandomEngine ( CLHEP::HepRandomEngine *  v)
inline

Definition at line 38 of file Py8InterfaceBase.h.

References p8RndmEngine_, and gen::P8RndmEngine::setRandomEngine().

Referenced by gen::Py8GunBase::setRandomEngine().

double v[5][pyjets_maxn]
void setRandomEngine(CLHEP::HepRandomEngine *v)
Definition: P8RndmEngine.h:35
P8RndmEngine p8RndmEngine_
P8RndmEngine& gen::Py8InterfaceBase::randomEngine ( )
inline
bool gen::Py8InterfaceBase::readSettings ( int  )

Definition at line 32 of file Py8InterfaceBase.cc.

References geometryCSVtoXML::line.

33 {
34 
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 }
ParameterCollector fParameters
std::auto_ptr< Pythia8::Pythia > fMasterGen
std::auto_ptr< Pythia8::Pythia > fDecayer
const_iterator end() const
const_iterator begin() const
void gen::Py8InterfaceBase::statistics ( )
virtual

Reimplemented in Pythia8Hadronizer, and gen::Py8GunBase.

Definition at line 110 of file Py8InterfaceBase.cc.

111 {
112 
113  fMasterGen->stat();
114  return;
115 
116 }
std::auto_ptr< Pythia8::Pythia > fMasterGen

Member Data Documentation

HepMC::IO_AsciiParticles* gen::Py8InterfaceBase::ascii_io
protected
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 46 of file Py8InterfaceBase.h.

Referenced by Pythia8Hadronizer::Pythia8Hadronizer().

unsigned int gen::Py8InterfaceBase::maxEventsToPrint
protected
P8RndmEngine gen::Py8InterfaceBase::p8RndmEngine_
private

Definition at line 56 of file Py8InterfaceBase.h.

Referenced by p8SetRandomEngine(), and randomEngine().

bool gen::Py8InterfaceBase::pythiaHepMCVerbosity
protected
bool gen::Py8InterfaceBase::pythiaHepMCVerbosityParticles
protected
unsigned int gen::Py8InterfaceBase::pythiaPylistVerbosity
protected
HepMC::Pythia8ToHepMC gen::Py8InterfaceBase::toHepMC
protected