CMS 3D CMS Logo

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