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 Types | Private Member Functions | Private Attributes | Static Private Attributes
Pythia8Hadronizer Class Reference
Inheritance diagram for Pythia8Hadronizer:
gen::BaseHadronizer gen::Py8InterfaceBase

Public Member Functions

const char * classname () const override
 
void finalizeEvent () override
 
bool generatePartonsAndHadronize () override
 
bool hadronize ()
 
bool initializeForExternalPartons ()
 
bool initializeForInternalPartons () override
 
 Pythia8Hadronizer (const edm::ParameterSet &params)
 
virtual bool residualDecay ()
 
void statistics () override
 
 ~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)
 
virtual bool select (HepMC::GenEvent *) const
 
void setEDMEvent (edm::Event &event)
 
void setLHEEvent (lhef::LHEEvent *event)
 
void setLHERunInfo (lhef::LHERunInfo *runInfo)
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
std::vector< std::string > const & sharedResources () const
 
 ~BaseHadronizer ()
 
- Public Member Functions inherited from gen::Py8InterfaceBase
bool decay ()
 
bool declareSpecialSettings (const std::vector< std::string > &)
 
bool declareStableParticles (const std::vector< int > &)
 
void p8SetRandomEngine (CLHEP::HepRandomEngine *v)
 
 Py8InterfaceBase (edm::ParameterSet const &ps)
 
P8RndmEnginerandomEngine ()
 
bool readSettings (int)
 
 ~Py8InterfaceBase ()
 

Private Types

enum  { PP, PPbar, ElectronPositron }
 

Private Member Functions

virtual void doSetRandomEngine (CLHEP::HepRandomEngine *v) override
 
virtual std::vector
< std::string > const & 
doSharedResources () const override
 

Private Attributes

double comEnergy
 Center-of-Mass energy. More...
 
vector< float > DJR
 
int EV1_emittedMode
 
int EV1_maxVetoCount
 
bool EV1_MPIvetoOn
 
int EV1_nFinal
 
int EV1_pTdefMode
 
int EV1_pTempMode
 
int EV1_pThardMode
 
bool EV1_vetoOn
 
double fBeam1PZ
 
double fBeam2PZ
 
PowhegHooks * fEmissionVetoHook
 
EmissionVetoHook1fEmissionVetoHook1
 
int fInitialState
 
JetMatchingHookfJetMatchingHook
 
Pythia8::JetMatchingMadgraph * fJetMatchingPy8InternalHook
 
Pythia8::amcnlo_unitarised_interface * fMergingHook
 
UserHooks * fReweightPtHatRapUserHook
 
UserHooks * fReweightRapUserHook
 
UserHooks * fReweightUserHook
 
std::auto_ptr< LHAupLesHoucheslhaUP
 
string LHEInputFileName
 
int nFSRveto
 
int NHooks
 
int nISRveto
 
int nME
 
int nMEFiltered
 
std::string slhafile_
 

Static Private Attributes

static const std::vector
< std::string > 
p8SharedResources = { edm::SharedResourceNames::kPythia8 }
 

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 ()
 
- Protected Attributes inherited from gen::Py8InterfaceBase
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
 

Detailed Description

Definition at line 55 of file Pythia8Hadronizer.cc.

Member Enumeration Documentation

anonymous enum
private

Constructor & Destructor Documentation

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

Definition at line 136 of file Pythia8Hadronizer.cc.

References gen::ParameterCollector::begin(), edm::errors::Configuration, ElectronPositron, gen::ParameterCollector::end(), EV1_emittedMode, EV1_maxVetoCount, EV1_MPIvetoOn, EV1_nFinal, EV1_pTdefMode, EV1_pTempMode, EV1_pThardMode, EV1_vetoOn, Exception, edm::ParameterSet::exists(), python.connectstrParser::f1, ztee::fd, fEmissionVetoHook1, fInitialState, fJetMatchingHook, gen::Py8InterfaceBase::fMasterGen, gen::Py8InterfaceBase::fParameters, fReweightPtHatRapUserHook, fReweightRapUserHook, fReweightUserHook, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), geometryCSVtoXML::line, NHooks, PP, PPbar, slhafile_, AlCaHLTBitMon_QueryRunRegistry::string, and TablePrint::write.

