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