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