136  :
137  BaseHadronizer(params), Py8InterfaceBase(params),
138  comEnergy(params.getParameter<double>("comEnergy")),
139  LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
140  fInitialState(PP),
144  NHooks(0)
145 {
146 
147  // J.Y.: the following 3 parameters are hacked "for a reason"
148  //
149  if ( params.exists( "PPbarInitialState" ) )
150  {
151  if ( fInitialState == PP )
152  {
154  edm::LogImportant("GeneratorInterface|Pythia8Interface")
155  << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
156  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
157  }
158  else
159  {
160  // probably need to throw on attempt to override ?
161  }
162  }
163  else if ( params.exists( "ElectronPositronInitialState" ) )
164  {
165  if ( fInitialState == PP )
166  {
168  edm::LogInfo("GeneratorInterface|Pythia8Interface")
169  << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
170  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
171  }
172  else
173  {
174  // probably need to throw on attempt to override ?
175  }
176  }
177  else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
178  {
179  // throw on unknown initial state !
180  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
181  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
182  }
183 
184  if( params.exists( "SLHAFileForPythia8" ) ) {
185  std::string slhafilenameshort = params.getParameter<string>("SLHAFileForPythia8");
186  edm::FileInPath f1( slhafilenameshort );
187 
188  fMasterGen->settings.mode("SLHA:readFrom", 2);
189  fMasterGen->settings.word("SLHA:file", f1.fullPath());
190 
192  line != fParameters.end(); ++line ) {
193  if (line->find("SLHA:file") != std::string::npos)
194  throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
195  << "using Pythia8 card SLHA:file and Pythia8Interface card SLHAFileForPythia8"
196  << std::endl;
197  }
198  }
199  else if( params.exists( "SLHATableForPythia8" ) ) {
200  std::string slhatable = params.getParameter<string>("SLHATableForPythia8");
201 
202  char tempslhaname[] = "pythia8SLHAtableXXXXXX";
203  int fd = mkstemp(tempslhaname);
204  write(fd,slhatable.c_str(),slhatable.size());
205  close(fd);
206 
207  slhafile_ = tempslhaname;
208 
209  fMasterGen->settings.mode("SLHA:readFrom", 2);
210  fMasterGen->settings.word("SLHA:file", slhafile_);
211 
213  line != fParameters.end(); ++line ) {
214  if (line->find("SLHA:file") != std::string::npos)
215  throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
216  << "using Pythia8 card SLHA:file and Pythia8Interface card SLHATableForPythia8"
217  << std::endl;
218  }
219  }
220 
221  // Reweight user hook
222  //
223  if( params.exists( "reweightGen" ) )
225  if( params.exists( "reweightGenRap" ) )
226  {
227  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
228  edm::ParameterSet rgrParams =
229  params.getParameter<edm::ParameterSet>("reweightGenRap");
231  new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
232  rgrParams.getParameter<double>("yLabPower"),
233  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
234  rgrParams.getParameter<double>("yCMPower"),
235  rgrParams.getParameter<double>("pTHatMin"),
236  rgrParams.getParameter<double>("pTHatMax"));
237  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
238  }
239  if( params.exists( "reweightGenPtHatRap" ) )
240  {
241  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
242  edm::ParameterSet rgrParams =
243  params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
245  new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
246  rgrParams.getParameter<double>("yLabPower"),
247  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
248  rgrParams.getParameter<double>("yCMPower"),
249  rgrParams.getParameter<double>("pTHatMin"),
250  rgrParams.getParameter<double>("pTHatMax"));
251  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
252  }
253 
254  if( params.exists( "useUserHook" ) )
255  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
256  <<" Obsolete parameter: useUserHook \n Please use the actual one instead \n";
257 
258  // PS matching prototype
259  //
260  if ( params.exists("jetMatching") )
261  {
262  edm::ParameterSet jmParams =
263  params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
264  std::string scheme = jmParams.getParameter<std::string>("scheme");
265  if ( scheme == "Madgraph" || scheme == "MadgraphFastJet" )
266  {
267  fJetMatchingHook = new JetMatchingHook( jmParams, &fMasterGen->info );
268  }
269  }
270 
271  // Pythia8Interface emission veto
272  //
273  if ( params.exists("emissionVeto1") )
274  {
275  EV1_nFinal = -1;
276  if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
277  EV1_vetoOn = true;
278  if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
279  EV1_maxVetoCount = 10;
280  if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
281  EV1_pThardMode = 1;
282  if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
283  EV1_pTempMode = 0;
284  if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
285  if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
286  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
287  <<" Wrong value for EV1_pTempMode code\n";
288  EV1_emittedMode = 0;
289  if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
290  EV1_pTdefMode = 1;
291  if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
292  EV1_MPIvetoOn = false;
293  if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
297  }
298 
302  if(fJetMatchingHook) NHooks++;
304  if(NHooks > 1)
305  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
306  <<" Too many User Hooks. \n Please choose one from: reweightGen, reweightGenRap, reweightGenPtHatRap, jetMatching, emissionVeto1 \n";
307  if(fReweightUserHook) fMasterGen->setUserHooksPtr(fReweightUserHook);
310  if(fJetMatchingHook) fMasterGen->setUserHooksPtr(fJetMatchingHook);
311  if(fEmissionVetoHook1) {
312  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
313  fMasterGen->setUserHooksPtr(fEmissionVetoHook1);
314  }
315 
316 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
ParameterCollector fParameters
double comEnergy
Center-of-Mass energy.
std::auto_ptr< Pythia8::Pythia > fMasterGen
EmissionVetoHook1 * fEmissionVetoHook1
UserHooks * fReweightUserHook
BaseHadronizer(edm::ParameterSet const &ps)
bool exists(std::string const &parameterName) const
checks if a parameter exists
Pythia8::JetMatchingMadgraph * fJetMatchingPy8InternalHook
Py8InterfaceBase(edm::ParameterSet const &ps)
PowhegHooks * fEmissionVetoHook
UserHooks * fReweightPtHatRapUserHook
tuple fd
Definition: ztee.py:136
UserHooks * fReweightRapUserHook
Pythia8::amcnlo_unitarised_interface * fMergingHook
JetMatchingHook * fJetMatchingHook
const_iterator end() const
const_iterator begin() const
Pythia8Hadronizer::~Pythia8Hadronizer ( )

