CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Pythia8Hadronizer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <string>
4 #include <memory>
5 #include <stdint.h>
6 #include <vector>
7 
8 #include "HepMC/GenEvent.h"
9 #include "HepMC/GenParticle.h"
10 
11 #include "Pythia8/Pythia.h"
12 #include "Pythia8Plugins/HepMC2.h"
13 
14 using namespace Pythia8;
15 
17 
19 
20 // PS matchning prototype
21 //
23 #include "Pythia8Plugins/JetMatching.h"
24 #include "Pythia8Plugins/aMCatNLOHooks.h"
25 
26 // Emission Veto Hooks
27 //
28 #include "Pythia8Plugins/PowhegHooks.h"
30 
31 // EvtGen plugin
32 //
33 #include "Pythia8Plugins/EvtGen.h"
34 
40 
43 
47 
49 
50 #include "HepPID/ParticleIDTranslations.hh"
51 
53 
54 namespace CLHEP {
55  class HepRandomEngine;
56 }
57 
58 using namespace gen;
59 
60 
62 
63  public:
64 
65  Pythia8Hadronizer(const edm::ParameterSet &params);
67 
68  bool initializeForInternalPartons() override;
69  bool initializeForExternalPartons();
70 
71  bool generatePartonsAndHadronize() override;
72  bool hadronize();
73 
74  virtual bool residualDecay();
75 
76  void finalizeEvent() override;
77 
78  void statistics() override;
79 
80  const char *classname() const override { return "Pythia8Hadronizer"; }
81 
82  private:
83 
84  virtual void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { p8SetRandomEngine(v); }
85  virtual std::vector<std::string> const& doSharedResources() const override { return p8SharedResources; }
86 
88  double comEnergy;
89 
91  std::auto_ptr<LHAupLesHouches> lhaUP;
92 
93  enum { PP, PPbar, ElectronPositron };
94  int fInitialState ; // pp, ppbar, or e-e+
95 
96  double fBeam1PZ;
97  double fBeam2PZ;
98 
99  // Reweight user hooks
100  //
101  UserHooks* fReweightUserHook;
102  UserHooks* fReweightRapUserHook;
104 
105  // PS matching prototype
106  //
108  Pythia8::JetMatchingMadgraph *fJetMatchingPy8InternalHook;
109  Pythia8::amcnlo_unitarised_interface *fMergingHook;
110 
111  // Emission Veto Hooks
112  //
113  PowhegHooks* fEmissionVetoHook;
115 
124 
125  static const std::vector<std::string> p8SharedResources;
126 
128 
129  vector<float> DJR;
130  int nME;
132 
133  int nISRveto;
134  int nFSRveto;
135 
136  int NHooks;
137 
138 };
139 
141 
143  BaseHadronizer(params), Py8InterfaceBase(params),
144  comEnergy(params.getParameter<double>("comEnergy")),
145  LHEInputFileName(params.getUntrackedParameter<std::string>("LHEInputFileName","")),
146  fInitialState(PP),
147  fReweightUserHook(0),fReweightRapUserHook(0),fReweightPtHatRapUserHook(0),
148  fJetMatchingHook(0),fJetMatchingPy8InternalHook(0), fMergingHook(0),
149  fEmissionVetoHook(0), fEmissionVetoHook1(0), nME(-1), nMEFiltered(-1), nISRveto(0), nFSRveto(0),
150  NHooks(0)
151 {
152 
153  // J.Y.: the following 3 parameters are hacked "for a reason"
154  //
155  if ( params.exists( "PPbarInitialState" ) )
156  {
157  if ( fInitialState == PP )
158  {
160  edm::LogImportant("GeneratorInterface|Pythia8Interface")
161  << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
162  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
163  }
164  else
165  {
166  // probably need to throw on attempt to override ?
167  }
168  }
169  else if ( params.exists( "ElectronPositronInitialState" ) )
170  {
171  if ( fInitialState == PP )
172  {
174  edm::LogInfo("GeneratorInterface|Pythia8Interface")
175  << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
176  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
177  }
178  else
179  {
180  // probably need to throw on attempt to override ?
181  }
182  }
183  else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
184  {
185  // throw on unknown initial state !
186  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
187  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
188  }
189 
190  if( params.exists( "SLHAFileForPythia8" ) ) {
191  std::string slhafilenameshort = params.getParameter<std::string>("SLHAFileForPythia8");
192  edm::FileInPath f1( slhafilenameshort );
193 
194  fMasterGen->settings.mode("SLHA:readFrom", 2);
195  fMasterGen->settings.word("SLHA:file", f1.fullPath());
196 
198  line != fParameters.end(); ++line ) {
199  if (line->find("SLHA:file") != std::string::npos)
200  throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
201  << "using Pythia8 card SLHA:file and Pythia8Interface card SLHAFileForPythia8"
202  << std::endl;
203  }
204  }
205  else if( params.exists( "SLHATableForPythia8" ) ) {
206  std::string slhatable = params.getParameter<std::string>("SLHATableForPythia8");
207 
208  char tempslhaname[] = "pythia8SLHAtableXXXXXX";
209  int fd = mkstemp(tempslhaname);
210  write(fd,slhatable.c_str(),slhatable.size());
211  close(fd);
212 
213  slhafile_ = tempslhaname;
214 
215  fMasterGen->settings.mode("SLHA:readFrom", 2);
216  fMasterGen->settings.word("SLHA:file", slhafile_);
217 
219  line != fParameters.end(); ++line ) {
220  if (line->find("SLHA:file") != std::string::npos)
221  throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
222  << "using Pythia8 card SLHA:file and Pythia8Interface card SLHATableForPythia8"
223  << std::endl;
224  }
225  }
226 
227  // Reweight user hook
228  //
229  if( params.exists( "reweightGen" ) )
231  if( params.exists( "reweightGenRap" ) )
232  {
233  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
234  edm::ParameterSet rgrParams =
235  params.getParameter<edm::ParameterSet>("reweightGenRap");
237  new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
238  rgrParams.getParameter<double>("yLabPower"),
239  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
240  rgrParams.getParameter<double>("yCMPower"),
241  rgrParams.getParameter<double>("pTHatMin"),
242  rgrParams.getParameter<double>("pTHatMax"));
243  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
244  }
245  if( params.exists( "reweightGenPtHatRap" ) )
246  {
247  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
248  edm::ParameterSet rgrParams =
249  params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
251  new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
252  rgrParams.getParameter<double>("yLabPower"),
253  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
254  rgrParams.getParameter<double>("yCMPower"),
255  rgrParams.getParameter<double>("pTHatMin"),
256  rgrParams.getParameter<double>("pTHatMax"));
257  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
258  }
259 
260  if( params.exists( "useUserHook" ) )
261  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
262  <<" Obsolete parameter: useUserHook \n Please use the actual one instead \n";
263 
264  // PS matching prototype
265  //
266  if ( params.exists("jetMatching") )
267  {
268  edm::ParameterSet jmParams =
269  params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
270  std::string scheme = jmParams.getParameter<std::string>("scheme");
271  if ( scheme == "Madgraph" || scheme == "MadgraphFastJet" )
272  {
273  fJetMatchingHook = new JetMatchingHook( jmParams, &fMasterGen->info );
274  }
275  }
276 
277  // Pythia8Interface emission veto
278  //
279  if ( params.exists("emissionVeto1") )
280  {
281  EV1_nFinal = -1;
282  if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
283  EV1_vetoOn = true;
284  if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
285  EV1_maxVetoCount = 10;
286  if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
287  EV1_pThardMode = 1;
288  if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
289  EV1_pTempMode = 0;
290  if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
291  if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
292  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
293  <<" Wrong value for EV1_pTempMode code\n";
294  EV1_emittedMode = 0;
295  if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
296  EV1_pTdefMode = 1;
297  if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
298  EV1_MPIvetoOn = false;
299  if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
303  }
304 
308  if(fJetMatchingHook) NHooks++;
310  if(NHooks > 1)
311  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
312  <<" Too many User Hooks. \n Please choose one from: reweightGen, reweightGenRap, reweightGenPtHatRap, jetMatching, emissionVeto1 \n";
313  if(fReweightUserHook) fMasterGen->setUserHooksPtr(fReweightUserHook);
316  if(fJetMatchingHook) fMasterGen->setUserHooksPtr(fJetMatchingHook);
317  if(fEmissionVetoHook1) {
318  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
319  fMasterGen->setUserHooksPtr(fEmissionVetoHook1);
320  }
321 
322 }
323 
324 
326 {
327 // do we need to delete UserHooks/JetMatchingHook here ???
330 
331  //clean up temp file
332  if (!slhafile_.empty()) {
333  std::remove(slhafile_.c_str());
334  }
335 
336 }
337 
339 {
340 
341  bool status = false, status1 = false;
342 
343  if ( fInitialState == PP ) // default
344  {
345  fMasterGen->settings.mode("Beams:idA", 2212);
346  fMasterGen->settings.mode("Beams:idB", 2212);
347  }
348  else if ( fInitialState == PPbar )
349  {
350  fMasterGen->settings.mode("Beams:idA", 2212);
351  fMasterGen->settings.mode("Beams:idB", -2212);
352  }
353  else if ( fInitialState == ElectronPositron )
354  {
355  fMasterGen->settings.mode("Beams:idA", 11);
356  fMasterGen->settings.mode("Beams:idB", -11);
357  }
358  else
359  {
360  // throw on unknown initial state !
361  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
362  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
363  }
364 
365  fMasterGen->settings.parm("Beams:eCM", comEnergy);
366  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
367  status = fMasterGen->init();
368 
369  if ( pythiaPylistVerbosity > 10 )
370  {
371  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
372  fMasterGen->settings.listAll();
373  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
374  fMasterGen->particleData.listAll();
375  }
376 
377  // init decayer
378  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
379  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
380  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
381  status1 = fDecayer->init();
382 
383  if (useEvtGen) {
384  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
385 
386  evtgenDecays = new EvtGenDecays(fMasterGen.get(), evtgenDecFile.c_str(), evtgenPdlFile.c_str());
387  evtgenDecays->readDecayFile("evtgen_userfile.dec");
388  }
389 
390  return (status&&status1);
391 }
392 
393 
395 {
396 
397  edm::LogInfo("Pythia8Interface") << "Initializing for external partons";
398 
399  bool status = false, status1 = false;
400 
401  if((fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) && !fEmissionVetoHook) {
402 
404  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
405  <<" Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible are : jetMatching, emissionVeto1 \n";
406 
407  fEmissionVetoHook = new PowhegHooks();
408 
409  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
410  fMasterGen->setUserHooksPtr(fEmissionVetoHook);
411 
412  }
413 
414  //adapted from main89.cc in pythia8 examples
415  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
416  bool internalMerging = !(fMasterGen->settings.word("Merging:Process").compare("void")==0);
417 
418  if (internalMatching && internalMerging) {
419  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
420  <<" Only one jet matching/merging scheme can be used at a time. \n";
421  }
422 
423  if (internalMatching && !fJetMatchingPy8InternalHook) {
424  fJetMatchingPy8InternalHook = new Pythia8::JetMatchingMadgraph;
425  fMasterGen->setUserHooksPtr(fJetMatchingPy8InternalHook);
426  }
427 
428  if (internalMerging && !fMergingHook) {
429  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
430  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
431  1 :
432  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
433  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
434  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
435  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
436  2 :
437  0 );
438  fMergingHook = new Pythia8::amcnlo_unitarised_interface(scheme);
439  fMasterGen->setUserHooksPtr(fMergingHook);
440  }
441 
442 
443  if(LHEInputFileName != std::string()) {
444 
445  edm::LogInfo("Pythia8Interface") << "Initialize direct pythia8 reading from LHE file "
446  << LHEInputFileName;
447  edm::LogInfo("Pythia8Interface") << "Some LHE information can be not stored";
448  fMasterGen->settings.mode("Beams:frameType", 4);
449  fMasterGen->settings.word("Beams:LHEF", LHEInputFileName);
450  status = fMasterGen->init();
451 
452  } else {
453 
454  lhaUP.reset(new LHAupLesHouches());
455  lhaUP->setScalesFromLHEF(fMasterGen->settings.flag("Beams:setProductionScalesFromLHEF"));
456  lhaUP->loadRunInfo(lheRunInfo());
457 
458  if ( fJetMatchingHook )
459  {
461  }
462 
463  fMasterGen->settings.mode("Beams:frameType", 5);
464  fMasterGen->setLHAupPtr(lhaUP.get());
465  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
466  status = fMasterGen->init();
467  }
468 
469  if ( pythiaPylistVerbosity > 10 )
470  {
471  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
472  fMasterGen->settings.listAll();
473  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
474  fMasterGen->particleData.listAll();
475  }
476 
477  // init decayer
478  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
479  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
480  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
481  status1 = fDecayer->init();
482 
483  if (useEvtGen) {
484  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
485 
486  std::string evtgenpath(getenv("EVTGENDATA"));
487  evtgenDecays = new EvtGenDecays(fMasterGen.get(), evtgenDecFile.c_str(), evtgenPdlFile.c_str());
488  evtgenDecays->readDecayFile("evtgen_userfile.dec");
489  }
490 
491  return (status&&status1);
492 }
493 
494 
496 {
497  fMasterGen->stat();
498 
499  if(fEmissionVetoHook) {
500  edm::LogPrint("Pythia8Interface") << "\n"
501  << "Number of ISR vetoed = " << nISRveto;
502  edm::LogPrint("Pythia8Interface")
503  << "Number of FSR vetoed = " << nFSRveto;
504  }
505 
506  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
507  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
508  double err = fMasterGen->info.sigmaErr(); // cross section err in mb
509  err *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
511 }
512 
513 
515 {
516 
517  if (!fMasterGen->next()) return false;
518 
519  if (evtgenDecays) evtgenDecays->decay();
520 
521  event().reset(new HepMC::GenEvent);
522  return toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
523 
524 }
525 
526 
528 {
529  DJR.resize(0);
530  nME = -1;
531  nMEFiltered = -1;
532  if(LHEInputFileName == std::string()) lhaUP->loadEvent(lheEvent());
533 
534  if ( fJetMatchingHook )
535  {
538  }
539 
540  bool py8next = fMasterGen->next();
541 
542  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
543  if (fMergingHook) {
544  mergeweight *= fMergingHook->getNormFactor();
545  }
546 
547 
548  //protect against 0-weight from ckkw or similar
549  if (!py8next || std::abs(mergeweight)==0.)
550  {
551  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
552  event().reset();
553  return false;
554  }
555 
557  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->GetDJR();
558  //cap size of djr vector to save storage space (keep only up to first 6 elements)
559  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
560  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
561  DJR.push_back(djrmatch[idjr]);
562  }
563 
564  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
565  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
566  }
567 
568  // update LHE matching statistics
569  //
570  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
571 
572  if (evtgenDecays) evtgenDecays->decay();
573 
574  event().reset(new HepMC::GenEvent);
575  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
576  if (!py8hepmc) {
577  return false;
578  }
579 
580  //add ckkw/umeps/unlops merging weight
581  if (mergeweight!=1.) {
582  event()->weights()[0] *= mergeweight;
583  }
584 
585  if (fEmissionVetoHook) {
586  nISRveto += fEmissionVetoHook->getNISRveto();
587  nFSRveto += fEmissionVetoHook->getNFSRveto();
588  }
589 
590  return true;
591 
592 }
593 
594 
596 {
597 
598  Event* pythiaEvent = &(fMasterGen->event);
599 
600  int NPartsBeforeDecays = pythiaEvent->size();
601  int NPartsAfterDecays = event().get()->particles_size();
602 
603  if(NPartsAfterDecays == NPartsBeforeDecays) return true;
604 
605  bool result = true;
606 
607  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
608  {
609 
610  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
611 
612  if ( part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id()) )
613  {
614  fDecayer->event.reset();
615  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
616  part->momentum().x(),
617  part->momentum().y(),
618  part->momentum().z(),
619  part->momentum().t(),
620  part->generated_mass() );
621  HepMC::GenVertex* ProdVtx = part->production_vertex();
622  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
623  ProdVtx->position().z(), ProdVtx->position().t() );
624  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
625  fDecayer->event.append( py8part );
626  int nentries = fDecayer->event.size();
627  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
628  fDecayer->next();
629  int nentries1 = fDecayer->event.size();
630  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
631 
632  part->set_status(2);
633 
634  result = toHepMC.fill_next_event( *(fDecayer.get()), event().get(), -1, true, part);
635 
636  }
637  }
638 
639  return result;
640 
641 }
642 
643 
645 {
646  bool lhe = lheEvent() != 0;
647 
648  // now create the GenEventInfo product from the GenEvent and fill
649  // the missing pieces
650  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
651 
652  // in pythia pthat is used to subdivide samples into different bins
653  // in LHE mode the binning is done by the external ME generator
654  // which is likely not pthat, so only filling it for Py6 internal mode
655  if (!lhe) {
656  eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
657  }
658 
659  eventInfo()->setDJR(DJR);
660  eventInfo()->setNMEPartons(nME);
661  eventInfo()->setNMEPartonsFiltered(nMEFiltered);
662 
663  //******** Verbosity ********
664 
665  if (maxEventsToPrint > 0 &&
669  if (pythiaPylistVerbosity) {
670  fMasterGen->info.list(std::cout);
671  fMasterGen->event.list(std::cout);
672  }
673 
674  if (pythiaHepMCVerbosity) {
675  std::cout << "Event process = "
676  << fMasterGen->info.code() << "\n"
677  << "----------------------" << std::endl;
678  event()->print();
679  }
681  std::cout << "Event process = "
682  << fMasterGen->info.code() << "\n"
683  << "----------------------" << std::endl;
684  ascii_io->write_event(event().get());
685  }
686  }
687 }
688 
691 
692 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual bool residualDecay()
ParameterCollector fParameters
double comEnergy
Center-of-Mass energy.
std::auto_ptr< Pythia8::Pythia > fMasterGen
EmissionVetoHook1 * fEmissionVetoHook1
edm::GeneratorFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8GeneratorFilter
bool initializeForInternalPartons() override
UserHooks * fReweightUserHook
HepMC::IO_AsciiParticles * ascii_io
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void statistics() override
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::string LHEInputFileName
Pythia8::JetMatchingMadgraph * fJetMatchingPy8InternalHook
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:207
void setInternalXSec(const XSec &xsec)
uint16_t size_type
PowhegHooks * fEmissionVetoHook
virtual void beforeHadronization(lhef::LHEEvent *lhee)
std::auto_ptr< HepMC::GenEvent > & event()
virtual std::vector< std::string > const & doSharedResources() const override
GenRunInfoProduct & runInfo()
std::auto_ptr< LHAupLesHouches > lhaUP
lhef::LHEEvent * lheEvent()
UserHooks * fReweightPtHatRapUserHook
tuple result
Definition: query.py:137
EvtGenDecays * evtgenDecays
static const std::vector< std::string > p8SharedResources
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int pythiaPylistVerbosity
const char * classname() const override
T min(T a, T b)
Definition: MathUtil.h:58
UserHooks * fReweightRapUserHook
std::auto_ptr< GenEventInfoProduct > & eventInfo()
lhef::LHERunInfo * lheRunInfo()
bool generatePartonsAndHadronize() override
virtual void init(lhef::LHERunInfo *runInfo)
Pythia8::amcnlo_unitarised_interface * fMergingHook
JetMatchingHook * fJetMatchingHook
part
Definition: HCALResponse.h:20
unsigned int maxEventsToPrint
std::auto_ptr< Pythia8::Pythia > fDecayer
Pythia8Hadronizer(const edm::ParameterSet &params)
vector< float > DJR
void resetMatchingStatus()
const_iterator end() const
const_iterator begin() const
HepMC::Pythia8ToHepMC toHepMC
tuple cout
Definition: gather_cfg.py:121
static const std::string kPythia8
void finalizeEvent() override
tuple status
Definition: ntuplemaker.py:245
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
edm::HadronizerFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8HadronizerFilter