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