Definition at line 319 of file Pythia8Hadronizer.cc.

References fEmissionVetoHook, fEmissionVetoHook1, python.multivaluedict::remove(), and slhafile_.

320 {
321 // do we need to delete UserHooks/JetMatchingHook here ???
324 
325  //clean up temp file
326  if (!slhafile_.empty()) {
327  std::remove(slhafile_.c_str());
328  }
329 
330 }
EmissionVetoHook1 * fEmissionVetoHook1
PowhegHooks * fEmissionVetoHook

Member Function Documentation

const char* Pythia8Hadronizer::classname ( ) const
inlineoverridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 74 of file Pythia8Hadronizer.cc.

74 { return "Pythia8Hadronizer"; }
virtual void Pythia8Hadronizer::doSetRandomEngine ( CLHEP::HepRandomEngine *  v)
inlineoverrideprivatevirtual

Reimplemented from gen::BaseHadronizer.

Definition at line 78 of file Pythia8Hadronizer.cc.

78 { p8SetRandomEngine(v); }
void p8SetRandomEngine(CLHEP::HepRandomEngine *v)
virtual std::vector<std::string> const& Pythia8Hadronizer::doSharedResources ( ) const
inlineoverrideprivatevirtual

Reimplemented from gen::BaseHadronizer.

Definition at line 79 of file Pythia8Hadronizer.cc.

79 { return p8SharedResources; }
static const std::vector< std::string > p8SharedResources
void Pythia8Hadronizer::finalizeEvent ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 620 of file Pythia8Hadronizer.cc.

References gen::Py8InterfaceBase::ascii_io, gather_cfg::cout, DJR, gen::BaseHadronizer::event(), gen::BaseHadronizer::eventInfo(), gen::Py8InterfaceBase::fMasterGen, gen::BaseHadronizer::lheEvent(), gen::Py8InterfaceBase::maxEventsToPrint, nME, nMEFiltered, gen::Py8InterfaceBase::pythiaHepMCVerbosity, gen::Py8InterfaceBase::pythiaHepMCVerbosityParticles, and gen::Py8InterfaceBase::pythiaPylistVerbosity.

