CMS 3D CMS Logo

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