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