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 
157 
158  static const std::vector<std::string> p8SharedResources;
159 
160  vector<float> DJR;
161  int nME;
163 
164  int nISRveto;
165  int nFSRveto;
167 
168 };
169 
171 
173  Py8InterfaceBase(params),
174  comEnergy(params.getParameter<double>("comEnergy")),
175  LHEInputFileName(params.getUntrackedParameter<std::string>("LHEInputFileName","")),
176  fInitialState(PP),
177  nME(-1), nMEFiltered(-1), nISRveto(0), nFSRveto(0), nFSRvetoBB4L(0)
178 {
179 
180  // J.Y.: the following 3 parameters are hacked "for a reason"
181  //
182  if ( params.exists( "PPbarInitialState" ) )
183  {
184  if ( fInitialState == PP )
185  {
187  edm::LogImportant("GeneratorInterface|Pythia8Interface")
188  << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
189  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
190  }
191  else
192  {
193  // probably need to throw on attempt to override ?
194  }
195  }
196  else if ( params.exists( "ElectronPositronInitialState" ) )
197  {
198  if ( fInitialState == PP )
199  {
201  edm::LogInfo("GeneratorInterface|Pythia8Interface")
202  << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
203  << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
204  }
205  else
206  {
207  // probably need to throw on attempt to override ?
208  }
209  }
210  else if ( params.exists( "ElectronProtonInitialState" ) || params.exists( "PositronProtonInitialState" ) )
211  {
212  // throw on unknown initial state !
213  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
214  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
215  }
216 
217  // Reweight user hook
218  //
219  if( params.exists( "reweightGen" ) )
220  {
221  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGen";
222  edm::ParameterSet rgParams =
223  params.getParameter<edm::ParameterSet>("reweightGen");
224  fReweightUserHook.reset(
225  new PtHatReweightUserHook(rgParams.getParameter<double>("pTRef"),
226  rgParams.getParameter<double>("power"))
227  );
228  edm::LogInfo("Pythia8Interface") << "End setup for reweightGen";
229  }
230  if( params.exists( "reweightGenEmp" ) )
231  {
232  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenEmp";
233  edm::ParameterSet rgeParams =
234  params.getParameter<edm::ParameterSet>("reweightGenEmp");
236  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenEmp";
237  }
238  if( params.exists( "reweightGenRap" ) )
239  {
240  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
241  edm::ParameterSet rgrParams =
242  params.getParameter<edm::ParameterSet>("reweightGenRap");
243  fReweightRapUserHook.reset(
244  new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
245  rgrParams.getParameter<double>("yLabPower"),
246  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
247  rgrParams.getParameter<double>("yCMPower"),
248  rgrParams.getParameter<double>("pTHatMin"),
249  rgrParams.getParameter<double>("pTHatMax"))
250  );
251  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
252  }
253  if( params.exists( "reweightGenPtHatRap" ) )
254  {
255  edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
256  edm::ParameterSet rgrParams =
257  params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
259  new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
260  rgrParams.getParameter<double>("yLabPower"),
261  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
262  rgrParams.getParameter<double>("yCMPower"),
263  rgrParams.getParameter<double>("pTHatMin"),
264  rgrParams.getParameter<double>("pTHatMax"))
265  );
266  edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
267  }
268 
269  if( params.exists( "useUserHook" ) )
270  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
271  <<" Obsolete parameter: useUserHook \n Please use the actual one instead \n";
272 
273  // PS matching prototype
274  //
275  if ( params.exists("jetMatching") )
276  {
277  edm::ParameterSet jmParams =
278  params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
279  std::string scheme = jmParams.getParameter<std::string>("scheme");
280  if ( scheme == "Madgraph" || scheme == "MadgraphFastJet" )
281  {
282  fJetMatchingHook.reset(new JetMatchingHook( jmParams, &fMasterGen->info ));
283  }
284  }
285 
286  // Pythia8Interface emission veto
287  //
288  if ( params.exists("emissionVeto1") )
289  {
290  EV1_nFinal = -1;
291  if(params.exists("EV1_nFinal")) EV1_nFinal = params.getParameter<int>("EV1_nFinal");
292  EV1_vetoOn = true;
293  if(params.exists("EV1_vetoOn")) EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
294  EV1_maxVetoCount = 10;
295  if(params.exists("EV1_maxVetoCount")) EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
296  EV1_pThardMode = 1;
297  if(params.exists("EV1_pThardMode")) EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
298  EV1_pTempMode = 0;
299  if(params.exists("EV1_pTempMode")) EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
300  if(EV1_pTempMode > 2 || EV1_pTempMode < 0)
301  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
302  <<" Wrong value for EV1_pTempMode code\n";
303  EV1_emittedMode = 0;
304  if(params.exists("EV1_emittedMode")) EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
305  EV1_pTdefMode = 1;
306  if(params.exists("EV1_pTdefMode")) EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
307  EV1_MPIvetoOn = false;
308  if(params.exists("EV1_MPIvetoOn")) EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
312  }
313 
314  if( params.exists( "VinciaPlugin" ) ) {
315  fMasterGen.reset(new Pythia);
316  fvincia.reset(new Vincia::VinciaPlugin(fMasterGen.get()));
317  }
318  if( params.exists( "DirePlugin" ) ) {
319  fMasterGen.reset(new Pythia);
320  fDire.reset(new Pythia8::Dire());
321  fDire->initSettings(*fMasterGen.get());
322  fDire->initShowersAndWeights(*fMasterGen.get(), nullptr, nullptr);
323  }
324 
325 }
326 
327 
329 {
330 
331 }
332 
334 {
335 
336  bool status = false, status1 = false;
337 
338  if (lheFile_.empty()) {
339  if ( fInitialState == PP ) // default
340  {
341  fMasterGen->settings.mode("Beams:idA", 2212);
342  fMasterGen->settings.mode("Beams:idB", 2212);
343  }
344  else if ( fInitialState == PPbar )
345  {
346  fMasterGen->settings.mode("Beams:idA", 2212);
347  fMasterGen->settings.mode("Beams:idB", -2212);
348  }
349  else if ( fInitialState == ElectronPositron )
350  {
351  fMasterGen->settings.mode("Beams:idA", 11);
352  fMasterGen->settings.mode("Beams:idB", -11);
353  }
354  else
355  {
356  // throw on unknown initial state !
357  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
358  <<" UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
359  }
360  fMasterGen->settings.parm("Beams:eCM", comEnergy);
361  }
362  else {
363  fMasterGen->settings.mode("Beams:frameType", 4);
364  fMasterGen->settings.word("Beams:LHEF", lheFile_);
365  }
366 
367  fMultiUserHook.reset(new MultiUserHook);
368 
369  if(fReweightUserHook.get()) fMultiUserHook->addHook(fReweightUserHook.get());
370  if(fReweightEmpUserHook.get()) fMultiUserHook->addHook(fReweightEmpUserHook.get());
371  if(fReweightRapUserHook.get()) fMultiUserHook->addHook(fReweightRapUserHook.get());
373  if(fJetMatchingHook.get()) fMultiUserHook->addHook(fJetMatchingHook.get());
374  if(fEmissionVetoHook1.get()) {
375  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
376  fMultiUserHook->addHook(fEmissionVetoHook1.get());
377  }
378 
379  if (fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) {
380 
381  if(fJetMatchingHook.get() || fEmissionVetoHook1.get())
382  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
383  <<" Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible are : jetMatching, emissionVeto1 \n";
384 
385  fEmissionVetoHook.reset(new PowhegHooks());
386 
387  edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
388  fMultiUserHook->addHook(fEmissionVetoHook.get());
389  }
390 
391  bool PowhegRes = fMasterGen->settings.flag("POWHEGres:calcScales");
392  if (PowhegRes) {
393  edm::LogInfo("Pythia8Interface") << "Turning on resonance scale setting from CMSSW Pythia8Interface";
394  fPowhegResHook.reset(new PowhegResHook());
395  fMultiUserHook->addHook(fPowhegResHook.get());
396  }
397 
398  bool PowhegBB4L = fMasterGen->settings.flag("POWHEG:bb4l");
399  if (PowhegBB4L) {
400  edm::LogInfo("Pythia8Interface") << "Turning on BB4l hook from CMSSW Pythia8Interface";
401  fPowhegHooksBB4L.reset(new PowhegHooksBB4L());
402  fMultiUserHook->addHook(fPowhegHooksBB4L.get());
403  }
404 
405  //adapted from main89.cc in pythia8 examples
406  bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
407  bool internalMerging = !(fMasterGen->settings.word("Merging:Process")=="void");
408 
409  if (internalMatching && internalMerging) {
410  throw edm::Exception(edm::errors::Configuration,"Pythia8Interface")
411  <<" Only one jet matching/merging scheme can be used at a time. \n";
412  }
413 
414  if (internalMatching) {
415  fJetMatchingPy8InternalHook.reset(new Pythia8::JetMatchingMadgraph);
417  }
418 
419  if (internalMerging) {
420  int scheme = ( fMasterGen->settings.flag("Merging:doUMEPSTree")
421  || fMasterGen->settings.flag("Merging:doUMEPSSubt")) ?
422  1 :
423  ( ( fMasterGen->settings.flag("Merging:doUNLOPSTree")
424  || fMasterGen->settings.flag("Merging:doUNLOPSSubt")
425  || fMasterGen->settings.flag("Merging:doUNLOPSLoop")
426  || fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO")) ?
427  2 :
428  0 );
429  fMergingHook.reset(new Pythia8::amcnlo_unitarised_interface(scheme));
430  fMultiUserHook->addHook(fMergingHook.get());
431  }
432 
433  bool resonanceDecayFilter = fMasterGen->settings.flag("ResonanceDecayFilter:filter");
434  if (resonanceDecayFilter) {
437  }
438 
439  bool PTFilter = fMasterGen->settings.flag("PTFilter:filter");
440  if (PTFilter) {
441  fPTFilterHook.reset(new PTFilterHook);
442  fMultiUserHook->addHook(fPTFilterHook.get());
443  }
444 
445  if (fMultiUserHook->nHooks()>0) {
446  fMasterGen->setUserHooksPtr(fMultiUserHook.get());
447  }
448 
449  edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
450  if( fvincia.get() ) {
451  fvincia->init(); status = true;
452  }
453  else if( fDire.get() ) {
454  //fDire->initTune(*fMasterGen.get());
455  fDire->weightsPtr->setup();
456  fMasterGen->init();
457  fDire->setup(*fMasterGen.get());
458  status = true;
459  }
460  else {
461  status = fMasterGen->init();
462  }
463 
464  //clean up temp file
465  if (!slhafile_.empty()) {
466  std::remove(slhafile_.c_str());
467  }
468 
469  if ( pythiaPylistVerbosity > 10 )
470  {
471  if ( pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13 )
472  fMasterGen->settings.listAll();
473  if ( pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13 )
474  fMasterGen->particleData.listAll();
475  }
476 
477  // init decayer
478  fDecayer->settings.flag("ProcessLevel:all", false ); // trick
479  fDecayer->settings.flag("ProcessLevel:resonanceDecays", true );
480  edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
481  status1 = fDecayer->init();
482 
483  if (useEvtGen) {
484  edm::LogInfo("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
485 
486  evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile, evtgenPdlFile));
487 
488  for (unsigned int i=0; i<evtgenUserFiles.size(); i++) {
489  edm::FileInPath evtgenUserFile(evtgenUserFiles.at(i));
490  evtgenDecays->readDecayFile(evtgenUserFile.fullPath());
491  }
492 
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("Pythia8Interface") << "Creating and initializing pythia8 EvtGen plugin";
635 
636  std::string evtgenpath(getenv("EVTGENDATA"));
637  evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile, evtgenPdlFile));
638 
639  for (unsigned int i=0; i<evtgenUserFiles.size(); i++) {
640  edm::FileInPath evtgenUserFile(evtgenUserFiles.at(i));
641  evtgenDecays->readDecayFile(evtgenUserFile.fullPath());
642  }
643 
644  }
645 
646  return (status&&status1);
647 }
648 
649 
651 {
652  fMasterGen->stat();
653 
654  if(fEmissionVetoHook.get()) {
655  edm::LogPrint("Pythia8Interface") << "\n"
656  << "Number of ISR vetoed = " << nISRveto;
657  edm::LogPrint("Pythia8Interface")
658  << "Number of FSR vetoed = " << nFSRveto;
659  }
660 
661  if(fPowhegHooksBB4L.get()) {
662  edm::LogInfo("Pythia8Interface") << "\n"
663  << "BB4L: Number of FSR vetoed = " << nFSRvetoBB4L;
664  }
665 
666  double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
667  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
668  double err = fMasterGen->info.sigmaErr(); // cross section err in mb
669  err *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
671 }
672 
673 
675 {
676 
677  DJR.resize(0);
678  nME = -1;
679  nMEFiltered = -1;
680 
681  if ( fJetMatchingHook.get() )
682  {
683  fJetMatchingHook->resetMatchingStatus();
684  fJetMatchingHook->beforeHadronization( lheEvent() );
685  }
686 
687  if (!fMasterGen->next()) return false;
688 
689  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
690  if (fMergingHook.get()) {
691  mergeweight *= fMergingHook->getNormFactor();
692  }
693 
694  //protect against 0-weight from ckkw or similar
695  if (std::abs(mergeweight)==0.)
696  {
697  event().reset();
698  return false;
699  }
700 
701  if (fJetMatchingPy8InternalHook.get()) {
702  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
703  //cap size of djr vector to save storage space (keep only up to first 6 elements)
704  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
705  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
706  DJR.push_back(djrmatch[idjr]);
707  }
708 
709  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
710  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
711  }
712 
713  if (evtgenDecays.get()) evtgenDecays->decay();
714 
715  event().reset(new HepMC::GenEvent);
716  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
717 
718  if (!py8hepmc) {
719  return false;
720  }
721 
722  //add ckkw/umeps/unlops merging weight
723  if (mergeweight!=1.) {
724  event()->weights()[0] *= mergeweight;
725  }
726 
727  if (fEmissionVetoHook.get()) {
728  nISRveto += fEmissionVetoHook->getNISRveto();
729  nFSRveto += fEmissionVetoHook->getNFSRveto();
730  }
731  if (fPowhegHooksBB4L.get()) {
732  nFSRvetoBB4L += fPowhegHooksBB4L->getNFSRveto();
733  }
734 
735  //fill additional weights for systematic uncertainties
736  if (fMasterGen->info.getWeightsDetailedSize() > 0) {
737  for (const string &key : fMasterGen->info.initrwgt->weightsKeys) {
738  double wgt = (*fMasterGen->info.weights_detailed)[key];
739  event()->weights().push_back(wgt);
740  }
741  }
742  else if (fMasterGen->info.getWeightsCompressedSize() > 0) {
743  for (unsigned int i = 0; i < fMasterGen->info.getWeightsCompressedSize(); i++) {
744  double wgt = fMasterGen->info.getWeightsCompressedValue(i);
745  event()->weights().push_back(wgt);
746  }
747  }
748 
749  // fill shower weights
750  // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
751  if( fMasterGen->info.nWeights() > 1 ){
752  for(int i = 0; i < fMasterGen->info.nWeights(); ++i) {
753  double wgt = fMasterGen->info.weight(i);
754  event()->weights().push_back(wgt);
755  }
756  }
757 
758  // VINCIA shower weights
759  // http://vincia.hepforge.org/current/share/Vincia/htmldoc/VinciaUncertainties.html
760  if( fvincia.get() ) {
761  event()->weights()[0] *= fvincia->weight(0);
762  for (int iVar=1; iVar < fvincia->nWeights(); iVar++) {
763  event()->weights().push_back(fvincia->weight(iVar));
764  }
765  }
766 
767  // Retrieve Dire shower weights
768  if( fDire.get() ) {
769  fDire->weightsPtr->calcWeight(0.);
770  fDire->weightsPtr->reset();
771 
772  //Make sure the base weight comes first
773  event()->weights()[0] *= fDire->weightsPtr->getShowerWeight("base");
774 
775  map<string, double>::iterator it;
776  for ( it = fDire->weightsPtr->getShowerWeights()->begin(); it != fDire->weightsPtr->getShowerWeights()->end(); it++ ) {
777  if (it->first == "base") continue;
778  event()->weights().push_back(it->second);
779  }
780  }
781 
782  return true;
783 
784 }
785 
786 
788 {
789  DJR.resize(0);
790  nME = -1;
791  nMEFiltered = -1;
792  if(LHEInputFileName.empty()) lhaUP->loadEvent(lheEvent());
793 
794  if ( fJetMatchingHook.get() )
795  {
796  fJetMatchingHook->resetMatchingStatus();
797  fJetMatchingHook->beforeHadronization( lheEvent() );
798  }
799 
800  bool py8next = fMasterGen->next();
801 
802  double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
803  if (fMergingHook.get()) {
804  mergeweight *= fMergingHook->getNormFactor();
805  }
806 
807 
808  //protect against 0-weight from ckkw or similar
809  if (!py8next || std::abs(mergeweight)==0.)
810  {
811  lheEvent()->count( lhef::LHERunInfo::kSelected, 1.0, mergeweight );
812  event().reset();
813  return false;
814  }
815 
816  if (fJetMatchingPy8InternalHook.get()) {
817  const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
818  //cap size of djr vector to save storage space (keep only up to first 6 elements)
819  unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
820  for (unsigned int idjr=0; idjr<ndjr; ++idjr) {
821  DJR.push_back(djrmatch[idjr]);
822  }
823 
824  nME=fJetMatchingPy8InternalHook->nMEpartons().first;
825  nMEFiltered=fJetMatchingPy8InternalHook->nMEpartons().second;
826  }
827 
828  // update LHE matching statistics
829  //
830  lheEvent()->count( lhef::LHERunInfo::kAccepted, 1.0, mergeweight );
831 
832  if (evtgenDecays.get()) evtgenDecays->decay();
833 
834  event().reset(new HepMC::GenEvent);
835  bool py8hepmc = toHepMC.fill_next_event( *(fMasterGen.get()), event().get());
836  if (!py8hepmc) {
837  return false;
838  }
839 
840  //add ckkw/umeps/unlops merging weight
841  if (mergeweight!=1.) {
842  event()->weights()[0] *= mergeweight;
843  }
844 
845  if (fEmissionVetoHook.get()) {
846  nISRveto += fEmissionVetoHook->getNISRveto();
847  nFSRveto += fEmissionVetoHook->getNFSRveto();
848  }
849 
850  // fill shower weights
851  // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
852  if( fMasterGen->info.nWeights() > 1 ){
853  for(int i = 0; i < fMasterGen->info.nWeights(); ++i) {
854  double wgt = fMasterGen->info.weight(i);
855  event()->weights().push_back(wgt);
856  }
857  }
858 
859  return true;
860 
861 }
862 
863 
865 {
866 
867  Event* pythiaEvent = &(fMasterGen->event);
868 
869  int NPartsBeforeDecays = pythiaEvent->size();
870  int NPartsAfterDecays = event().get()->particles_size();
871 
872  if(NPartsAfterDecays == NPartsBeforeDecays) return true;
873 
874  bool result = true;
875 
876  for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
877  {
878 
879  HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
880 
881  if ( part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id()) )
882  {
883  fDecayer->event.reset();
884  Particle py8part( part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
885  part->momentum().x(),
886  part->momentum().y(),
887  part->momentum().z(),
888  part->momentum().t(),
889  part->generated_mass() );
890  HepMC::GenVertex* ProdVtx = part->production_vertex();
891  py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(),
892  ProdVtx->position().z(), ProdVtx->position().t() );
893  py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
894  fDecayer->event.append( py8part );
895  int nentries = fDecayer->event.size();
896  if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
897  fDecayer->next();
898  int nentries1 = fDecayer->event.size();
899  if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
900 
901  part->set_status(2);
902 
903  result = toHepMC.fill_next_event( *(fDecayer.get()), event().get(), -1, true, part);
904 
905  }
906  }
907 
908  return result;
909 
910 }
911 
912 
914 {
915  bool lhe = lheEvent() != nullptr;
916 
917  // now create the GenEventInfo product from the GenEvent and fill
918  // the missing pieces
919  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
920 
921  // in pythia pthat is used to subdivide samples into different bins
922  // in LHE mode the binning is done by the external ME generator
923  // which is likely not pthat, so only filling it for Py6 internal mode
924  if (!lhe) {
925  eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
926  }
927 
928  eventInfo()->setDJR(DJR);
929  eventInfo()->setNMEPartons(nME);
930  eventInfo()->setNMEPartonsFiltered(nMEFiltered);
931 
932  //******** Verbosity ********
933 
934  if (maxEventsToPrint > 0 &&
938  if (pythiaPylistVerbosity) {
939  fMasterGen->info.list();
940  fMasterGen->event.list();
941  }
942 
943  if (pythiaHepMCVerbosity) {
944  std::cout << "Event process = "
945  << fMasterGen->info.code() << "\n"
946  << "----------------------" << std::endl;
947  event()->print();
948  }
950  std::cout << "Event process = "
951  << fMasterGen->info.code() << "\n"
952  << "----------------------" << std::endl;
953  ascii_io->write_event(event().get());
954  }
955  }
956 }
957 
959  GenLumiInfoHeader *genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
960 
961  //fill lhe headers
962  //*FIXME* initrwgt header is corrupt due to pythia bug
963  for (const std::string &key : fMasterGen->info.headerKeys()) {
964  genLumiInfoHeader->lheHeaders().emplace_back(key,fMasterGen->info.header(key));
965  }
966 
967  //check, if it is not only nominal weight
968  int weights_number = fMasterGen->info.nWeights();
969  if (fMasterGen->info.initrwgt) weights_number += fMasterGen->info.initrwgt->weightsKeys.size();
970  if(weights_number > 1){
971  genLumiInfoHeader->weightNames().reserve(weights_number + 1);
972  genLumiInfoHeader->weightNames().push_back("nominal");
973  }
974 
975  //fill weight names
976  if (fMasterGen->info.initrwgt) {
977  for (const std::string &key : fMasterGen->info.initrwgt->weightsKeys) {
978  std::string weightgroupname;
979  for (const auto &wgtgrp : fMasterGen->info.initrwgt->weightgroups) {
980  const auto &wgtgrpwgt = wgtgrp.second.weights.find(key);
981  if (wgtgrpwgt != wgtgrp.second.weights.end()) {
982  weightgroupname = wgtgrp.first;
983  }
984  }
985 
986  std::ostringstream weightname;
987  weightname << "LHE, id = " << key << ", ";
988  if (!weightgroupname.empty()) {
989  weightname << "group = " << weightgroupname << ", ";
990  }
991  weightname<< fMasterGen->info.initrwgt->weights[key].contents;
992  genLumiInfoHeader->weightNames().push_back(weightname.str());
993  }
994  }
995 
996  //fill shower labels
997  // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
998  // http://home.thep.lu.se/~torbjorn/doxygen/classPythia8_1_1Info.html
999  if( fMasterGen->info.nWeights() > 1 ){
1000  for(int i = 0; i < fMasterGen->info.nWeights(); ++i) {
1001  genLumiInfoHeader->weightNames().push_back( fMasterGen->info.weightLabel(i) );
1002  }
1003  }
1004 
1005  // VINCIA shower weights
1006  // http://vincia.hepforge.org/current/share/Vincia/htmldoc/VinciaUncertainties.html
1007  if( fvincia.get() ) {
1008  for (int iVar=0; iVar < fvincia->nWeights(); iVar++) {
1009  genLumiInfoHeader->weightNames().push_back( fvincia->weightLabel(iVar) );
1010  }
1011  }
1012 
1013  if( fDire.get() ) {
1014  //Make sure the base weight comes first
1015  genLumiInfoHeader->weightNames().push_back("base");
1016 
1017  map<string, double>::iterator it;
1018  for ( it = fDire->weightsPtr->getShowerWeights()->begin(); it != fDire->weightsPtr->getShowerWeights()->end(); it++ ) {
1019  if (it->first == "base") continue;
1020  genLumiInfoHeader->weightNames().push_back(it->first);
1021  }
1022  }
1023 
1024  return genLumiInfoHeader;
1025 }
1026 
1029 
1030 
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