621 {
622  bool lhe = lheEvent() != 0;
623 
624  // now create the GenEventInfo product from the GenEvent and fill
625  // the missing pieces
626  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
627 
628  // in pythia pthat is used to subdivide samples into different bins
629  // in LHE mode the binning is done by the external ME generator
630  // which is likely not pthat, so only filling it for Py6 internal mode
631  if (!lhe) {
632  eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
633  }
634 
635  eventInfo()->setDJR(DJR);
636  eventInfo()->setNMEPartons(nME);
637  eventInfo()->setNMEPartonsFiltered(nMEFiltered);
638 
639  //******** Verbosity ********
640 
641  if (maxEventsToPrint > 0 &&
645  if (pythiaPylistVerbosity) {
646  fMasterGen->info.list(std::cout);
647  fMasterGen->event.list(std::cout);
648  }
649 
650  if (pythiaHepMCVerbosity) {
651  std::cout << "Event process = "
652  << fMasterGen->info.code() << "\n"
653  << "----------------------" << std::endl;
654  event()->print();
655  }
657  std::cout << "Event process = "
658  << fMasterGen->info.code() << "\n"
659  << "----------------------" << std::endl;
660  ascii_io->write_event(event().get());
661  }
662  }
663 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
HepMC::IO_AsciiParticles * ascii_io
std::auto_ptr< HepMC::GenEvent > & event()
lhef::LHEEvent * lheEvent()
unsigned int pythiaPylistVerbosity
std::auto_ptr< GenEventInfoProduct > & eventInfo()
unsigned int maxEventsToPrint
vector< float > DJR
tuple cout
Definition: gather_cfg.py:121
bool Pythia8Hadronizer::generatePartonsAndHadronize ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 493 of file Pythia8Hadronizer.cc.

References gen::BaseHadronizer::event(), gen::Py8InterfaceBase::fMasterGen, and gen::Py8InterfaceBase::toHepMC.

494 {
495 
496  if (!fMasterGen->next()) return false;
497 
498  event().reset(new HepMC::GenEvent);
499  return toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
500 
501 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
std::auto_ptr< HepMC::GenEvent > & event()
HepMC::Pythia8ToHepMC toHepMC
bool Pythia8Hadronizer::hadronize ( )

Definition at line 504 of file Pythia8Hadronizer.cc.

References funct::abs(), JetMatchingHook::beforeHadronization(), lhef::LHEEvent::count(), DJR, gen::BaseHadronizer::event(), fEmissionVetoHook, fJetMatchingHook, fJetMatchingPy8InternalHook, gen::Py8InterfaceBase::fMasterGen, fMergingHook, lhef::LHERunInfo::kAccepted, lhef::LHERunInfo::kSelected, lhaUP, gen::BaseHadronizer::lheEvent(), LHEInputFileName, min(), nFSRveto, nISRveto, nME, nMEFiltered, JetMatchingHook::resetMatchingStatus(), and gen::Py8InterfaceBase::toHepMC.

505 {
506  DJR.resize(0);
507  nME = -1;
508  nMEFiltered = -1;
509  if(LHEInputFileName == string()) lhaUP->loadEvent(lheEvent());
510 
511  if ( fJetMatchingHook )
512  {
515  }
516 
517  bool py8next = fMasterGen->next();
518 
519  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
520  if (fMergingHook) {
521  mergeweight *= fMergingHook->getNormFactor();
522  }
523 
524 
525  //protect against 0-weight from ckkw or similar
526  if (!py8next || std::abs(mergeweight)==0.)
527  {
528  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
529  event().reset();
530  return false;
531  }
532 
534  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->GetDJR();
535  //cap size of djr vector to save storage space (keep only up to first 6 elements)
536  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
537  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
538  DJR.push_back(djrmatch[idjr]);
539  }
540 
541  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
542  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
543  }
544 
545  // update LHE matching statistics
546  //
547  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
548 
549  event().reset(new HepMC::GenEvent);
550  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
551  if (!py8hepmc) {
552  return false;
553  }
554 
555  //add ckkw/umeps/unlops merging weight
556  if (mergeweight!=1.) {
557  event()->weights()[0] *= mergeweight;
558  }
559 
560  if (fEmissionVetoHook) {
561  nISRveto += fEmissionVetoHook->getNISRveto();
562  nFSRveto += fEmissionVetoHook->getNFSRveto();
563  }
564 
565  return true;
566 
567 
568 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
Pythia8::JetMatchingMadgraph * fJetMatchingPy8InternalHook
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:207
uint16_t size_type
PowhegHooks * fEmissionVetoHook
virtual void beforeHadronization(lhef::LHEEvent *lhee)
std::auto_ptr< HepMC::GenEvent > & event()
std::auto_ptr< LHAupLesHouches > lhaUP
lhef::LHEEvent * lheEvent()
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
Pythia8::amcnlo_unitarised_interface * fMergingHook
JetMatchingHook * fJetMatchingHook
vector< float > DJR
void resetMatchingStatus()
HepMC::Pythia8ToHepMC toHepMC
bool Pythia8Hadronizer::initializeForExternalPartons ( )

