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