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