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