CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
Pythia8Hadronizer Class Reference
Inheritance diagram for Pythia8Hadronizer:
gen::BaseHadronizer

Public Member Functions

const char * classname () const
 
bool decay ()
 
bool declareSpecialSettings (const std::vector< std::string >)
 
bool declareStableParticles (const std::vector< int > &pdgIds)
 
void finalizeEvent ()
 
bool generatePartonsAndHadronize ()
 
bool hadronize ()
 
bool initializeForExternalPartons ()
 
bool initializeForInternalPartons ()
 
 Pythia8Hadronizer (const edm::ParameterSet &params)
 
bool residualDecay ()
 
void statistics ()
 
 ~Pythia8Hadronizer ()
 
- Public Member Functions inherited from gen::BaseHadronizer
 BaseHadronizer (edm::ParameterSet const &ps)
 
edm::EventgetEDMEvent () const
 
HepMC::GenEvent * getGenEvent ()
 
GenEventInfoProductgetGenEventInfo ()
 
GenRunInfoProductgetGenRunInfo ()
 
const boost::shared_ptr
< lhef::LHERunInfo > & 
getLHERunInfo () const
 
void resetEvent (HepMC::GenEvent *event)
 
void resetEventInfo (GenEventInfoProduct *eventInfo)
 
void setEDMEvent (edm::Event &event)
 
void setLHEEvent (lhef::LHEEvent *event)
 
void setLHERunInfo (lhef::LHERunInfo *runInfo)
 
 ~BaseHadronizer ()
 

Private Attributes

double comEnergy
 Center-of-Mass energy. More...
 
std::auto_ptr< LHAupLesHoucheslhaUP
 
string LHEInputFileName
 
unsigned int maxEventsToPrint
 Events to print if verbosity. More...
 
ParameterCollector parameters
 
std::auto_ptr< Pythia > pythia
 
EventpythiaEvent
 
bool pythiaHepMCVerbosity
 HepMC verbosity flag. More...
 
unsigned int pythiaPylistVerbosity
 Pythia PYLIST Verbosity flag. More...
 
HepMC::I_Pythia8 toHepMC
 
UserHooks * UserHookPointer
 
bool useUserHook
 Switch User Hook flag. More...
 

Additional Inherited Members

- Protected Member Functions inherited from gen::BaseHadronizer
std::auto_ptr< HepMC::GenEvent > & event ()
 
std::auto_ptr
< GenEventInfoProduct > & 
eventInfo ()
 
lhef::LHEEventlheEvent ()
 
lhef::LHERunInfolheRunInfo ()
 
GenRunInfoProductrunInfo ()
 

Detailed Description

Definition at line 39 of file Pythia8Hadronizer.cc.

Constructor & Destructor Documentation

Pythia8Hadronizer::Pythia8Hadronizer ( const edm::ParameterSet params)

Definition at line 87 of file Pythia8Hadronizer.cc.

References edm::ParameterSet::exists(), gen::getEngineReference(), edm::ParameterSet::getParameter(), randomEngine, and useUserHook.