Definition at line 381 of file Pythia8Hadronizer.cc.

References edm::errors::Configuration, Exception, gen::Py8InterfaceBase::fDecayer, fEmissionVetoHook, fEmissionVetoHook1, fJetMatchingHook, fJetMatchingPy8InternalHook, gen::Py8InterfaceBase::fMasterGen, fMergingHook, JetMatchingHook::init(), lhaUP, LHEInputFileName, gen::BaseHadronizer::lheRunInfo(), gen::Py8InterfaceBase::pythiaPylistVerbosity, and ntuplemaker::status.

382 {
383 
384  edm::LogInfo("Pythia8Interface") << "Initializing for external partons";
385 
386  bool status = false, status1 = false;
387 
388  if((fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) && !fEmissionVetoHook) {
389 
391  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
392  <<" Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible are : jetMatching, emissionVeto1 \n";
393 
394  fEmissionVetoHook = new PowhegHooks();
395 
396  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
397  fMasterGen->setUserHooksPtr(fEmissionVetoHook);
398 
399  }
400 
401  //adapted from main89.cc in pythia8 examples
402  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
403  bool internalMerging = !(fMasterGen->settings.word("Merging:Process").compare("void")==0);
404 
405  if (internalMatching && internalMerging) {
406  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
407  <<" Only one jet matching/merging scheme can be used at a time. \n";
408  }
409 
410  if (internalMatching && !fJetMatchingPy8InternalHook) {
411  fJetMatchingPy8InternalHook = new Pythia8::JetMatchingMadgraph;
412  fMasterGen->setUserHooksPtr(fJetMatchingPy8InternalHook);
413  }
414 
415  if (internalMerging && !fMergingHook) {
416  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
417  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
418  1 :
419  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
420  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
421  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
422  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
423  2 :
424  0 );
425  fMergingHook = new Pythia8::amcnlo_unitarised_interface(scheme);
426  fMasterGen->setUserHooksPtr(fMergingHook);
427  }
428 
429 
430  if(LHEInputFileName != string()) {
431 
432  edm::LogInfo("Pythia8Interface") << "Initialize direct pythia8 reading from LHE file "
433  << LHEInputFileName;
434  edm::LogInfo("Pythia8Interface") << "Some LHE information can be not stored";
435  fMasterGen->settings.mode("Beams:frameType", 4);
436  fMasterGen->settings.word("Beams:LHEF", LHEInputFileName);
437  status = fMasterGen->init();
438 
439  } else {
440 
441  lhaUP.reset(new LHAupLesHouches());
442  lhaUP->setScalesFromLHEF(fMasterGen->settings.flag("Beams:setProductionScalesFromLHEF"));
443  lhaUP->loadRunInfo(lheRunInfo());
444 
445  if ( fJetMatchingHook )
446  {
448  }
449 
450  fMasterGen->settings.mode("Beams:frameType", 5);
451  fMasterGen->setLHAupPtr(lhaUP.get());
452  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
453  status = fMasterGen->init();
454  }
455 
456  if ( pythiaPylistVerbosity > 10 )
457  {
458  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
459  fMasterGen->settings.listAll();
460  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
461  fMasterGen->particleData.listAll();
462  }
463 
464  // init decayer
465  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
466  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
467  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
468  status1 = fDecayer->init();
469 
470  return (status&&status1);
471 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
EmissionVetoHook1 * fEmissionVetoHook1
Pythia8::JetMatchingMadgraph * fJetMatchingPy8InternalHook
PowhegHooks * fEmissionVetoHook
std::auto_ptr< LHAupLesHouches > lhaUP
unsigned int pythiaPylistVerbosity
lhef::LHERunInfo * lheRunInfo()
virtual void init(lhef::LHERunInfo *runInfo)
Pythia8::amcnlo_unitarised_interface * fMergingHook
JetMatchingHook * fJetMatchingHook
std::auto_ptr< Pythia8::Pythia > fDecayer
tuple status
Definition: ntuplemaker.py:245
bool Pythia8Hadronizer::initializeForInternalPartons ( )
overridevirtual

