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