87  :
88  BaseHadronizer(params),
89  parameters(params.getParameter<edm::ParameterSet>("PythiaParameters")),
90  comEnergy(params.getParameter<double>("comEnergy")),
91  pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
92  pythiaHepMCVerbosity(params.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
93  maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
94  LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
95  useUserHook(false),
97 {
98  if( params.exists( "useUserHook" ) )
99  useUserHook = params.getParameter<bool>("useUserHook");
101 }
bool pythiaHepMCVerbosity
HepMC verbosity flag.
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double comEnergy
Center-of-Mass energy.
BaseHadronizer(edm::ParameterSet const &ps)
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool useUserHook
Switch User Hook flag.
CLHEP::HepRandomEngine & getEngineReference()
CLHEP::HepRandomEngine * randomEngine
Definition: PYR.cc:4
unsigned int pythiaPylistVerbosity
Pythia PYLIST Verbosity flag.
ParameterCollector parameters
unsigned int maxEventsToPrint
Events to print if verbosity.
UserHooks * UserHookPointer
Pythia8Hadronizer::~Pythia8Hadronizer ( )

Definition at line 103 of file Pythia8Hadronizer.cc.

104 {
105 }

Member Function Documentation

const char* Pythia8Hadronizer::classname ( ) const
inline

Definition at line 58 of file Pythia8Hadronizer.cc.

58 { return "Pythia8Hadronizer"; }
bool Pythia8Hadronizer::decay ( )

Definition at line 292 of file Pythia8Hadronizer.cc.

293 {
294  return true;
295 }
bool Pythia8Hadronizer::declareSpecialSettings ( const std::vector< std::string >  )

Definition at line 245 of file Pythia8Hadronizer.cc.

246 {
247  return true;
248 }
bool Pythia8Hadronizer::declareStableParticles ( const std::vector< int > &  pdgIds)

Definition at line 224 of file Pythia8Hadronizer.cc.

References i, and pythia.

225 {
226 
227  for ( size_t i=0; i<pdgIds.size(); i++ )
228  {
229  // FIXME: need to double check if PID's are the same in Py6 & Py8,
230  // because the HepPDT translation tool is actually for **Py6**
231  //
232  int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
233  std::ostringstream pyCard ;
234  pyCard << PyID <<":mayDecay=false";
235  pythia->readString( pyCard.str() );
236  // alternative:
237  // set the 2nd input argument warn=false
238  // - this way Py8 will NOT print warnings about unknown particle code(s)
239  // pythia->readString( pyCard.str(), false )
240  }
241 
242  return true;
243 
244 }
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< Pythia > pythia
void Pythia8Hadronizer::finalizeEvent ( )

Definition at line 302 of file Pythia8Hadronizer.cc.

References abs, gather_cfg::cout, gen::BaseHadronizer::event(), gen::BaseHadronizer::eventInfo(), gen::BaseHadronizer::lheEvent(), gen::BaseHadronizer::lheRunInfo(), maxEventsToPrint, pythia, pythiaHepMCVerbosity, pythiaPylistVerbosity, UserHookPointer, and useUserHook.

303 {
304  bool lhe = lheEvent() != 0;
305 
306  event()->set_signal_process_id(pythia->info.code());
307  event()->set_event_scale(pythia->info.pTHat()); //FIXME
308 
309  // Putting pdf info into the HepMC record
310  // There is the overloaded pythia8 HepMCInterface method fill_next_event
311  // that does this, but CMSSW GeneratorInterface does not fill HepMC
312  // record according to the HepMC convention (stores f(x) instead of x*f(x)
313  // and converts gluon PDG ID to zero). For this reason we use the
314  // method fill_next_event (above) that does NOT this and fill pdf info here
315  //
316  int id1 = pythia->info.id1();
317  int id2 = pythia->info.id2();
318  if (id1 == 21) id1 = 0;
319  if (id2 == 21) id2 = 0;
320  double x1 = pythia->info.x1();
321  double x2 = pythia->info.x2();
322  double Q = pythia->info.QRen();
323  double pdf1 = pythia->info.pdf1() / pythia->info.x1();
324  double pdf2 = pythia->info.pdf2() / pythia->info.x2();
325  event()->set_pdf_info(HepMC::PdfInfo(id1,id2,x1,x2,Q,pdf1,pdf2));
326 
327  // Storing weights. Will be moved to pythia8 HepMCInterface
328  //
329  if (lhe && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
330  // translate mb to pb (CMS/Gen "convention" as of May 2009)
331  event()->weights().push_back( pythia->info.weight() * 1.0e9 );
332  else
333  event()->weights().push_back( pythia->info.weight() );
334  //
335  // Below is the weight to be used in the case of reweighting by UserHook
336  // This is a temporary solution. It requires the explicit conversion of
337  // the pointer. In future versions (>150) the additional factor introduced
338  // in the UserHook will be stored by pythia. Also, putting the weights
339  // into the HepMC record will be made in the pythia8 HepMCInterface
340  //
341  if(useUserHook)
342  event()->weights().push_back(1./((PtHatReweightUserHook*)UserHookPointer)->getFactor());
343 
344  // now create the GenEventInfo product from the GenEvent and fill
345  // the missing pieces
346  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
347 
348  // in pythia pthat is used to subdivide samples into different bins
349  // in LHE mode the binning is done by the external ME generator
350  // which is likely not pthat, so only filling it for Py6 internal mode
351  if (!lhe) {
352  eventInfo()->setBinningValues(std::vector<double>(1, pythia->info.pTHat()));
353  }
354 
355  //******** Verbosity ********
356 
357  if (maxEventsToPrint > 0 &&
360  if (pythiaPylistVerbosity) {
361  pythia->info.list(std::cout);
362  pythia->event.list(std::cout);
363  }
364 
365  if (pythiaHepMCVerbosity) {
366  std::cout << "Event process = "
367  << pythia->info.code() << "\n"
368  << "----------------------" << std::endl;
369  event()->print();
370  }
371  }
372 }
bool pythiaHepMCVerbosity
HepMC verbosity flag.
std::auto_ptr< Pythia > pythia
#define abs(x)
Definition: mlp_lapack.h:159
bool useUserHook
Switch User Hook flag.
std::auto_ptr< HepMC::GenEvent > & event()
lhef::LHEEvent * lheEvent()
unsigned int pythiaPylistVerbosity
Pythia PYLIST Verbosity flag.
std::auto_ptr< GenEventInfoProduct > & eventInfo()
lhef::LHERunInfo * lheRunInfo()
unsigned int maxEventsToPrint
Events to print if verbosity.
tuple cout
Definition: gather_cfg.py:41
UserHooks * UserHookPointer
bool Pythia8Hadronizer::generatePartonsAndHadronize ( )

Definition at line 260 of file Pythia8Hadronizer.cc.

References gen::BaseHadronizer::event(), pythia, pythiaEvent, and toHepMC.

261 {
262  if (!pythia->next())
263  return false;
264 
265  event().reset(new HepMC::GenEvent);
266  toHepMC.fill_next_event(*pythiaEvent, event().get());
267 
268  return true;
269 }
std::auto_ptr< Pythia > pythia
std::auto_ptr< HepMC::GenEvent > & event()
HepMC::I_Pythia8 toHepMC
bool Pythia8Hadronizer::hadronize ( )

Definition at line 271 of file Pythia8Hadronizer.cc.

References lhef::LHEEvent::count(), gen::BaseHadronizer::event(), lhef::LHERunInfo::kAccepted, lhaUP, gen::BaseHadronizer::lheEvent(), LHEInputFileName, pythia, pythiaEvent, and toHepMC.

272 {
273  if(LHEInputFileName == string()) {
274  //cout << "start loading event" << endl;
275  lhaUP->loadEvent(lheEvent());
276  //cout << "finish loading event" << endl;
277  }
278 
279  if (!pythia->next())
280  return false;
281 
282  // update LHE matching statistics
283  //
285 
286  event().reset(new HepMC::GenEvent);
287  toHepMC.fill_next_event(*pythiaEvent, event().get());
288 
289  return true;
290 }
std::auto_ptr< Pythia > pythia
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:198
std::auto_ptr< HepMC::GenEvent > & event()
std::auto_ptr< LHAupLesHouches > lhaUP
lhef::LHEEvent * lheEvent()
HepMC::I_Pythia8 toHepMC
bool Pythia8Hadronizer::initializeForExternalPartons ( )

Definition at line 156 of file Pythia8Hadronizer.cc.

References gen::ParameterCollector::begin(), gather_cfg::cout, gen::ParameterCollector::end(), lhaUP, LHEInputFileName, gen::BaseHadronizer::lheRunInfo(), geometryCSVtoXML::line, parameters, pythia, pythiaEvent, and pythiaPylistVerbosity.

157 {
158 
159  std::cout << "Initializing for external partons" << std::endl;
160 
161  RandomP8* RP8 = new RandomP8();
162 
163  pythia.reset(new Pythia);
164 
165  pythia->setRndmEnginePtr(RP8);
166 
167  pythiaEvent = &pythia->event;
168 
170  line != parameters.end(); ++line) {
171  if (line->find("Random:") != std::string::npos)
172  throw cms::Exception("PythiaError")
173  << "Attempted to set random number "
174  "using Pythia commands. Please use "
175  "the RandomNumberGeneratorService."
176  << std::endl;
177 
178  if (!pythia->readString(*line))
179  throw cms::Exception("PythiaError")
180  << "Pythia 8 did not accept \""
181  << *line << "\"." << std::endl;
182  }
183 
184  if(pythiaPylistVerbosity > 10) {
186  pythia->settings.listAll();
188  pythia->particleData.listAll();
189  }
190 
191  if(LHEInputFileName != string()) {
192 
193  cout << endl;
194  cout << "LHE Input File Name = " << LHEInputFileName << endl;
195  cout << endl;
196  pythia->init(LHEInputFileName);
197 
198  } else {
199 
200  lhaUP.reset(new LHAupLesHouches());
201  lhaUP->loadRunInfo(lheRunInfo());
202  pythia->init(lhaUP.get());
203 
204  }
205 
206  return true;
207 }
std::auto_ptr< Pythia > pythia
std::auto_ptr< LHAupLesHouches > lhaUP
unsigned int pythiaPylistVerbosity
Pythia PYLIST Verbosity flag.
lhef::LHERunInfo * lheRunInfo()
ParameterCollector parameters
const_iterator end() const
const_iterator begin() const
tuple cout
Definition: gather_cfg.py:41
bool Pythia8Hadronizer::initializeForInternalPartons ( )

Definition at line 107 of file Pythia8Hadronizer.cc.

References gen::ParameterCollector::begin(), comEnergy, gen::ParameterCollector::end(), geometryCSVtoXML::line, parameters, pythia, pythiaEvent, pythiaPylistVerbosity, UserHookPointer, and useUserHook.

108 {
109  //Old code that used Pythia8 own random engine
110  //edm::Service<edm::RandomNumberGenerator> rng;
111  //uint32_t seed = rng->mySeed();
112  //Pythia8::Rndm::init(seed);
113 
114  RandomP8* RP8 = new RandomP8();
115 
116  pythia.reset(new Pythia);
117 
118  pythia->setRndmEnginePtr(RP8);
119  if(useUserHook){
121  pythia->setUserHooksPtr(UserHookPointer);
122  }
123 
124  pythiaEvent = &pythia->event;
125 
127  line != parameters.end(); ++line) {
128  if (line->find("Random:") != std::string::npos)
129  throw cms::Exception("PythiaError")
130  << "Attempted to set random number "
131  "using Pythia commands. Please use "
132  "the RandomNumberGeneratorService."
133  << std::endl;
134 
135  if (!pythia->readString(*line))
136  throw cms::Exception("PythiaError")
137  << "Pythia 8 did not accept \""
138  << *line << "\"." << std::endl;
139  }
140 
141  if(pythiaPylistVerbosity > 10) {
143  pythia->settings.listAll();
145  pythia->particleData.listAll();
146  }
147 
148  pythia->init(2212, 2212, comEnergy);
149 
150  pythia->settings.listChanged();
151 
152  return true;
153 }
double comEnergy
Center-of-Mass energy.
std::auto_ptr< Pythia > pythia
bool useUserHook
Switch User Hook flag.
unsigned int pythiaPylistVerbosity
Pythia PYLIST Verbosity flag.
ParameterCollector parameters
const_iterator end() const
const_iterator begin() const
UserHooks * UserHookPointer
bool Pythia8Hadronizer::residualDecay ( )