Implements gen::Py8InterfaceBase.

Definition at line 332 of file Pythia8Hadronizer.cc.

References comEnergy, edm::errors::Configuration, ElectronPositron, Exception, gen::Py8InterfaceBase::fDecayer, fInitialState, gen::Py8InterfaceBase::fMasterGen, PP, PPbar, gen::Py8InterfaceBase::pythiaPylistVerbosity, and ntuplemaker::status.

333 {
334 
335  bool status = false, status1 = false;
336 
337  if ( fInitialState == PP ) // default
338  {
339  fMasterGen->settings.mode("Beams:idA", 2212);
340  fMasterGen->settings.mode("Beams:idB", 2212);
341  }
342  else if ( fInitialState == PPbar )
343  {
344  fMasterGen->settings.mode("Beams:idA", 2212);
345  fMasterGen->settings.mode("Beams:idB", -2212);
346  }
347  else if ( fInitialState == ElectronPositron )
348  {
349  fMasterGen->settings.mode("Beams:idA", 11);
350  fMasterGen->settings.mode("Beams:idB", -11);
351  }
352  else
353  {
354  // throw on unknown initial state !
355  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
356  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
357  }
358 
359  fMasterGen->settings.parm("Beams:eCM", comEnergy);
360  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
361  status = fMasterGen->init();
362 
363  if ( pythiaPylistVerbosity > 10 )
364  {
365  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
366  fMasterGen->settings.listAll();
367  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
368  fMasterGen->particleData.listAll();
369  }
370 
371  // init decayer
372  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
373  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
374  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
375  status1 = fDecayer->init();
376 
377  return (status&&status1);
378 }
double comEnergy
Center-of-Mass energy.
std::auto_ptr< Pythia8::Pythia > fMasterGen
unsigned int pythiaPylistVerbosity
std::auto_ptr< Pythia8::Pythia > fDecayer
tuple status
Definition: ntuplemaker.py:245
bool Pythia8Hadronizer::residualDecay ( )
virtual

Definition at line 571 of file Pythia8Hadronizer.cc.

References gen::BaseHadronizer::event(), gen::Py8InterfaceBase::fDecayer, gen::Py8InterfaceBase::fMasterGen, GenParticle::GenParticle, query::result, and gen::Py8InterfaceBase::toHepMC.

572 {
573 
574  Event* pythiaEvent = &(fMasterGen->event);
575 
576  int NPartsBeforeDecays = pythiaEvent->size();
577  int NPartsAfterDecays = event().get()->particles_size();
578 
579  if(NPartsAfterDecays == NPartsBeforeDecays) return true;
580 
581  bool result = true;
582 
583  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
584  {
585 
586  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
587 
588  if ( part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id()) )
589  {
590  fDecayer->event.reset();
591  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
592  part->momentum().x(),
593  part->momentum().y(),
594  part->momentum().z(),
595  part->momentum().t(),
596  part->generated_mass() );
597  HepMC::GenVertex* ProdVtx = part->production_vertex();
598  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
599  ProdVtx->position().z(), ProdVtx->position().t() );
600  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
601  fDecayer->event.append( py8part );
602  int nentries = fDecayer->event.size();
603  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
604  fDecayer->next();
605  int nentries1 = fDecayer->event.size();
606  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
607 
608  part->set_status(2);
609 
610  result = toHepMC.fill_next_event( *(fDecayer.get()), event().get(), -1, true, part);
611 
612  }
613  }
614 
615  return result;
616 
617 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
std::auto_ptr< HepMC::GenEvent > & event()
tuple result
Definition: query.py:137
part
Definition: HCALResponse.h:20
std::auto_ptr< Pythia8::Pythia > fDecayer
HepMC::Pythia8ToHepMC toHepMC
void Pythia8Hadronizer::statistics ( )
overridevirtual

Reimplemented from gen::Py8InterfaceBase.

