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 
15 
17 
18 // PS matchning prototype
19 //
21 #include "Pythia8Plugins/JetMatching.h"
22 #include "Pythia8Plugins/aMCatNLOHooks.h"
23 
24 // Emission Veto Hooks
25 //
28 
34 
37 
41 
43 
44 #include "HepPID/ParticleIDTranslations.hh"
45 
47 
48 namespace CLHEP {
49  class HepRandomEngine;
50 }
51 
52 using namespace gen;
53 using namespace Pythia8;
54 
56 
57  public:
58 
59  Pythia8Hadronizer(const edm::ParameterSet &params);
61 
62  bool initializeForInternalPartons() override;
63  bool initializeForExternalPartons();
64 
65  bool generatePartonsAndHadronize() override;
66  bool hadronize();
67 
68  virtual bool residualDecay();
69 
70  void finalizeEvent() override;
71 
72  void statistics() override;
73 
74  const char *classname() const override { return "Pythia8Hadronizer"; }
75 
76  private:
77 
78  virtual void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { p8SetRandomEngine(v); }
79  virtual std::vector<std::string> const& doSharedResources() const override { return p8SharedResources; }
80 
82  double comEnergy;
83 
85  std::auto_ptr<LHAupLesHouches> lhaUP;
86 
87  enum { PP, PPbar, ElectronPositron };
88  int fInitialState ; // pp, ppbar, or e-e+
89 
90  double fBeam1PZ;
91  double fBeam2PZ;
92 
93  // Reweight user hooks
94  //
95  UserHooks* fReweightUserHook;
96  UserHooks* fReweightRapUserHook;
98 
99  // PS matching prototype
100  //
102  Pythia8::JetMatchingMadgraph *fJetMatchingPy8InternalHook;
103  Pythia8::amcnlo_unitarised_interface *fMergingHook;
104 
105  // Emission Veto Hooks
106  //
109 
118 
119  static const std::vector<std::string> p8SharedResources;
120 
122 
123  vector<float> DJR;
124  int nME;
126 };
127 
129 
131  BaseHadronizer(params), Py8InterfaceBase(params),
132  comEnergy(params.getParameter<double>("comEnergy")),
133  LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
134  fInitialState(PP),
135  fReweightUserHook(0),fReweightRapUserHook(0),fReweightPtHatRapUserHook(0),
136  fJetMatchingHook(0),fJetMatchingPy8InternalHook(0), fMergingHook(0),
137  fEmissionVetoHook(0),fEmissionVetoHook1(0), nME(-1), nMEFiltered(-1)
138 {
139 
140  // J.Y.: the following 3 parameters are hacked "for a reason"
141  //
142  if ( params.exists( "PPbarInitialState" ) )
143  {
144  if ( fInitialState == PP )
145  {
147  edm::LogImportant("GeneratorInterface|Pythia8Interface")
148  << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
149  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
150  }
151  else
152  {
153  // probably need to throw on attempt to override ?
154  }
155  }
156  else if ( params.exists( "ElectronPositronInitialState" ) )
157  {
158  if ( fInitialState == PP )
159  {
161  edm::LogInfo("GeneratorInterface|Pythia8Interface")
162  << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
163  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
164  }
165  else
166  {
167  // probably need to throw on attempt to override ?
168  }
169  }
170  else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
171  {
172  // throw on unknown initial state !
173  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
174  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
175  }
176 
177  if( params.exists( "SLHAFileForPythia8" ) ) {
178  std::string slhafilenameshort = params.getParameter<string>("SLHAFileForPythia8");
179  edm::FileInPath f1( slhafilenameshort );
180  slhafile_ = f1.fullPath();
181 
182  fMasterGen->settings.mode("SLHA:readFrom", 2);
183  fMasterGen->settings.word("SLHA:file", slhafile_);
184 
186  line != fParameters.end(); ++line ) {
187  if (line->find("SLHA:file") != std::string::npos)
188  throw cms::Exception("PythiaError") << "Attempted to set SLHA file name twice, "
189  << "using Pythia8 card SLHA:file and Pythia8Interface card SLHAFileForPythia8"
190  << std::endl;
191  }
192  }
193 
194  // Reweight user hook
195  //
196  if( params.exists( "reweightGen" ) )
198  if( params.exists( "reweightGenRap" ) )
199  {
200  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
201  edm::ParameterSet rgrParams =
202  params.getParameter<edm::ParameterSet>("reweightGenRap");
204  new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
205  rgrParams.getParameter<double>("yLabPower"),
206  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
207  rgrParams.getParameter<double>("yCMPower"),
208  rgrParams.getParameter<double>("pTHatMin"),
209  rgrParams.getParameter<double>("pTHatMax"));
210  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
211  }
212  if( params.exists( "reweightGenPtHatRap" ) )
213  {
214  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
215  edm::ParameterSet rgrParams =
216  params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
218  new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
219  rgrParams.getParameter<double>("yLabPower"),
220  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
221  rgrParams.getParameter<double>("yCMPower"),
222  rgrParams.getParameter<double>("pTHatMin"),
223  rgrParams.getParameter<double>("pTHatMax"));
224  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
225  }
226 
227  if( params.exists( "useUserHook" ) )
228  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
229  <<" Obsolete parameter: useUserHook \n Please use the actual one instead \n";
230 
231  // PS matching prototype
232  //
233  if ( params.exists("jetMatching") )
234  {
235  edm::ParameterSet jmParams =
236  params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
237  std::string scheme = jmParams.getParameter<std::string>("scheme");
238  if ( scheme == "Madgraph" || scheme == "MadgraphFastJet" )
239  {
240  fJetMatchingHook = new JetMatchingHook( jmParams, &fMasterGen->info );
241  }
242  }
243 
244  // Emission vetos
245  //
246  if ( params.exists("emissionVeto") )
247  {
249  }
250 
251  if ( params.exists("emissionVeto1") )
252  {
253  EV1_nFinal = -1;
254  if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
255  EV1_vetoOn = true;
256  if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
257  EV1_maxVetoCount = 10;
258  if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
259  EV1_pThardMode = 1;
260  if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
261  EV1_pTempMode = 0;
262  if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
263  if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
264  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
265  <<" Wrong value for EV1_pTempMode code\n";
266  EV1_emittedMode = 0;
267  if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
268  EV1_pTdefMode = 1;
269  if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
270  EV1_MPIvetoOn = false;
271  if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
275  }
276 
277  int NHooks=0;
278  if(fReweightUserHook) NHooks++;
279  if(fReweightRapUserHook) NHooks++;
280  if(fReweightPtHatRapUserHook) NHooks++;
281  if(fJetMatchingHook) NHooks++;
282  if(fEmissionVetoHook) NHooks++;
283  if(fEmissionVetoHook1) NHooks++;
284  if(NHooks > 1)
285  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
286  <<" Too many User Hooks. \n Please choose one from: reweightGen, reweightGenRap, reweightGenPtHatRap, jetMatching, emissionVeto, emissionVeto1 \n";
287  if(fReweightUserHook) fMasterGen->setUserHooksPtr(fReweightUserHook);
290  if(fJetMatchingHook) fMasterGen->setUserHooksPtr(fJetMatchingHook);
292  std::cout << "Turning on Emission Veto Hook";
293  if(fEmissionVetoHook1) std::cout << " 1";
294  std::cout << std::endl;
295  if(fEmissionVetoHook) fMasterGen->setUserHooksPtr(fEmissionVetoHook);
296  if(fEmissionVetoHook1) fMasterGen->setUserHooksPtr(fEmissionVetoHook1);
297  }
298 }
299 
300 
302 {
303 // do we need to delete UserHooks/JetMatchingHook here ???
304 
307 }
308 
310 {
311  if ( fInitialState == PP ) // default
312  {
313  //fMasterGen->init(2212, 2212, comEnergy);
314  fMasterGen->settings.mode("Beams:idA", 2212);
315  fMasterGen->settings.mode("Beams:idB", 2212);
316  fMasterGen->settings.parm("Beams:eCM", comEnergy);
317  fMasterGen->init();
318  }
319  else if ( fInitialState == PPbar )
320  {
321  //fMasterGen->init(2212, -2212, comEnergy);
322  fMasterGen->settings.mode("Beams:idA", 2212);
323  fMasterGen->settings.mode("Beams:idB", -2212);
324  fMasterGen->settings.parm("Beams:eCM", comEnergy);
325  fMasterGen->init();
326  }
327  else if ( fInitialState == ElectronPositron )
328  {
329  //fMasterGen->init(11, -11, comEnergy);
330  fMasterGen->settings.mode("Beams:idA", 11);
331  fMasterGen->settings.mode("Beams:idB", -11);
332  fMasterGen->settings.parm("Beams:eCM", comEnergy);
333  fMasterGen->init();
334  }
335  else
336  {
337  // throw on unknown initial state !
338  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
339  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
340  }
341 
342  fMasterGen->settings.listChanged();
343 
344  if ( pythiaPylistVerbosity > 10 )
345  {
346  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
347  fMasterGen->settings.listAll();
348  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
349  fMasterGen->particleData.listAll();
350  }
351 
352  // init decayer
353  //fDecayer->readString("ProcessLevel:all = off"); // trick
354  //fDecayer->readString("ProcessLevel::resonanceDecays=on");
355  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
356  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
357  fDecayer->init();
358 
359  return true;
360 }
361 
362 
364 {
365 
366  std::cout << "Initializing for external partons" << std::endl;
367 
368  //adapted from main89.cc in pythia8 examples
369  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
370  bool internalMerging = !(fMasterGen->settings.word("Merging:Process").compare("void")==0);
371 
372  if (internalMatching && internalMerging) {
373  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
374  <<" Only one jet matching/merging scheme can be used at a time. \n";
375  }
376 
377  if (internalMatching && !fJetMatchingPy8InternalHook) {
378  fJetMatchingPy8InternalHook = new Pythia8::JetMatchingMadgraph;
379  fMasterGen->setUserHooksPtr(fJetMatchingPy8InternalHook);
380  }
381 
382  if (internalMerging && !fMergingHook) {
383  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
384  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
385  1 :
386  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
387  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
388  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
389  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
390  2 :
391  0 );
392  fMergingHook = new Pythia8::amcnlo_unitarised_interface(scheme);
393  fMasterGen->setUserHooksPtr(fMergingHook);
394  }
395 
396 
397  if(LHEInputFileName != string()) {
398 
399  edm::LogInfo("Pythia8Interface") << "Initialize direct pythia8 reading from LHE file "
400  << LHEInputFileName;
401  edm::LogInfo("Pythia8Interface") << "Some LHE information can be not stored";
402  //fMasterGen->init(LHEInputFileName);
403  fMasterGen->settings.mode("Beams:frameType", 4);
404  fMasterGen->settings.word("Beams:LHEF", LHEInputFileName);
405  fMasterGen->init();
406 
407  } else {
408 
409  lhaUP.reset(new LHAupLesHouches());
410  lhaUP->loadRunInfo(lheRunInfo());
411 
412  if ( fJetMatchingHook )
413  {
415  }
416 
417  //fMasterGen->init(lhaUP.get());
418  fMasterGen->settings.mode("Beams:frameType", 5);
419  fMasterGen->setLHAupPtr(lhaUP.get());
420  fMasterGen->init();
421  }
422 
423  if ( pythiaPylistVerbosity > 10 )
424  {
425  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
426  fMasterGen->settings.listAll();
427  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
428  fMasterGen->particleData.listAll();
429  }
430 
431  // init decayer
432  //fDecayer->readString("ProcessLevel:all = off"); // trick
433  //fDecayer->readString("ProcessLevel::resonanceDecays=on");
434  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
435  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
436  fDecayer->init();
437 
438  return true;
439 }
440 
441 
443 {
444  fMasterGen->stat();
445 
446  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
447  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
448  runInfo().setInternalXSec(xsec);
449 }
450 
451 
453 {
454 
455  if (!fMasterGen->next()) return false;
456 
457  event().reset(new HepMC::GenEvent);
458  return toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
459 
460 }
461 
462 
464 {
465  DJR.resize(0);
466  nME = -1;
467  nMEFiltered = -1;
468  if(LHEInputFileName == string()) lhaUP->loadEvent(lheEvent());
469 
470  if ( fJetMatchingHook )
471  {
474  }
475 
476  bool py8next = fMasterGen->next();
477 
478  double mergeweight = fMasterGen.get()->info.mergingWeight();
479 
480  //protect against 0-weight from ckkw or similar
481  if (!py8next || std::abs(mergeweight)==0.)
482  {
483  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
484  event().reset();
485  return false;
486  }
487 
489  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->GetDJR();
490  //cap size of djr vector to save storage space (keep only up to first 6 elements)
491  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
492  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
493  DJR.push_back(djrmatch[idjr]);
494  }
495 
496  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
497  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
498  }
499 
500  // update LHE matching statistics
501  //
502  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
503 
504  event().reset(new HepMC::GenEvent);
505  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
506  if (!py8hepmc) {
507  return false;
508  }
509 
510  //add ckkw merging weight
511  if (mergeweight!=1.) {
512  event()->weights().push_back(mergeweight);
513  }
514 
515  return true;
516 
517 
518 }
519 
520 
522 {
523 
524  Event* pythiaEvent = &(fMasterGen->event);
525 
526  int NPartsBeforeDecays = pythiaEvent->size();
527  int NPartsAfterDecays = event().get()->particles_size();
528  int NewBarcode = NPartsAfterDecays;
529 
530  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
531  {
532 
533  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
534 
535  if ( part->status() == 1 )
536  {
537  fDecayer->event.reset();
538  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
539  part->momentum().x(),
540  part->momentum().y(),
541  part->momentum().z(),
542  part->momentum().t(),
543  part->generated_mass() );
544  HepMC::GenVertex* ProdVtx = part->production_vertex();
545  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
546  ProdVtx->position().z(), ProdVtx->position().t() );
547  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
548  fDecayer->event.append( py8part );
549  int nentries = fDecayer->event.size();
550  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
551  fDecayer->next();
552  int nentries1 = fDecayer->event.size();
553  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
554 
555  part->set_status(2);
556 
557  Particle& py8daughter = fDecayer->event[nentries]; // the 1st daughter
558  HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
559  py8daughter.yProd(),
560  py8daughter.zProd(),
561  py8daughter.tProd()) );
562 
563  DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
564  // I presume (vtx) barcode will be given automatically
565 
566  HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
567 
568  HepMC::GenParticle* daughter =
569  new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
570 
571  NewBarcode++;
572  daughter->suggest_barcode( NewBarcode );
573  DecVtx->add_particle_out( daughter );
574 
575  for ( int ipart1=nentries+1; ipart1<nentries1; ipart1++ )
576  {
577  py8daughter = fDecayer->event[ipart1];
578  HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
579  HepMC::GenParticle* daughterN =
580  new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
581  NewBarcode++;
582  daughterN->suggest_barcode( NewBarcode );
583  DecVtx->add_particle_out( daughterN );
584  }
585 
586  event().get()->add_vertex( DecVtx );
587 
588  }
589  }
590  return true;
591 
592 }
593 
594 
596 {
597  bool lhe = lheEvent() != 0;
598 
599  // now create the GenEventInfo product from the GenEvent and fill
600  // the missing pieces
601  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
602 
603  // in pythia pthat is used to subdivide samples into different bins
604  // in LHE mode the binning is done by the external ME generator
605  // which is likely not pthat, so only filling it for Py6 internal mode
606  if (!lhe) {
607  eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
608  }
609 
610  eventInfo()->setDJR(DJR);
611  eventInfo()->setNMEPartons(nME);
612  eventInfo()->setNMEPartonsFiltered(nMEFiltered);
613 
614  //******** Verbosity ********
615 
616  if (maxEventsToPrint > 0 &&
619  if (pythiaPylistVerbosity) {
620  fMasterGen->info.list(std::cout);
621  fMasterGen->event.list(std::cout);
622  }
623 
624  if (pythiaHepMCVerbosity) {
625  std::cout << "Event process = "
626  << fMasterGen->info.code() << "\n"
627  << "----------------------" << std::endl;
628  event()->print();
629  }
630  }
631 }
632 
635 
636 
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
EmissionVetoHook * fEmissionVetoHook
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void statistics() override
bool exists(std::string const &parameterName) const
checks if a parameter exists
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
virtual void beforeHadronization(lhef::LHEEvent *lhee)
double v[5][pyjets_maxn]
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
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
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
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
edm::HadronizerFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8HadronizerFilter