Definition at line 297 of file Pythia8Hadronizer.cc.

298 {
299  return true;
300 }
void Pythia8Hadronizer::statistics ( )

Definition at line 251 of file Pythia8Hadronizer.cc.

References pythia, gen::BaseHadronizer::runInfo(), and GenRunInfoProduct::setInternalXSec().

252 {
253  pythia->statistics();
254 
255  double xsec = pythia->info.sigmaGen(); // cross section in mb
256  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
257  runInfo().setInternalXSec(xsec);
258 }
std::auto_ptr< Pythia > pythia
void setInternalXSec(const XSec &xsec)
GenRunInfoProduct & runInfo()

Member Data Documentation

double Pythia8Hadronizer::comEnergy
private

Center-of-Mass energy.

Definition at line 64 of file Pythia8Hadronizer.cc.

Referenced by initializeForInternalPartons().

std::auto_ptr<LHAupLesHouches> Pythia8Hadronizer::lhaUP
private

Definition at line 77 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

string Pythia8Hadronizer::LHEInputFileName
private

Definition at line 72 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

unsigned int Pythia8Hadronizer::maxEventsToPrint
private

Events to print if verbosity.

Definition at line 70 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent().

ParameterCollector Pythia8Hadronizer::parameters
private
std::auto_ptr<Pythia> Pythia8Hadronizer::pythia
private
Event* Pythia8Hadronizer::pythiaEvent
private
bool Pythia8Hadronizer::pythiaHepMCVerbosity
private

HepMC verbosity flag.

Definition at line 68 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent().

unsigned int Pythia8Hadronizer::pythiaPylistVerbosity
private

Pythia PYLIST Verbosity flag.

Definition at line 66 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), initializeForExternalPartons(), and initializeForInternalPartons().

HepMC::I_Pythia8 Pythia8Hadronizer::toHepMC
private

Definition at line 81 of file Pythia8Hadronizer.cc.

Referenced by generatePartonsAndHadronize(), and hadronize().

UserHooks* Pythia8Hadronizer::UserHookPointer
private

Definition at line 83 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), and initializeForInternalPartons().

bool Pythia8Hadronizer::useUserHook
private

Switch User Hook flag.

Definition at line 75 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), initializeForInternalPartons(), and Pythia8Hadronizer().