Definition at line 474 of file Pythia8Hadronizer.cc.

References fEmissionVetoHook, gen::Py8InterfaceBase::fMasterGen, nFSRveto, nISRveto, gen::BaseHadronizer::runInfo(), and GenRunInfoProduct::setInternalXSec().

475 {
476  fMasterGen->stat();
477 
478  if(fEmissionVetoHook) {
479  edm::LogPrint("Pythia8Interface") << "\n"
480  << "Number of ISR vetoed = " << nISRveto;
481  edm::LogPrint("Pythia8Interface")
482  << "Number of FSR vetoed = " << nFSRveto;
483  }
484 
485  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
486  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
487  double err = fMasterGen->info.sigmaErr(); // cross section err in mb
488  err *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
490 }
std::auto_ptr< Pythia8::Pythia > fMasterGen
void setInternalXSec(const XSec &xsec)
PowhegHooks * fEmissionVetoHook
GenRunInfoProduct & runInfo()

Member Data Documentation

double Pythia8Hadronizer::comEnergy
private

Center-of-Mass energy.

Definition at line 82 of file Pythia8Hadronizer.cc.

Referenced by initializeForInternalPartons().

vector<float> Pythia8Hadronizer::DJR
private

Definition at line 123 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), and hadronize().

int Pythia8Hadronizer::EV1_emittedMode
private

Definition at line 115 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::EV1_maxVetoCount
private

Definition at line 112 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

bool Pythia8Hadronizer::EV1_MPIvetoOn
private

Definition at line 117 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::EV1_nFinal
private

Definition at line 110 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::EV1_pTdefMode
private

Definition at line 116 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::EV1_pTempMode
private

Definition at line 114 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::EV1_pThardMode
private

Definition at line 113 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

bool Pythia8Hadronizer::EV1_vetoOn
private

Definition at line 111 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

double Pythia8Hadronizer::fBeam1PZ
private

Definition at line 90 of file Pythia8Hadronizer.cc.

double Pythia8Hadronizer::fBeam2PZ
private

Definition at line 91 of file Pythia8Hadronizer.cc.

PowhegHooks* Pythia8Hadronizer::fEmissionVetoHook
private
EmissionVetoHook1* Pythia8Hadronizer::fEmissionVetoHook1
private
int Pythia8Hadronizer::fInitialState
private

Definition at line 88 of file Pythia8Hadronizer.cc.

Referenced by initializeForInternalPartons(), and Pythia8Hadronizer().

JetMatchingHook* Pythia8Hadronizer::fJetMatchingHook
private
Pythia8::JetMatchingMadgraph* Pythia8Hadronizer::fJetMatchingPy8InternalHook
private

Definition at line 102 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

Pythia8::amcnlo_unitarised_interface* Pythia8Hadronizer::fMergingHook
private

Definition at line 103 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

UserHooks* Pythia8Hadronizer::fReweightPtHatRapUserHook
private

Definition at line 97 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

UserHooks* Pythia8Hadronizer::fReweightRapUserHook
private

Definition at line 96 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

UserHooks* Pythia8Hadronizer::fReweightUserHook
private

Definition at line 95 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

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

Definition at line 85 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

string Pythia8Hadronizer::LHEInputFileName
private

Definition at line 84 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and initializeForExternalPartons().

int Pythia8Hadronizer::nFSRveto
private

Definition at line 128 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and statistics().

int Pythia8Hadronizer::NHooks
private

Definition at line 130 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer().

int Pythia8Hadronizer::nISRveto
private

Definition at line 127 of file Pythia8Hadronizer.cc.

Referenced by hadronize(), and statistics().

int Pythia8Hadronizer::nME
private

Definition at line 124 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), and hadronize().

int Pythia8Hadronizer::nMEFiltered
private

Definition at line 125 of file Pythia8Hadronizer.cc.

Referenced by finalizeEvent(), and hadronize().

const std::vector< std::string > Pythia8Hadronizer::p8SharedResources = { edm::SharedResourceNames::kPythia8 }
staticprivate

Definition at line 119 of file Pythia8Hadronizer.cc.

std::string Pythia8Hadronizer::slhafile_
private

Definition at line 121 of file Pythia8Hadronizer.cc.

Referenced by Pythia8Hadronizer(), and ~Pythia8Hadronizer().