CMS 3D CMS Logo

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 
27 
28 // Emission Veto Hooks
29 //
30 #include "Pythia8Plugins/PowhegHooks.h"
32 
33 // Resonance scale hook
36 
37 //decay filter hook
39 
40 //decay filter hook
42 
43 // EvtGen plugin
44 //
45 #include "Pythia8Plugins/EvtGen.h"
46 
52 
55 
58 
60 
61 #include "HepPID/ParticleIDTranslations.hh"
62 
64 
65 namespace CLHEP {
66  class HepRandomEngine;
67 }
68 
69 using namespace gen;
70 
71 
73 
74  public:
75 
76  Pythia8Hadronizer(const edm::ParameterSet &params);
78 
79  bool initializeForInternalPartons() override;
80  bool initializeForExternalPartons();
81 
82  bool generatePartonsAndHadronize() override;
83  bool hadronize();
84 
85  virtual bool residualDecay();
86 
87  void finalizeEvent() override;
88 
89  void statistics() override;
90 
91  const char *classname() const override { return "Pythia8Hadronizer"; }
92 
93  GenLumiInfoHeader *getGenLumiInfoHeader() const override;
94 
95  private:
96 
97  virtual void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { p8SetRandomEngine(v); }
98  virtual std::vector<std::string> const& doSharedResources() const override { return p8SharedResources; }
99 
101  double comEnergy;
102 
104  std::auto_ptr<LHAupLesHouches> lhaUP;
105 
106  enum { PP, PPbar, ElectronPositron };
107  int fInitialState ; // pp, ppbar, or e-e+
108 
109  double fBeam1PZ;
110  double fBeam2PZ;
111 
112  //helper class to allow multiple user hooks simultaneously
113  std::auto_ptr<MultiUserHook> fMultiUserHook;
114 
115  // Reweight user hooks
116  //
117  std::auto_ptr<UserHooks> fReweightUserHook;
118  std::auto_ptr<UserHooks> fReweightRapUserHook;
119  std::auto_ptr<UserHooks> fReweightPtHatRapUserHook;
120 
121  // PS matching prototype
122  //
123  std::auto_ptr<JetMatchingHook> fJetMatchingHook;
124  std::auto_ptr<Pythia8::JetMatchingMadgraph> fJetMatchingPy8InternalHook;
125  std::auto_ptr<Pythia8::amcnlo_unitarised_interface> fMergingHook;
126 
127  // Emission Veto Hooks
128  //
129  std::auto_ptr<PowhegHooks> fEmissionVetoHook;
130  std::auto_ptr<EmissionVetoHook1> fEmissionVetoHook1;
131 
132  // Resonance scale hook
133  std::auto_ptr<PowhegResHook> fPowhegResHook;
134  std::auto_ptr<PowhegHooksBB4L> fPowhegHooksBB4L;
135 
136  //resonance decay filter hook
137  std::auto_ptr<ResonanceDecayFilterHook> fResonanceDecayFilterHook;
138 
139  //PT filter hook
140  std::auto_ptr<PTFilterHook> fPTFilterHook;
141 
150 
151  static const std::vector<std::string> p8SharedResources;
152 
153  vector<float> DJR;
154  int nME;
156 
157  int nISRveto;
158  int nFSRveto;
160 
161 };
162 
164 
166  Py8InterfaceBase(params),
167  comEnergy(params.getParameter<double>("comEnergy")),
168  LHEInputFileName(params.getUntrackedParameter<std::string>("LHEInputFileName","")),
169  fInitialState(PP),
170  nME(-1), nMEFiltered(-1), nISRveto(0), nFSRveto(0), nFSRvetoBB4L(0)
171 {
172 
173  // J.Y.: the following 3 parameters are hacked "for a reason"
174  //
175  if ( params.exists( "PPbarInitialState" ) )
176  {
177  if ( fInitialState == PP )
178  {
180  edm::LogImportant("GeneratorInterface|Pythia8Interface")
181  << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
182  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
183  }
184  else
185  {
186  // probably need to throw on attempt to override ?
187  }
188  }
189  else if ( params.exists( "ElectronPositronInitialState" ) )
190  {
191  if ( fInitialState == PP )
192  {
194  edm::LogInfo("GeneratorInterface|Pythia8Interface")
195  << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
196  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
197  }
198  else
199  {
200  // probably need to throw on attempt to override ?
201  }
202  }
203  else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
204  {
205  // throw on unknown initial state !
206  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
207  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
208  }
209 
210  // Reweight user hook
211  //
212  if( params.exists( "reweightGen" ) )
214  if( params.exists( "reweightGenRap" ) )
215  {
216  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
217  edm::ParameterSet rgrParams =
218  params.getParameter<edm::ParameterSet>("reweightGenRap");
219  fReweightRapUserHook.reset(
220  new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
221  rgrParams.getParameter<double>("yLabPower"),
222  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
223  rgrParams.getParameter<double>("yCMPower"),
224  rgrParams.getParameter<double>("pTHatMin"),
225  rgrParams.getParameter<double>("pTHatMax"))
226  );
227  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
228  }
229  if( params.exists( "reweightGenPtHatRap" ) )
230  {
231  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
232  edm::ParameterSet rgrParams =
233  params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
235  new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
236  rgrParams.getParameter<double>("yLabPower"),
237  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
238  rgrParams.getParameter<double>("yCMPower"),
239  rgrParams.getParameter<double>("pTHatMin"),
240  rgrParams.getParameter<double>("pTHatMax"))
241  );
242  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
243  }
244 
245  if( params.exists( "useUserHook" ) )
246  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
247  <<" Obsolete parameter: useUserHook \n Please use the actual one instead \n";
248 
249  // PS matching prototype
250  //
251  if ( params.exists("jetMatching") )
252  {
253  edm::ParameterSet jmParams =
254  params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
255  std::string scheme = jmParams.getParameter<std::string>("scheme");
256  if ( scheme == "Madgraph" || scheme == "MadgraphFastJet" )
257  {
258  fJetMatchingHook.reset(new JetMatchingHook( jmParams, &fMasterGen->info ));
259  }
260  }
261 
262  // Pythia8Interface emission veto
263  //
264  if ( params.exists("emissionVeto1") )
265  {
266  EV1_nFinal = -1;
267  if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
268  EV1_vetoOn = true;
269  if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
270  EV1_maxVetoCount = 10;
271  if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
272  EV1_pThardMode = 1;
273  if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
274  EV1_pTempMode = 0;
275  if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
276  if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
277  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
278  <<" Wrong value for EV1_pTempMode code\n";
279  EV1_emittedMode = 0;
280  if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
281  EV1_pTdefMode = 1;
282  if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
283  EV1_MPIvetoOn = false;
284  if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
288  }
289 
290 }
291 
292 
294 {
295 
296 }
297 
299 {
300 
301  bool status = false, status1 = false;
302 
303  if (lheFile_.empty()) {
304  if ( fInitialState == PP ) // default
305  {
306  fMasterGen->settings.mode("Beams:idA", 2212);
307  fMasterGen->settings.mode("Beams:idB", 2212);
308  }
309  else if ( fInitialState == PPbar )
310  {
311  fMasterGen->settings.mode("Beams:idA", 2212);
312  fMasterGen->settings.mode("Beams:idB", -2212);
313  }
314  else if ( fInitialState == ElectronPositron )
315  {
316  fMasterGen->settings.mode("Beams:idA", 11);
317  fMasterGen->settings.mode("Beams:idB", -11);
318  }
319  else
320  {
321  // throw on unknown initial state !
322  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
323  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
324  }
325  fMasterGen->settings.parm("Beams:eCM", comEnergy);
326  }
327  else {
328  fMasterGen->settings.mode("Beams:frameType", 4);
329  fMasterGen->settings.word("Beams:LHEF", lheFile_);
330  }
331 
332  fMultiUserHook.reset(new MultiUserHook);
333 
334  if(fReweightUserHook.get()) fMultiUserHook->addHook(fReweightUserHook.get());
335  if(fReweightRapUserHook.get()) fMultiUserHook->addHook(fReweightRapUserHook.get());
337  if(fJetMatchingHook.get()) fMultiUserHook->addHook(fJetMatchingHook.get());
338  if(fEmissionVetoHook1.get()) {
339  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
340  fMultiUserHook->addHook(fEmissionVetoHook1.get());
341  }
342 
343  if (fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) {
344 
345  if(fJetMatchingHook.get() || fEmissionVetoHook1.get())
346  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
347  <<" Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible are : jetMatching, emissionVeto1 \n";
348 
349  fEmissionVetoHook.reset(new PowhegHooks());
350 
351  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
352  fMultiUserHook->addHook(fEmissionVetoHook.get());
353  }
354 
355  bool PowhegRes = fMasterGen->settings.flag("POWHEGres:calcScales");
356  if (PowhegRes) {
357  edm::LogInfo("Pythia8Interface") << "Turning on resonance scale setting from CMSSW Pythia8Interface";
358  fPowhegResHook.reset(new PowhegResHook());
359  fMultiUserHook->addHook(fPowhegResHook.get());
360  }
361 
362  bool PowhegBB4L = fMasterGen->settings.flag("POWHEG:bb4l");
363  if (PowhegBB4L) {
364  edm::LogInfo("Pythia8Interface") << "Turning on BB4l hook from CMSSW Pythia8Interface";
365  fPowhegHooksBB4L.reset(new PowhegHooksBB4L());
366  fMultiUserHook->addHook(fPowhegHooksBB4L.get());
367  }
368 
369  //adapted from main89.cc in pythia8 examples
370  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
371  bool internalMerging = !(fMasterGen->settings.word("Merging:Process").compare("void")==0);
372 
373  if (internalMatching && internalMerging) {
374  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
375  <<" Only one jet matching/merging scheme can be used at a time. \n";
376  }
377 
378  if (internalMatching) {
379  fJetMatchingPy8InternalHook.reset(new Pythia8::JetMatchingMadgraph);
381  }
382 
383  if (internalMerging) {
384  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
385  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
386  1 :
387  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
388  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
389  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
390  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
391  2 :
392  0 );
393  fMergingHook.reset(new Pythia8::amcnlo_unitarised_interface(scheme));
394  fMultiUserHook->addHook(fMergingHook.get());
395  }
396 
397  bool resonanceDecayFilter = fMasterGen->settings.flag("ResonanceDecayFilter:filter");
398  if (resonanceDecayFilter) {
401  }
402 
403  bool PTFilter = fMasterGen->settings.flag("PTFilter:filter");
404  if (PTFilter) {
405  fPTFilterHook.reset(new PTFilterHook);
406  fMultiUserHook->addHook(fPTFilterHook.get());
407  }
408 
409  if (fMultiUserHook->nHooks()>0) {
410  fMasterGen->setUserHooksPtr(fMultiUserHook.get());
411  }
412 
413  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
414  status = fMasterGen->init();
415 
416  //clean up temp file
417  if (!slhafile_.empty()) {
418  std::remove(slhafile_.c_str());
419  }
420 
421  if ( pythiaPylistVerbosity > 10 )
422  {
423  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
424  fMasterGen->settings.listAll();
425  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
426  fMasterGen->particleData.listAll();
427  }
428 
429  // init decayer
430  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
431  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
432  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
433  status1 = fDecayer->init();
434 
435  if (useEvtGen) {
436  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
437 
438  evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile.c_str(), evtgenPdlFile.c_str()));
439 
440  for (unsigned int i=0; i<evtgenUserFiles.size(); i++) {
441  edm::FileInPath evtgenUserFile(evtgenUserFiles.at(i));
442  evtgenDecays->readDecayFile(evtgenUserFile.fullPath().c_str());
443  }
444 
445  }
446 
447  return (status&&status1);
448 }
449 
450 
452 {
453 
454  edm::LogInfo("Pythia8Interface") << "Initializing for external partons";
455 
456  bool status = false, status1 = false;
457 
458  fMultiUserHook.reset(new MultiUserHook);
459 
460  if(fReweightUserHook.get()) fMultiUserHook->addHook(fReweightUserHook.get());
461  if(fReweightRapUserHook.get()) fMultiUserHook->addHook(fReweightRapUserHook.get());
463  if(fJetMatchingHook.get()) fMultiUserHook->addHook(fJetMatchingHook.get());
464  if(fEmissionVetoHook1.get()) {
465  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
466  fMultiUserHook->addHook(fEmissionVetoHook1.get());
467  }
468 
469  if (fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) {
470 
471  if(fJetMatchingHook.get() || fEmissionVetoHook1.get())
472  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
473  <<" Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible are : jetMatching, emissionVeto1 \n";
474 
475  fEmissionVetoHook.reset(new PowhegHooks());
476 
477  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
478  fMultiUserHook->addHook(fEmissionVetoHook.get());
479  }
480 
481  bool PowhegRes = fMasterGen->settings.flag("POWHEGres:calcScales");
482  if (PowhegRes) {
483  edm::LogInfo("Pythia8Interface") << "Turning on resonance scale setting from CMSSW Pythia8Interface";
484  fPowhegResHook.reset(new PowhegResHook());
485  fMultiUserHook->addHook(fPowhegResHook.get());
486  }
487 
488  bool PowhegBB4L = fMasterGen->settings.flag("POWHEG:bb4l");
489  if (PowhegBB4L) {
490  edm::LogInfo("Pythia8Interface") << "Turning on BB4l hook from CMSSW Pythia8Interface";
491  fPowhegHooksBB4L.reset(new PowhegHooksBB4L());
492  fMultiUserHook->addHook(fPowhegHooksBB4L.get());
493  }
494 
495  //adapted from main89.cc in pythia8 examples
496  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
497  bool internalMerging = !(fMasterGen->settings.word("Merging:Process").compare("void")==0);
498 
499  if (internalMatching && internalMerging) {
500  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
501  <<" Only one jet matching/merging scheme can be used at a time. \n";
502  }
503 
504  if (internalMatching) {
505  fJetMatchingPy8InternalHook.reset(new Pythia8::JetMatchingMadgraph);
507  }
508 
509  if (internalMerging) {
510  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
511  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
512  1 :
513  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
514  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
515  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
516  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
517  2 :
518  0 );
519  fMergingHook.reset(new Pythia8::amcnlo_unitarised_interface(scheme));
520  fMultiUserHook->addHook(fMergingHook.get());
521  }
522 
523  bool resonanceDecayFilter = fMasterGen->settings.flag("ResonanceDecayFilter:filter");
524  if (resonanceDecayFilter) {
527  }
528 
529  bool PTFilter = fMasterGen->settings.flag("PTFilter:filter");
530  if (PTFilter) {
531  fPTFilterHook.reset(new PTFilterHook);
532  fMultiUserHook->addHook(fPTFilterHook.get());
533  }
534 
535  if (fMultiUserHook->nHooks()>0) {
536  fMasterGen->setUserHooksPtr(fMultiUserHook.get());
537  }
538 
539  if(LHEInputFileName != std::string()) {
540 
541  edm::LogInfo("Pythia8Interface") << "Initialize direct pythia8 reading from LHE file "
542  << LHEInputFileName;
543  edm::LogInfo("Pythia8Interface") << "Some LHE information can be not stored";
544  fMasterGen->settings.mode("Beams:frameType", 4);
545  fMasterGen->settings.word("Beams:LHEF", LHEInputFileName);
546  status = fMasterGen->init();
547 
548  } else {
549 
550  lhaUP.reset(new LHAupLesHouches());
551  lhaUP->setScalesFromLHEF(fMasterGen->settings.flag("Beams:setProductionScalesFromLHEF"));
552  lhaUP->loadRunInfo(lheRunInfo());
553 
554  if ( fJetMatchingHook.get() )
555  {
556  fJetMatchingHook->init ( lheRunInfo() );
557  }
558 
559  fMasterGen->settings.mode("Beams:frameType", 5);
560  fMasterGen->setLHAupPtr(lhaUP.get());
561  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
562  status = fMasterGen->init();
563  }
564 
565  //clean up temp file
566  if (!slhafile_.empty()) {
567  std::remove(slhafile_.c_str());
568  }
569 
570  if ( pythiaPylistVerbosity > 10 )
571  {
572  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
573  fMasterGen->settings.listAll();
574  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
575  fMasterGen->particleData.listAll();
576  }
577 
578  // init decayer
579  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
580  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
581  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
582  status1 = fDecayer->init();
583 
584  if (useEvtGen) {
585  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
586 
587  std::string evtgenpath(getenv("EVTGENDATA"));
588  evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile.c_str(), evtgenPdlFile.c_str()));
589 
590  for (unsigned int i=0; i<evtgenUserFiles.size(); i++) {
591  edm::FileInPath evtgenUserFile(evtgenUserFiles.at(i));
592  evtgenDecays->readDecayFile(evtgenUserFile.fullPath().c_str());
593  }
594 
595  }
596 
597  return (status&&status1);
598 }
599 
600 
602 {
603  fMasterGen->stat();
604 
605  if(fEmissionVetoHook.get()) {
606  edm::LogPrint("Pythia8Interface") << "\n"
607  << "Number of ISR vetoed = " << nISRveto;
608  edm::LogPrint("Pythia8Interface")
609  << "Number of FSR vetoed = " << nFSRveto;
610  }
611 
612  if(fPowhegHooksBB4L.get()) {
613  edm::LogInfo("Pythia8Interface") << "\n"
614  << "BB4L: Number of FSR vetoed = " << nFSRvetoBB4L;
615  }
616 
617  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
618  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
619  double err = fMasterGen->info.sigmaErr(); // cross section err in mb
620  err *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
622 }
623 
624 
626 {
627 
628  DJR.resize(0);
629  nME = -1;
630  nMEFiltered = -1;
631 
632  if ( fJetMatchingHook.get() )
633  {
634  fJetMatchingHook->resetMatchingStatus();
635  fJetMatchingHook->beforeHadronization( lheEvent() );
636  }
637 
638  if (!fMasterGen->next()) return false;
639 
640  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
641  if (fMergingHook.get()) {
642  mergeweight *= fMergingHook->getNormFactor();
643  }
644 
645  //protect against 0-weight from ckkw or similar
646  if (std::abs(mergeweight)==0.)
647  {
648  event().reset();
649  return false;
650  }
651 
652  if (fJetMatchingPy8InternalHook.get()) {
653  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
654  //cap size of djr vector to save storage space (keep only up to first 6 elements)
655  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
656  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
657  DJR.push_back(djrmatch[idjr]);
658  }
659 
660  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
661  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
662  }
663 
664  if (evtgenDecays.get()) evtgenDecays->decay();
665 
666  event().reset(new HepMC::GenEvent);
667  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
668 
669  if (!py8hepmc) {
670  return false;
671  }
672 
673  //add ckkw/umeps/unlops merging weight
674  if (mergeweight!=1.) {
675  event()->weights()[0] *= mergeweight;
676  }
677 
678  if (fEmissionVetoHook.get()) {
679  nISRveto += fEmissionVetoHook->getNISRveto();
680  nFSRveto += fEmissionVetoHook->getNFSRveto();
681  }
682  if (fPowhegHooksBB4L.get()) {
683  nFSRvetoBB4L += fPowhegHooksBB4L->getNFSRveto();
684  }
685 
686  //fill additional weights for systematic uncertainties
687  if (fMasterGen->info.getWeightsDetailedSize() > 0) {
688  for (const string &key : fMasterGen->info.initrwgt->weightsKeys) {
689  double wgt = (*fMasterGen->info.weights_detailed)[key];
690  event()->weights().push_back(wgt);
691  }
692  }
693  else if (fMasterGen->info.getWeightsCompressedSize() > 0) {
694  for (unsigned int i = 0; i < fMasterGen->info.getWeightsCompressedSize(); i++) {
695  double wgt = fMasterGen->info.getWeightsCompressedValue(i);
696  event()->weights().push_back(wgt);
697  }
698  }
699 
700  // fill shower weights
701  // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
702  if( fMasterGen->info.nWeights() > 1 ){
703  for(int i = 0; i < fMasterGen->info.nWeights(); ++i) {
704  double wgt = fMasterGen->info.weight(i);
705  event()->weights().push_back(wgt);
706  }
707  }
708 
709  return true;
710 
711 }
712 
713 
715 {
716  DJR.resize(0);
717  nME = -1;
718  nMEFiltered = -1;
719  if(LHEInputFileName == std::string()) lhaUP->loadEvent(lheEvent());
720 
721  if ( fJetMatchingHook.get() )
722  {
723  fJetMatchingHook->resetMatchingStatus();
724  fJetMatchingHook->beforeHadronization( lheEvent() );
725  }
726 
727  bool py8next = fMasterGen->next();
728 
729  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
730  if (fMergingHook.get()) {
731  mergeweight *= fMergingHook->getNormFactor();
732  }
733 
734 
735  //protect against 0-weight from ckkw or similar
736  if (!py8next || std::abs(mergeweight)==0.)
737  {
738  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
739  event().reset();
740  return false;
741  }
742 
743  if (fJetMatchingPy8InternalHook.get()) {
744  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
745  //cap size of djr vector to save storage space (keep only up to first 6 elements)
746  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
747  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
748  DJR.push_back(djrmatch[idjr]);
749  }
750 
751  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
752  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
753  }
754 
755  // update LHE matching statistics
756  //
757  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
758 
759  if (evtgenDecays.get()) evtgenDecays->decay();
760 
761  event().reset(new HepMC::GenEvent);
762  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
763  if (!py8hepmc) {
764  return false;
765  }
766 
767  //add ckkw/umeps/unlops merging weight
768  if (mergeweight!=1.) {
769  event()->weights()[0] *= mergeweight;
770  }
771 
772  if (fEmissionVetoHook.get()) {
773  nISRveto += fEmissionVetoHook->getNISRveto();
774  nFSRveto += fEmissionVetoHook->getNFSRveto();
775  }
776 
777  return true;
778 
779 }
780 
781 
783 {
784 
785  Event* pythiaEvent = &(fMasterGen->event);
786 
787  int NPartsBeforeDecays = pythiaEvent->size();
788  int NPartsAfterDecays = event().get()->particles_size();
789 
790  if(NPartsAfterDecays == NPartsBeforeDecays) return true;
791 
792  bool result = true;
793 
794  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
795  {
796 
797  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
798 
799  if ( part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id()) )
800  {
801  fDecayer->event.reset();
802  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
803  part->momentum().x(),
804  part->momentum().y(),
805  part->momentum().z(),
806  part->momentum().t(),
807  part->generated_mass() );
808  HepMC::GenVertex* ProdVtx = part->production_vertex();
809  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
810  ProdVtx->position().z(), ProdVtx->position().t() );
811  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
812  fDecayer->event.append( py8part );
813  int nentries = fDecayer->event.size();
814  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
815  fDecayer->next();
816  int nentries1 = fDecayer->event.size();
817  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
818 
819  part->set_status(2);
820 
821  result = toHepMC.fill_next_event( *(fDecayer.get()), event().get(), -1, true, part);
822 
823  }
824  }
825 
826  return result;
827 
828 }
829 
830 
832 {
833  bool lhe = lheEvent() != 0;
834 
835  // now create the GenEventInfo product from the GenEvent and fill
836  // the missing pieces
837  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
838 
839  // in pythia pthat is used to subdivide samples into different bins
840  // in LHE mode the binning is done by the external ME generator
841  // which is likely not pthat, so only filling it for Py6 internal mode
842  if (!lhe) {
843  eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
844  }
845 
846  eventInfo()->setDJR(DJR);
847  eventInfo()->setNMEPartons(nME);
848  eventInfo()->setNMEPartonsFiltered(nMEFiltered);
849 
850  //******** Verbosity ********
851 
852  if (maxEventsToPrint > 0 &&
856  if (pythiaPylistVerbosity) {
857  fMasterGen->info.list();
858  fMasterGen->event.list();
859  }
860 
861  if (pythiaHepMCVerbosity) {
862  std::cout << "Event process = "
863  << fMasterGen->info.code() << "\n"
864  << "----------------------" << std::endl;
865  event()->print();
866  }
868  std::cout << "Event process = "
869  << fMasterGen->info.code() << "\n"
870  << "----------------------" << std::endl;
871  ascii_io->write_event(event().get());
872  }
873  }
874 }
875 
877  GenLumiInfoHeader *genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
878 
879  //fill lhe headers
880  //*FIXME* initrwgt header is corrupt due to pythia bug
881  for (const std::string &key : fMasterGen->info.headerKeys()) {
882  genLumiInfoHeader->lheHeaders().emplace_back(key,fMasterGen->info.header(key));
883  }
884 
885  //check, if it is not only nominal weight
886  int weights_number = fMasterGen->info.nWeights();
887  if (fMasterGen->info.initrwgt) weights_number += fMasterGen->info.initrwgt->weightsKeys.size();
888  if(weights_number > 1){
889  genLumiInfoHeader->weightNames().reserve(weights_number + 1);
890  genLumiInfoHeader->weightNames().push_back("nominal");
891  }
892 
893  //fill weight names
894  if (fMasterGen->info.initrwgt) {
895  for (const std::string &key : fMasterGen->info.initrwgt->weightsKeys) {
896  std::string weightgroupname;
897  for (const auto &wgtgrp : fMasterGen->info.initrwgt->weightgroups) {
898  const auto &wgtgrpwgt = wgtgrp.second.weights.find(key);
899  if (wgtgrpwgt != wgtgrp.second.weights.end()) {
900  weightgroupname = wgtgrp.first;
901  }
902  }
903 
904  std::ostringstream weightname;
905  weightname << "LHE, id = " << key << ", ";
906  if (!weightgroupname.empty()) {
907  weightname << "group = " << weightgroupname << ", ";
908  }
909  weightname<< fMasterGen->info.initrwgt->weights[key].contents;
910  genLumiInfoHeader->weightNames().push_back(weightname.str());
911  }
912  }
913 
914  //fill shower labels
915  // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
916  // http://home.thep.lu.se/~torbjorn/doxygen/classPythia8_1_1Info.html
917  if( fMasterGen->info.nWeights() > 1 ){
918  for(int i = 0; i < fMasterGen->info.nWeights(); ++i) {
919  genLumiInfoHeader->weightNames().push_back( fMasterGen->info.weightLabel(i) );
920  }
921  }
922 
923  return genLumiInfoHeader;
924 }
925 
928 
929 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::auto_ptr< PowhegHooks > fEmissionVetoHook
virtual bool residualDecay()
double comEnergy
Center-of-Mass energy.
std::auto_ptr< Pythia8::Pythia > fMasterGen
edm::GeneratorFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8GeneratorFilter
std::auto_ptr< ResonanceDecayFilterHook > fResonanceDecayFilterHook
bool initializeForInternalPartons() override
std::auto_ptr< EvtGenDecays > evtgenDecays
std::auto_ptr< PowhegHooksBB4L > fPowhegHooksBB4L
HepMC::IO_AsciiParticles * ascii_io
std::auto_ptr< PTFilterHook > fPTFilterHook
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
GenLumiInfoHeader * getGenLumiInfoHeader() const override
void statistics() override
const std::vector< std::string > & weightNames() const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::string LHEInputFileName
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:207
std::auto_ptr< JetMatchingHook > fJetMatchingHook
std::string lheFile_
void setInternalXSec(const XSec &xsec)
uint16_t size_type
std::auto_ptr< HepMC::GenEvent > & event()
std::vector< std::string > evtgenUserFiles
const char * classname() const override
GenRunInfoProduct & runInfo()
std::auto_ptr< LHAupLesHouches > lhaUP
lhef::LHEEvent * lheEvent()
std::auto_ptr< MultiUserHook > fMultiUserHook
static const std::vector< std::string > p8SharedResources
const std::vector< std::pair< std::string, std::string > > & lheHeaders() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::auto_ptr< UserHooks > fReweightPtHatRapUserHook
std::auto_ptr< PowhegResHook > fPowhegResHook
unsigned int pythiaPylistVerbosity
T min(T a, T b)
Definition: MathUtil.h:58
std::auto_ptr< GenEventInfoProduct > & eventInfo()
lhef::LHERunInfo * lheRunInfo()
bool generatePartonsAndHadronize() override
virtual std::vector< std::string > const & doSharedResources() const override
std::auto_ptr< UserHooks > fReweightUserHook
part
Definition: HCALResponse.h:20
unsigned int maxEventsToPrint
def remove(d, key, TELL=False)
Definition: MatrixUtil.py:209
std::auto_ptr< UserHooks > fReweightRapUserHook
std::auto_ptr< Pythia8::Pythia > fDecayer
Pythia8Hadronizer(const edm::ParameterSet &params)
vector< float > DJR
HepMC::Pythia8ToHepMC toHepMC
static const std::string kPythia8
void finalizeEvent() override
std::auto_ptr< Pythia8::JetMatchingMadgraph > fJetMatchingPy8InternalHook
virtual void doSetRandomEngine(CLHEP::HepRandomEngine *v) override
edm::HadronizerFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8HadronizerFilter
std::auto_ptr< EmissionVetoHook1 > fEmissionVetoHook1
std::auto_ptr< Pythia8::amcnlo_unitarised_interface > fMergingHook