Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-05-24 02:14:50

0001 #include <iostream>
0002 #include <memory>
0003 
0004 #include <cstdint>
0005 #include <memory>
0006 #include <sstream>
0007 #include <string>
0008 #include <vector>
0009 
0010 #include "HepMC/GenEvent.h"
0011 #include "HepMC/GenParticle.h"
0012 
0013 #include "Pythia8/Pythia.h"
0014 #include "Pythia8Plugins/HepMC2.h"
0015 
0016 using namespace Pythia8;
0017 
0018 #include "GeneratorInterface/Pythia8Interface/interface/Py8InterfaceBase.h"
0019 
0020 #include "GeneratorInterface/Pythia8Interface/plugins/ReweightUserHooks.h"
0021 #include "GeneratorInterface/Pythia8Interface/interface/CustomHook.h"
0022 
0023 // PS matchning prototype
0024 //
0025 #include "GeneratorInterface/Pythia8Interface/plugins/JetMatchingHook.h"
0026 #include "Pythia8Plugins/JetMatching.h"
0027 #include "Pythia8Plugins/aMCatNLOHooks.h"
0028 
0029 // Emission Veto Hooks
0030 //
0031 #include "Pythia8Plugins/PowhegHooks.h"
0032 #include "GeneratorInterface/Pythia8Interface/plugins/EmissionVetoHook1.h"
0033 
0034 // Resonance scale hook
0035 #include "GeneratorInterface/Pythia8Interface/plugins/PowhegResHook.h"
0036 #include "GeneratorInterface/Pythia8Interface/plugins/PowhegHooksBB4L.h"
0037 
0038 //biased tau decayer
0039 #include "GeneratorInterface/Pythia8Interface/interface/BiasedTauDecayer.h"
0040 
0041 //decay filter hook
0042 #include "GeneratorInterface/Pythia8Interface/interface/ResonanceDecayFilterHook.h"
0043 
0044 //decay filter hook
0045 #include "GeneratorInterface/Pythia8Interface/interface/PTFilterHook.h"
0046 
0047 // EvtGen plugin
0048 //
0049 #include "Pythia8Plugins/EvtGen.h"
0050 
0051 #include "FWCore/Concurrency/interface/SharedResourceNames.h"
0052 #include "FWCore/ServiceRegistry/interface/Service.h"
0053 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0054 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0055 #include "FWCore/ParameterSet/interface/FileInPath.h"
0056 
0057 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0058 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
0059 
0060 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
0061 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
0062 #include "GeneratorInterface/Core/interface/ConcurrentGeneratorFilter.h"
0063 #include "GeneratorInterface/Core/interface/ConcurrentHadronizerFilter.h"
0064 
0065 #include "GeneratorInterface/Pythia8Interface/plugins/LHAupLesHouches.h"
0066 
0067 #include "HepPID/ParticleIDTranslations.hh"
0068 
0069 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
0070 #include "GeneratorInterface/ExternalDecays/interface/ConcurrentExternalDecayDriver.h"
0071 
0072 namespace CLHEP {
0073   class HepRandomEngine;
0074 }
0075 
0076 using namespace gen;
0077 
0078 class Pythia8Hadronizer : public Py8InterfaceBase {
0079 public:
0080   Pythia8Hadronizer(const edm::ParameterSet &params);
0081   ~Pythia8Hadronizer() override;
0082 
0083   bool initializeForInternalPartons() override;
0084   bool initializeForExternalPartons();
0085 
0086   bool generatePartonsAndHadronize() override;
0087   bool hadronize();
0088 
0089   virtual bool residualDecay();
0090 
0091   void finalizeEvent() override;
0092 
0093   void statistics() override;
0094 
0095   const char *classname() const override { return "Pythia8Hadronizer"; }
0096 
0097   std::unique_ptr<GenLumiInfoHeader> getGenLumiInfoHeader() const override;
0098 
0099 private:
0100   void doSetRandomEngine(CLHEP::HepRandomEngine *v) override { p8SetRandomEngine(v); }
0101   std::vector<std::string> const &doSharedResources() const override { return p8SharedResources; }
0102 
0103   /// Center-of-Mass energy
0104   double comEnergy;
0105 
0106   std::string LHEInputFileName;
0107   std::shared_ptr<LHAupLesHouches> lhaUP;
0108 
0109   enum { PP, PPbar, ElectronPositron };
0110   int fInitialState;  // pp, ppbar, or e-e+
0111 
0112   double fBeam1PZ;
0113   double fBeam2PZ;
0114 
0115   //helper class to allow multiple user hooks simultaneously
0116   std::shared_ptr<UserHooksVector> fUserHooksVector;
0117   bool UserHooksSet;
0118 
0119   // Reweight user hooks
0120   //
0121   std::shared_ptr<UserHooks> fReweightUserHook;
0122   std::shared_ptr<UserHooks> fReweightEmpUserHook;
0123   std::shared_ptr<UserHooks> fReweightRapUserHook;
0124   std::shared_ptr<UserHooks> fReweightPtHatRapUserHook;
0125 
0126   // PS matching prototype
0127   //
0128   std::shared_ptr<JetMatchingHook> fJetMatchingHook;
0129   std::shared_ptr<Pythia8::JetMatchingMadgraph> fJetMatchingPy8InternalHook;
0130   std::shared_ptr<Pythia8::amcnlo_unitarised_interface> fMergingHook;
0131 
0132   // Emission Veto Hooks
0133   //
0134   std::shared_ptr<PowhegHooks> fEmissionVetoHook;
0135   std::shared_ptr<EmissionVetoHook1> fEmissionVetoHook1;
0136 
0137   // Resonance scale hook
0138   std::shared_ptr<PowhegResHook> fPowhegResHook;
0139   std::shared_ptr<PowhegHooksBB4L> fPowhegHooksBB4L;
0140 
0141   // biased tau decayer
0142   std::shared_ptr<BiasedTauDecayer> fBiasedTauDecayer;
0143 
0144   //resonance decay filter hook
0145   std::shared_ptr<ResonanceDecayFilterHook> fResonanceDecayFilterHook;
0146 
0147   //PT filter hook
0148   std::shared_ptr<PTFilterHook> fPTFilterHook;
0149 
0150   //Generic customized hooks vector
0151   std::shared_ptr<UserHooksVector> fCustomHooksVector;
0152 
0153   int EV1_nFinal;
0154   bool EV1_vetoOn;
0155   int EV1_maxVetoCount;
0156   int EV1_pThardMode;
0157   int EV1_pTempMode;
0158   int EV1_emittedMode;
0159   int EV1_pTdefMode;
0160   bool EV1_MPIvetoOn;
0161   int EV1_QEDvetoMode;
0162   int EV1_nFinalMode;
0163 
0164   static const std::vector<std::string> p8SharedResources;
0165 
0166   vector<float> DJR;
0167   int nME;
0168   int nMEFiltered;
0169 
0170   int nISRveto;
0171   int nFSRveto;
0172 };
0173 
0174 const std::vector<std::string> Pythia8Hadronizer::p8SharedResources = {edm::SharedResourceNames::kPythia8};
0175 
0176 Pythia8Hadronizer::Pythia8Hadronizer(const edm::ParameterSet &params)
0177     : Py8InterfaceBase(params),
0178       comEnergy(params.getParameter<double>("comEnergy")),
0179       LHEInputFileName(params.getUntrackedParameter<std::string>("LHEInputFileName", "")),
0180       fInitialState(PP),
0181       UserHooksSet(false),
0182       nME(-1),
0183       nMEFiltered(-1),
0184       nISRveto(0),
0185       nFSRveto(0) {
0186   // J.Y.: the following 3 parameters are hacked "for a reason"
0187   //
0188   if (params.exists("PPbarInitialState")) {
0189     if (fInitialState == PP) {
0190       fInitialState = PPbar;
0191       edm::LogImportant("GeneratorInterface|Pythia8Interface")
0192           << "Pythia8 will be initialized for PROTON-ANTIPROTON INITIAL STATE. "
0193           << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
0194     } else {
0195       // probably need to throw on attempt to override ?
0196     }
0197   } else if (params.exists("ElectronPositronInitialState")) {
0198     if (fInitialState == PP) {
0199       fInitialState = ElectronPositron;
0200       edm::LogInfo("GeneratorInterface|Pythia8Interface")
0201           << "Pythia8 will be initialized for ELECTRON-POSITRON INITIAL STATE. "
0202           << "This is a user-request change from the DEFAULT PROTON-PROTON initial state.";
0203     } else {
0204       // probably need to throw on attempt to override ?
0205     }
0206   } else if (params.exists("ElectronProtonInitialState") || params.exists("PositronProtonInitialState")) {
0207     // throw on unknown initial state !
0208     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0209         << " UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
0210   }
0211 
0212   // avoid filling weights twice (from v8.30x)
0213   toHepMC.set_store_weights(false);
0214 
0215   // Reweight user hook
0216   //
0217   if (params.exists("reweightGen")) {
0218     edm::LogInfo("Pythia8Interface") << "Start setup for reweightGen";
0219     edm::ParameterSet rgParams = params.getParameter<edm::ParameterSet>("reweightGen");
0220     fReweightUserHook.reset(
0221         new PtHatReweightUserHook(rgParams.getParameter<double>("pTRef"), rgParams.getParameter<double>("power")));
0222     edm::LogInfo("Pythia8Interface") << "End setup for reweightGen";
0223   }
0224   if (params.exists("reweightGenEmp")) {
0225     edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenEmp";
0226     edm::ParameterSet rgeParams = params.getParameter<edm::ParameterSet>("reweightGenEmp");
0227 
0228     std::string tuneName = "";
0229     if (rgeParams.exists("tune"))
0230       tuneName = rgeParams.getParameter<std::string>("tune");
0231     fReweightEmpUserHook.reset(new PtHatEmpReweightUserHook(tuneName));
0232     edm::LogInfo("Pythia8Interface") << "End setup for reweightGenEmp";
0233   }
0234   if (params.exists("reweightGenRap")) {
0235     edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenRap";
0236     edm::ParameterSet rgrParams = params.getParameter<edm::ParameterSet>("reweightGenRap");
0237     fReweightRapUserHook.reset(new RapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
0238                                                        rgrParams.getParameter<double>("yLabPower"),
0239                                                        rgrParams.getParameter<std::string>("yCMSigmaFunc"),
0240                                                        rgrParams.getParameter<double>("yCMPower"),
0241                                                        rgrParams.getParameter<double>("pTHatMin"),
0242                                                        rgrParams.getParameter<double>("pTHatMax")));
0243     edm::LogInfo("Pythia8Interface") << "End setup for reweightGenRap";
0244   }
0245   if (params.exists("reweightGenPtHatRap")) {
0246     edm::LogInfo("Pythia8Interface") << "Start setup for reweightGenPtHatRap";
0247     edm::ParameterSet rgrParams = params.getParameter<edm::ParameterSet>("reweightGenPtHatRap");
0248     fReweightPtHatRapUserHook.reset(new PtHatRapReweightUserHook(rgrParams.getParameter<std::string>("yLabSigmaFunc"),
0249                                                                  rgrParams.getParameter<double>("yLabPower"),
0250                                                                  rgrParams.getParameter<std::string>("yCMSigmaFunc"),
0251                                                                  rgrParams.getParameter<double>("yCMPower"),
0252                                                                  rgrParams.getParameter<double>("pTHatMin"),
0253                                                                  rgrParams.getParameter<double>("pTHatMax")));
0254     edm::LogInfo("Pythia8Interface") << "End setup for reweightGenPtHatRap";
0255   }
0256 
0257   if (params.exists("useUserHook"))
0258     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0259         << " Obsolete parameter: useUserHook \n Please use the actual one instead \n";
0260 
0261   // PS matching prototype
0262   //
0263   if (params.exists("jetMatching")) {
0264     edm::ParameterSet jmParams = params.getUntrackedParameter<edm::ParameterSet>("jetMatching");
0265     std::string scheme = jmParams.getParameter<std::string>("scheme");
0266     if (scheme == "Madgraph" || scheme == "MadgraphFastJet") {
0267       fJetMatchingHook.reset(new JetMatchingHook(jmParams, &fMasterGen->info));
0268     }
0269   }
0270 
0271   // Pythia8Interface emission veto
0272   //
0273   if (params.exists("emissionVeto1")) {
0274     EV1_nFinal = -1;
0275     if (params.exists("EV1_nFinal"))
0276       EV1_nFinal = params.getParameter<int>("EV1_nFinal");
0277     EV1_vetoOn = true;
0278     if (params.exists("EV1_vetoOn"))
0279       EV1_vetoOn = params.getParameter<bool>("EV1_vetoOn");
0280     EV1_maxVetoCount = 10;
0281     if (params.exists("EV1_maxVetoCount"))
0282       EV1_maxVetoCount = params.getParameter<int>("EV1_maxVetoCount");
0283     EV1_pThardMode = 1;
0284     if (params.exists("EV1_pThardMode"))
0285       EV1_pThardMode = params.getParameter<int>("EV1_pThardMode");
0286     EV1_pTempMode = 0;
0287     if (params.exists("EV1_pTempMode"))
0288       EV1_pTempMode = params.getParameter<int>("EV1_pTempMode");
0289     if (EV1_pTempMode > 2 || EV1_pTempMode < 0)
0290       throw edm::Exception(edm::errors::Configuration, "Pythia8Interface") << " Wrong value for EV1_pTempMode code\n";
0291     EV1_emittedMode = 0;
0292     if (params.exists("EV1_emittedMode"))
0293       EV1_emittedMode = params.getParameter<int>("EV1_emittedMode");
0294     EV1_pTdefMode = 1;
0295     if (params.exists("EV1_pTdefMode"))
0296       EV1_pTdefMode = params.getParameter<int>("EV1_pTdefMode");
0297     EV1_MPIvetoOn = false;
0298     if (params.exists("EV1_MPIvetoOn"))
0299       EV1_MPIvetoOn = params.getParameter<bool>("EV1_MPIvetoOn");
0300     EV1_QEDvetoMode = 0;
0301     if (params.exists("EV1_QEDvetoMode"))
0302       EV1_QEDvetoMode = params.getParameter<int>("EV1_QEDvetoMode");
0303     EV1_nFinalMode = 0;
0304     if (params.exists("EV1_nFinalMode"))
0305       EV1_nFinalMode = params.getParameter<int>("EV1_nFinalMode");
0306     fEmissionVetoHook1.reset(new EmissionVetoHook1(EV1_nFinal,
0307                                                    EV1_vetoOn,
0308                                                    EV1_maxVetoCount,
0309                                                    EV1_pThardMode,
0310                                                    EV1_pTempMode,
0311                                                    EV1_emittedMode,
0312                                                    EV1_pTdefMode,
0313                                                    EV1_MPIvetoOn,
0314                                                    EV1_QEDvetoMode,
0315                                                    EV1_nFinalMode,
0316                                                    0));
0317   }
0318 
0319   if (params.exists("UserCustomization")) {
0320     fCustomHooksVector = std::make_shared<UserHooksVector>();
0321     const std::vector<edm::ParameterSet> userParams =
0322         params.getParameter<std::vector<edm::ParameterSet>>("UserCustomization");
0323     for (const auto &pluginParams : userParams) {
0324       (fCustomHooksVector->hooks)
0325           .push_back(
0326               CustomHookFactory::get()->create(pluginParams.getParameter<std::string>("pluginName"), pluginParams));
0327     }
0328   }
0329 
0330   if (params.exists("VinciaPlugin")) {
0331     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0332         << " Obsolete parameter: VinciaPlugin \n Please use the parameter PartonShowers:model instead \n";
0333   }
0334   if (params.exists("DirePlugin")) {
0335     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0336         << " Obsolete parameter: DirePlugin \n Please use the parameter PartonShowers:model instead \n";
0337   }
0338 }
0339 
0340 Pythia8Hadronizer::~Pythia8Hadronizer() {}
0341 
0342 bool Pythia8Hadronizer::initializeForInternalPartons() {
0343   bool status = false, status1 = false;
0344 
0345   if (lheFile_.empty()) {
0346     if (fInitialState == PP)  // default
0347     {
0348       fMasterGen->settings.mode("Beams:idA", 2212);
0349       fMasterGen->settings.mode("Beams:idB", 2212);
0350     } else if (fInitialState == PPbar) {
0351       fMasterGen->settings.mode("Beams:idA", 2212);
0352       fMasterGen->settings.mode("Beams:idB", -2212);
0353     } else if (fInitialState == ElectronPositron) {
0354       fMasterGen->settings.mode("Beams:idA", 11);
0355       fMasterGen->settings.mode("Beams:idB", -11);
0356     } else {
0357       // throw on unknown initial state !
0358       throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0359           << " UNKNOWN INITIAL STATE. \n The allowed initial states are: PP, PPbar, ElectronPositron \n";
0360     }
0361     fMasterGen->settings.parm("Beams:eCM", comEnergy);
0362   } else {
0363     fMasterGen->settings.mode("Beams:frameType", 4);
0364     fMasterGen->settings.word("Beams:LHEF", lheFile_);
0365   }
0366 
0367   if (!fUserHooksVector.get())
0368     fUserHooksVector.reset(new UserHooksVector);
0369   (fUserHooksVector->hooks).clear();
0370 
0371   if (fReweightUserHook.get())
0372     (fUserHooksVector->hooks).push_back(fReweightUserHook);
0373   if (fReweightEmpUserHook.get())
0374     (fUserHooksVector->hooks).push_back(fReweightEmpUserHook);
0375   if (fReweightRapUserHook.get())
0376     (fUserHooksVector->hooks).push_back(fReweightRapUserHook);
0377   if (fReweightPtHatRapUserHook.get())
0378     (fUserHooksVector->hooks).push_back(fReweightPtHatRapUserHook);
0379   if (fJetMatchingHook.get())
0380     (fUserHooksVector->hooks).push_back(fJetMatchingHook);
0381   if (fEmissionVetoHook1.get()) {
0382     edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
0383     (fUserHooksVector->hooks).push_back(fEmissionVetoHook1);
0384   }
0385 
0386   if (fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) {
0387     if (fJetMatchingHook.get() || fEmissionVetoHook1.get())
0388       throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0389           << " Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible "
0390              "are : jetMatching, emissionVeto1 \n";
0391 
0392     if (!fEmissionVetoHook.get())
0393       fEmissionVetoHook.reset(new PowhegHooks());
0394 
0395     edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
0396     (fUserHooksVector->hooks).push_back(fEmissionVetoHook);
0397   }
0398 
0399   bool PowhegRes = fMasterGen->settings.flag("POWHEGres:calcScales");
0400   if (PowhegRes) {
0401     edm::LogInfo("Pythia8Interface") << "Turning on resonance scale setting from CMSSW Pythia8Interface";
0402     if (!fPowhegResHook.get())
0403       fPowhegResHook.reset(new PowhegResHook());
0404     (fUserHooksVector->hooks).push_back(fPowhegResHook);
0405   }
0406 
0407   bool PowhegBB4L = fMasterGen->settings.flag("POWHEG:bb4l");
0408   if (PowhegBB4L) {
0409     edm::LogInfo("Pythia8Interface") << "Turning on BB4l hook from CMSSW Pythia8Interface";
0410     if (!fPowhegHooksBB4L.get())
0411       fPowhegHooksBB4L.reset(new PowhegHooksBB4L());
0412     (fUserHooksVector->hooks).push_back(fPowhegHooksBB4L);
0413   }
0414 
0415   //adapted from main89.cc in pythia8 examples
0416   bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
0417   bool internalMerging = !(fMasterGen->settings.word("Merging:Process") == "void");
0418 
0419   if (internalMatching && internalMerging) {
0420     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0421         << " Only one jet matching/merging scheme can be used at a time. \n";
0422   }
0423 
0424   if (internalMatching) {
0425     if (!fJetMatchingPy8InternalHook.get())
0426       fJetMatchingPy8InternalHook.reset(new Pythia8::JetMatchingMadgraph);
0427     (fUserHooksVector->hooks).push_back(fJetMatchingPy8InternalHook);
0428   }
0429 
0430   if (internalMerging) {
0431     int scheme = (fMasterGen->settings.flag("Merging:doUMEPSTree") || fMasterGen->settings.flag("Merging:doUMEPSSubt"))
0432                      ? 1
0433                      : ((fMasterGen->settings.flag("Merging:doUNLOPSTree") ||
0434                          fMasterGen->settings.flag("Merging:doUNLOPSSubt") ||
0435                          fMasterGen->settings.flag("Merging:doUNLOPSLoop") ||
0436                          fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO"))
0437                             ? 2
0438                             : 0);
0439     if (!fMergingHook.get())
0440       fMergingHook.reset(new Pythia8::amcnlo_unitarised_interface(scheme));
0441     (fUserHooksVector->hooks).push_back(fMergingHook);
0442   }
0443 
0444   bool biasedTauDecayer = fMasterGen->settings.flag("BiasedTauDecayer:filter");
0445   if (biasedTauDecayer) {
0446     if (!fBiasedTauDecayer.get()) {
0447       Pythia8::Info localInfo = fMasterGen->info;
0448       fBiasedTauDecayer.reset(new BiasedTauDecayer(&localInfo, &(fMasterGen->settings)));
0449     }
0450     std::vector<int> handledParticles;
0451     handledParticles.push_back(15);
0452     fMasterGen->setDecayPtr(fBiasedTauDecayer, handledParticles);
0453   }
0454 
0455   bool resonanceDecayFilter = fMasterGen->settings.flag("ResonanceDecayFilter:filter");
0456   if (resonanceDecayFilter) {
0457     fResonanceDecayFilterHook.reset(new ResonanceDecayFilterHook);
0458     (fUserHooksVector->hooks).push_back(fResonanceDecayFilterHook);
0459   }
0460 
0461   bool PTFilter = fMasterGen->settings.flag("PTFilter:filter");
0462   if (PTFilter) {
0463     fPTFilterHook.reset(new PTFilterHook);
0464     (fUserHooksVector->hooks).push_back(fPTFilterHook);
0465   }
0466 
0467   if (!(fUserHooksVector->hooks).empty() && !UserHooksSet) {
0468     for (auto &fUserHook : fUserHooksVector->hooks) {
0469       fMasterGen->addUserHooksPtr(fUserHook);
0470     }
0471     UserHooksSet = true;
0472   }
0473 
0474   if (fCustomHooksVector.get()) {
0475     edm::LogInfo("Pythia8Interface") << "Adding customized user hooks";
0476     for (const auto &fUserHook : fCustomHooksVector->hooks) {
0477       fMasterGen->addUserHooksPtr(fUserHook);
0478     }
0479   }
0480 
0481   edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
0482   status = fMasterGen->init();
0483 
0484   //clean up temp file
0485   if (!slhafile_.empty()) {
0486     std::remove(slhafile_.c_str());
0487   }
0488 
0489   if (pythiaPylistVerbosity > 10) {
0490     if (pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13)
0491       fMasterGen->settings.listAll();
0492     if (pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13)
0493       fMasterGen->particleData.listAll();
0494   }
0495 
0496   // init decayer
0497   fDecayer->settings.flag("ProcessLevel:all", false);  // trick
0498   fDecayer->settings.flag("ProcessLevel:resonanceDecays", true);
0499   edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
0500   status1 = fDecayer->init();
0501 
0502   if (useEvtGen) {
0503     edm::LogInfo("Pythia8Hadronizer") << "Creating and initializing pythia8 EvtGen plugin";
0504     if (!evtgenDecays.get()) {
0505       evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile, evtgenPdlFile));
0506       for (unsigned int i = 0; i < evtgenUserFiles.size(); i++)
0507         evtgenDecays->readDecayFile(evtgenUserFiles.at(i));
0508     }
0509   }
0510 
0511   return (status && status1);
0512 }
0513 
0514 bool Pythia8Hadronizer::initializeForExternalPartons() {
0515   edm::LogInfo("Pythia8Interface") << "Initializing for external partons";
0516 
0517   bool status = false, status1 = false;
0518 
0519   if (!fUserHooksVector.get())
0520     fUserHooksVector.reset(new UserHooksVector);
0521   (fUserHooksVector->hooks).clear();
0522 
0523   if (fReweightUserHook.get())
0524     (fUserHooksVector->hooks).push_back(fReweightUserHook);
0525   if (fReweightEmpUserHook.get())
0526     (fUserHooksVector->hooks).push_back(fReweightEmpUserHook);
0527   if (fReweightRapUserHook.get())
0528     (fUserHooksVector->hooks).push_back(fReweightRapUserHook);
0529   if (fReweightPtHatRapUserHook.get())
0530     (fUserHooksVector->hooks).push_back(fReweightPtHatRapUserHook);
0531   if (fJetMatchingHook.get())
0532     (fUserHooksVector->hooks).push_back(fJetMatchingHook);
0533   if (fEmissionVetoHook1.get()) {
0534     edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook 1 from CMSSW Pythia8Interface";
0535     (fUserHooksVector->hooks).push_back(fEmissionVetoHook1);
0536   }
0537 
0538   if (fCustomHooksVector.get()) {
0539     edm::LogInfo("Pythia8Interface") << "Adding customized user hook";
0540     for (const auto &fUserHook : fCustomHooksVector->hooks) {
0541       (fUserHooksVector->hooks).push_back(fUserHook);
0542     }
0543   }
0544 
0545   if (fMasterGen->settings.mode("POWHEG:veto") > 0 || fMasterGen->settings.mode("POWHEG:MPIveto") > 0) {
0546     if (fJetMatchingHook.get() || fEmissionVetoHook1.get())
0547       throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0548           << " Attempt to turn on PowhegHooks by pythia8 settings but there are incompatible hooks on \n Incompatible "
0549              "are : jetMatching, emissionVeto1 \n";
0550 
0551     if (!fEmissionVetoHook.get())
0552       fEmissionVetoHook.reset(new PowhegHooks());
0553 
0554     edm::LogInfo("Pythia8Interface") << "Turning on Emission Veto Hook from pythia8 code";
0555     (fUserHooksVector->hooks).push_back(fEmissionVetoHook);
0556   }
0557 
0558   bool PowhegRes = fMasterGen->settings.flag("POWHEGres:calcScales");
0559   if (PowhegRes) {
0560     edm::LogInfo("Pythia8Interface") << "Turning on resonance scale setting from CMSSW Pythia8Interface";
0561     if (!fPowhegResHook.get())
0562       fPowhegResHook.reset(new PowhegResHook());
0563     (fUserHooksVector->hooks).push_back(fPowhegResHook);
0564   }
0565 
0566   bool PowhegBB4L = fMasterGen->settings.flag("POWHEG:bb4l");
0567   if (PowhegBB4L) {
0568     edm::LogInfo("Pythia8Interface") << "Turning on BB4l hook from CMSSW Pythia8Interface";
0569     if (!fPowhegHooksBB4L.get())
0570       fPowhegHooksBB4L.reset(new PowhegHooksBB4L());
0571     (fUserHooksVector->hooks).push_back(fPowhegHooksBB4L);
0572   }
0573 
0574   //adapted from main89.cc in pythia8 examples
0575   bool internalMatching = fMasterGen->settings.flag("JetMatching:merge");
0576   bool internalMerging = !(fMasterGen->settings.word("Merging:Process") == "void");
0577 
0578   if (internalMatching && internalMerging) {
0579     throw edm::Exception(edm::errors::Configuration, "Pythia8Interface")
0580         << " Only one jet matching/merging scheme can be used at a time. \n";
0581   }
0582 
0583   if (internalMatching) {
0584     if (!fJetMatchingPy8InternalHook.get())
0585       fJetMatchingPy8InternalHook.reset(new Pythia8::JetMatchingMadgraph);
0586     (fUserHooksVector->hooks).push_back(fJetMatchingPy8InternalHook);
0587   }
0588 
0589   if (internalMerging) {
0590     int scheme = (fMasterGen->settings.flag("Merging:doUMEPSTree") || fMasterGen->settings.flag("Merging:doUMEPSSubt"))
0591                      ? 1
0592                      : ((fMasterGen->settings.flag("Merging:doUNLOPSTree") ||
0593                          fMasterGen->settings.flag("Merging:doUNLOPSSubt") ||
0594                          fMasterGen->settings.flag("Merging:doUNLOPSLoop") ||
0595                          fMasterGen->settings.flag("Merging:doUNLOPSSubtNLO"))
0596                             ? 2
0597                             : 0);
0598     if (!fMergingHook.get())
0599       fMergingHook.reset(new Pythia8::amcnlo_unitarised_interface(scheme));
0600     (fUserHooksVector->hooks).push_back(fMergingHook);
0601   }
0602 
0603   bool biasedTauDecayer = fMasterGen->settings.flag("BiasedTauDecayer:filter");
0604   if (biasedTauDecayer) {
0605     if (!fBiasedTauDecayer.get()) {
0606       Pythia8::Info localInfo = fMasterGen->info;
0607       fBiasedTauDecayer.reset(new BiasedTauDecayer(&localInfo, &(fMasterGen->settings)));
0608     }
0609     std::vector<int> handledParticles;
0610     handledParticles.push_back(15);
0611     fMasterGen->setDecayPtr(fBiasedTauDecayer, handledParticles);
0612   }
0613 
0614   bool resonanceDecayFilter = fMasterGen->settings.flag("ResonanceDecayFilter:filter");
0615   if (resonanceDecayFilter) {
0616     fResonanceDecayFilterHook.reset(new ResonanceDecayFilterHook);
0617     (fUserHooksVector->hooks).push_back(fResonanceDecayFilterHook);
0618   }
0619 
0620   bool PTFilter = fMasterGen->settings.flag("PTFilter:filter");
0621   if (PTFilter) {
0622     fPTFilterHook.reset(new PTFilterHook);
0623     (fUserHooksVector->hooks).push_back(fPTFilterHook);
0624   }
0625 
0626   if (!(fUserHooksVector->hooks).empty() && !UserHooksSet) {
0627     for (auto &fUserHook : fUserHooksVector->hooks) {
0628       fMasterGen->addUserHooksPtr(fUserHook);
0629     }
0630     UserHooksSet = true;
0631   }
0632 
0633   if (!LHEInputFileName.empty()) {
0634     edm::LogInfo("Pythia8Interface") << "Initialize direct pythia8 reading from LHE file " << LHEInputFileName;
0635     edm::LogInfo("Pythia8Interface") << "Some LHE information can be not stored";
0636     fMasterGen->settings.mode("Beams:frameType", 4);
0637     fMasterGen->settings.word("Beams:LHEF", LHEInputFileName);
0638     status = fMasterGen->init();
0639 
0640   } else {
0641     lhaUP.reset(new LHAupLesHouches());
0642     lhaUP->setScalesFromLHEF(fMasterGen->settings.flag("Beams:setProductionScalesFromLHEF"));
0643     lhaUP->loadRunInfo(lheRunInfo());
0644 
0645     if (fJetMatchingHook.get()) {
0646       fJetMatchingHook->init(lheRunInfo());
0647     }
0648 
0649     fMasterGen->settings.mode("Beams:frameType", 5);
0650     fMasterGen->setLHAupPtr(lhaUP);
0651     edm::LogInfo("Pythia8Interface") << "Initializing MasterGen";
0652     status = fMasterGen->init();
0653   }
0654 
0655   //clean up temp file
0656   if (!slhafile_.empty()) {
0657     std::remove(slhafile_.c_str());
0658   }
0659 
0660   if (pythiaPylistVerbosity > 10) {
0661     if (pythiaPylistVerbosity == 11 || pythiaPylistVerbosity == 13)
0662       fMasterGen->settings.listAll();
0663     if (pythiaPylistVerbosity == 12 || pythiaPylistVerbosity == 13)
0664       fMasterGen->particleData.listAll();
0665   }
0666 
0667   // init decayer
0668   fDecayer->settings.flag("ProcessLevel:all", false);  // trick
0669   fDecayer->settings.flag("ProcessLevel:resonanceDecays", true);
0670   edm::LogInfo("Pythia8Interface") << "Initializing Decayer";
0671   status1 = fDecayer->init();
0672 
0673   if (useEvtGen) {
0674     edm::LogInfo("Pythia8Hadronizer") << "Creating and initializing pythia8 EvtGen plugin";
0675     if (!evtgenDecays.get()) {
0676       evtgenDecays.reset(new EvtGenDecays(fMasterGen.get(), evtgenDecFile, evtgenPdlFile));
0677       for (unsigned int i = 0; i < evtgenUserFiles.size(); i++)
0678         evtgenDecays->readDecayFile(evtgenUserFiles.at(i));
0679     }
0680   }
0681 
0682   return (status && status1);
0683 }
0684 
0685 void Pythia8Hadronizer::statistics() {
0686   fMasterGen->stat();
0687 
0688   if (fEmissionVetoHook.get()) {
0689     edm::LogPrint("Pythia8Interface") << "\n"
0690                                       << "Number of ISR vetoed = " << nISRveto;
0691     edm::LogPrint("Pythia8Interface") << "Number of FSR vetoed = " << nFSRveto;
0692   }
0693 
0694   double xsec = fMasterGen->info.sigmaGen();  // cross section in mb
0695   xsec *= 1.0e9;                              // translate to pb (CMS/Gen "convention" as of May 2009)
0696   double err = fMasterGen->info.sigmaErr();   // cross section err in mb
0697   err *= 1.0e9;                               // translate to pb (CMS/Gen "convention" as of May 2009)
0698   runInfo().setInternalXSec(GenRunInfoProduct::XSec(xsec, err));
0699 }
0700 
0701 bool Pythia8Hadronizer::generatePartonsAndHadronize() {
0702   DJR.resize(0);
0703   nME = -1;
0704   nMEFiltered = -1;
0705 
0706   if (fJetMatchingHook.get()) {
0707     fJetMatchingHook->resetMatchingStatus();
0708     fJetMatchingHook->beforeHadronization(lheEvent());
0709   }
0710 
0711   if (!fMasterGen->next())
0712     return false;
0713 
0714   double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
0715   if (fMergingHook.get()) {
0716     mergeweight *= fMergingHook->getNormFactor();
0717   }
0718 
0719   //protect against 0-weight from ckkw or similar
0720   if (std::abs(mergeweight) == 0.) {
0721     event().reset();
0722     return false;
0723   }
0724 
0725   if (fJetMatchingPy8InternalHook.get()) {
0726     const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
0727     //cap size of djr vector to save storage space (keep only up to first 6 elements)
0728     unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
0729     for (unsigned int idjr = 0; idjr < ndjr; ++idjr) {
0730       DJR.push_back(djrmatch[idjr]);
0731     }
0732 
0733     nME = fJetMatchingPy8InternalHook->nMEpartons().first;
0734     nMEFiltered = fJetMatchingPy8InternalHook->nMEpartons().second;
0735   }
0736 
0737   if (evtgenDecays.get())
0738     evtgenDecays->decay();
0739 
0740   event() = std::make_unique<HepMC::GenEvent>();
0741   bool py8hepmc = toHepMC.fill_next_event(*(fMasterGen.get()), event().get());
0742 
0743   if (!py8hepmc) {
0744     return false;
0745   }
0746 
0747   // 0th weight is not filled anymore since v8.30x, pushback manually
0748   event()->weights().push_back(fMasterGen->info.weight());
0749 
0750   //add ckkw/umeps/unlops merging weight
0751   if (mergeweight != 1.) {
0752     event()->weights()[0] *= mergeweight;
0753   }
0754 
0755   if (fEmissionVetoHook.get()) {
0756     nISRveto += fEmissionVetoHook->getNISRveto();
0757     nFSRveto += fEmissionVetoHook->getNFSRveto();
0758   }
0759 
0760   //fill additional weights for systematic uncertainties
0761   if (fMasterGen->info.getWeightsDetailedSize() > 0) {
0762     for (const string &key : fMasterGen->info.initrwgt->weightsKeys) {
0763       double wgt = (*fMasterGen->info.weights_detailed)[key];
0764       event()->weights().push_back(wgt);
0765     }
0766   } else if (fMasterGen->info.getWeightsCompressedSize() > 0) {
0767     for (unsigned int i = 0; i < fMasterGen->info.getWeightsCompressedSize(); i++) {
0768       double wgt = fMasterGen->info.getWeightsCompressedValue(i);
0769       event()->weights().push_back(wgt);
0770     }
0771   }
0772 
0773   // fill shower weights
0774   // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
0775   if (fMasterGen->info.nWeights() > 1) {
0776     for (int i = 0; i < fMasterGen->info.nWeights(); ++i) {
0777       double wgt = fMasterGen->info.weight(i);
0778       event()->weights().push_back(wgt);
0779     }
0780   }
0781 
0782 #if 0
0783   // VINCIA shower weights
0784   // http://vincia.hepforge.org/current/share/Vincia/htmldoc/VinciaUncertainties.html
0785   if (fvincia.get()) {
0786     event()->weights()[0] *= fvincia->weight(0);
0787     for (int iVar = 1; iVar < fvincia->nWeights(); iVar++) {
0788       event()->weights().push_back(fvincia->weight(iVar));
0789     }
0790   }
0791 
0792   // Retrieve Dire shower weights
0793   if (fDire.get()) {
0794     fDire->weightsPtr->calcWeight(0.);
0795     fDire->weightsPtr->reset();
0796 
0797     //Make sure the base weight comes first
0798     event()->weights()[0] *= fDire->weightsPtr->getShowerWeight("base");
0799 
0800     unordered_map<string, double>::iterator it;
0801     for (it = fDire->weightsPtr->getShowerWeights()->begin(); it != fDire->weightsPtr->getShowerWeights()->end();
0802          it++) {
0803       if (it->first == "base")
0804         continue;
0805       event()->weights().push_back(it->second);
0806     }
0807   }
0808 #endif
0809 
0810   return true;
0811 }
0812 
0813 bool Pythia8Hadronizer::hadronize() {
0814   DJR.resize(0);
0815   nME = -1;
0816   nMEFiltered = -1;
0817   if (LHEInputFileName.empty())
0818     lhaUP->loadEvent(lheEvent());
0819 
0820   if (fJetMatchingHook.get()) {
0821     fJetMatchingHook->resetMatchingStatus();
0822     fJetMatchingHook->beforeHadronization(lheEvent());
0823   }
0824 
0825   bool py8next = fMasterGen->next();
0826 
0827   double mergeweight = fMasterGen.get()->info.mergingWeightNLO();
0828   if (fMergingHook.get()) {
0829     mergeweight *= fMergingHook->getNormFactor();
0830   }
0831 
0832   //protect against 0-weight from ckkw or similar
0833   if (!py8next || std::abs(mergeweight) == 0.) {
0834     lheEvent()->count(lhef::LHERunInfo::kSelected, 1.0, mergeweight);
0835     event().reset();
0836     return false;
0837   }
0838 
0839   if (fJetMatchingPy8InternalHook.get()) {
0840     const std::vector<double> djrmatch = fJetMatchingPy8InternalHook->getDJR();
0841     //cap size of djr vector to save storage space (keep only up to first 6 elements)
0842     unsigned int ndjr = std::min(djrmatch.size(), std::vector<double>::size_type(6));
0843     for (unsigned int idjr = 0; idjr < ndjr; ++idjr) {
0844       DJR.push_back(djrmatch[idjr]);
0845     }
0846 
0847     nME = fJetMatchingPy8InternalHook->nMEpartons().first;
0848     nMEFiltered = fJetMatchingPy8InternalHook->nMEpartons().second;
0849   }
0850 
0851   // update LHE matching statistics
0852   //
0853   lheEvent()->count(lhef::LHERunInfo::kAccepted, 1.0, mergeweight);
0854 
0855   if (evtgenDecays.get())
0856     evtgenDecays->decay();
0857 
0858   event() = std::make_unique<HepMC::GenEvent>();
0859 
0860   bool py8hepmc = toHepMC.fill_next_event(*(fMasterGen.get()), event().get());
0861 
0862   if (!py8hepmc) {
0863     return false;
0864   }
0865 
0866   // 0th weight is not filled anymore since v8.30x, pushback manually
0867   event()->weights().push_back(fMasterGen->info.weight());
0868 
0869   //add ckkw/umeps/unlops merging weight
0870   if (mergeweight != 1.) {
0871     event()->weights()[0] *= mergeweight;
0872   }
0873 
0874   if (fEmissionVetoHook.get()) {
0875     nISRveto += fEmissionVetoHook->getNISRveto();
0876     nFSRveto += fEmissionVetoHook->getNFSRveto();
0877   }
0878 
0879   // fill shower weights
0880   // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
0881   if (fMasterGen->info.nWeights() > 1) {
0882     for (int i = 0; i < fMasterGen->info.nWeights(); ++i) {
0883       double wgt = fMasterGen->info.weight(i);
0884       event()->weights().push_back(wgt);
0885     }
0886   }
0887 
0888   return true;
0889 }
0890 
0891 bool Pythia8Hadronizer::residualDecay() {
0892   Event *pythiaEvent = &(fMasterGen->event);
0893 
0894   int NPartsBeforeDecays = pythiaEvent->size();
0895   int NPartsAfterDecays = event().get()->particles_size();
0896 
0897   if (NPartsAfterDecays == NPartsBeforeDecays)
0898     return true;
0899 
0900   bool result = true;
0901 
0902   for (int ipart = NPartsAfterDecays; ipart > NPartsBeforeDecays; ipart--) {
0903     HepMC::GenParticle *part = event().get()->barcode_to_particle(ipart);
0904 
0905     if (part->status() == 1 && (fDecayer->particleData).canDecay(part->pdg_id())) {
0906       fDecayer->event.reset();
0907       Particle py8part(part->pdg_id(),
0908                        93,
0909                        0,
0910                        0,
0911                        0,
0912                        0,
0913                        0,
0914                        0,
0915                        part->momentum().x(),
0916                        part->momentum().y(),
0917                        part->momentum().z(),
0918                        part->momentum().t(),
0919                        part->generated_mass());
0920       HepMC::GenVertex *ProdVtx = part->production_vertex();
0921       py8part.vProd(ProdVtx->position().x(), ProdVtx->position().y(), ProdVtx->position().z(), ProdVtx->position().t());
0922       py8part.tau((fDecayer->particleData).tau0(part->pdg_id()));
0923       fDecayer->event.append(py8part);
0924       int nentries = fDecayer->event.size();
0925       if (!fDecayer->event[nentries - 1].mayDecay())
0926         continue;
0927       fDecayer->next();
0928       int nentries1 = fDecayer->event.size();
0929       if (nentries1 <= nentries)
0930         continue;  //same number of particles, no decays...
0931 
0932       part->set_status(2);
0933 
0934       result = toHepMC.fill_next_event(*(fDecayer.get()), event().get(), -1, true, part);
0935     }
0936   }
0937 
0938   return result;
0939 }
0940 
0941 void Pythia8Hadronizer::finalizeEvent() {
0942   bool lhe = lheEvent() != nullptr;
0943 
0944   // protection against empty weight container
0945   if ((event()->weights()).empty())
0946     (event()->weights()).push_back(1.);
0947 
0948   // now create the GenEventInfo product from the GenEvent and fill
0949   // the missing pieces
0950   eventInfo() = std::make_unique<GenEventInfoProduct>(event().get());
0951 
0952   // in pythia pthat is used to subdivide samples into different bins
0953   // in LHE mode the binning is done by the external ME generator
0954   // which is likely not pthat, so only filling it for Py6 internal mode
0955   if (!lhe) {
0956     eventInfo()->setBinningValues(std::vector<double>(1, fMasterGen->info.pTHat()));
0957   }
0958 
0959   eventInfo()->setDJR(DJR);
0960   eventInfo()->setNMEPartons(nME);
0961   eventInfo()->setNMEPartonsFiltered(nMEFiltered);
0962 
0963   //******** Verbosity ********
0964 
0965   if (maxEventsToPrint > 0 && (pythiaPylistVerbosity || pythiaHepMCVerbosity || pythiaHepMCVerbosityParticles)) {
0966     maxEventsToPrint--;
0967     if (pythiaPylistVerbosity) {
0968       fMasterGen->info.list();
0969       fMasterGen->event.list();
0970     }
0971 
0972     if (pythiaHepMCVerbosity) {
0973       std::cout << "Event process = " << fMasterGen->info.code() << "\n"
0974                 << "----------------------" << std::endl;
0975       event()->print();
0976     }
0977     if (pythiaHepMCVerbosityParticles) {
0978       std::cout << "Event process = " << fMasterGen->info.code() << "\n"
0979                 << "----------------------" << std::endl;
0980       ascii_io->write_event(event().get());
0981     }
0982   }
0983 }
0984 
0985 std::unique_ptr<GenLumiInfoHeader> Pythia8Hadronizer::getGenLumiInfoHeader() const {
0986   auto genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();
0987 
0988   //fill lhe headers
0989   //*FIXME* initrwgt header is corrupt due to pythia bug
0990   for (const std::string &key : fMasterGen->info.headerKeys()) {
0991     genLumiInfoHeader->lheHeaders().emplace_back(key, fMasterGen->info.header(key));
0992   }
0993 
0994   //check, if it is not only nominal weight
0995   int weights_number = fMasterGen->info.nWeights();
0996   if (fMasterGen->info.initrwgt)
0997     weights_number += fMasterGen->info.initrwgt->weightsKeys.size();
0998   if (weights_number > 1) {
0999     genLumiInfoHeader->weightNames().reserve(weights_number + 1);
1000     genLumiInfoHeader->weightNames().push_back("nominal");
1001   }
1002 
1003   //fill weight names
1004   if (fMasterGen->info.initrwgt) {
1005     for (const std::string &key : fMasterGen->info.initrwgt->weightsKeys) {
1006       std::string weightgroupname;
1007       for (const auto &wgtgrp : fMasterGen->info.initrwgt->weightgroups) {
1008         const auto &wgtgrpwgt = wgtgrp.second.weights.find(key);
1009         if (wgtgrpwgt != wgtgrp.second.weights.end()) {
1010           weightgroupname = wgtgrp.first;
1011         }
1012       }
1013 
1014       std::ostringstream weightname;
1015       weightname << "LHE, id = " << key << ", ";
1016       if (!weightgroupname.empty()) {
1017         weightname << "group = " << weightgroupname << ", ";
1018       }
1019       weightname << fMasterGen->info.initrwgt->weights[key].contents;
1020       genLumiInfoHeader->weightNames().push_back(weightname.str());
1021     }
1022   }
1023 
1024   //fill shower labels
1025   // http://home.thep.lu.se/~torbjorn/pythia82html/Variations.html
1026   // http://home.thep.lu.se/~torbjorn/doxygen/classPythia8_1_1Info.html
1027   if (fMasterGen->info.nWeights() > 1) {
1028     for (int i = 0; i < fMasterGen->info.nWeights(); ++i) {
1029       genLumiInfoHeader->weightNames().push_back(fMasterGen->info.weightLabel(i));
1030     }
1031   }
1032 
1033 #if 0
1034   // VINCIA shower weights
1035   // http://vincia.hepforge.org/current/share/Vincia/htmldoc/VinciaUncertainties.html
1036   if (fvincia.get()) {
1037     for (int iVar = 0; iVar < fvincia->nWeights(); iVar++) {
1038       genLumiInfoHeader->weightNames().push_back(fvincia->weightLabel(iVar));
1039     }
1040   }
1041 
1042   if (fDire.get()) {
1043     //Make sure the base weight comes first
1044     genLumiInfoHeader->weightNames().push_back("base");
1045 
1046     unordered_map<string, double>::iterator it;
1047     for (it = fDire->weightsPtr->getShowerWeights()->begin(); it != fDire->weightsPtr->getShowerWeights()->end();
1048          it++) {
1049       if (it->first == "base")
1050         continue;
1051       genLumiInfoHeader->weightNames().push_back(it->first);
1052     }
1053   }
1054 #endif
1055 
1056   return genLumiInfoHeader;
1057 }
1058 
1059 typedef edm::GeneratorFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8GeneratorFilter;
1060 DEFINE_FWK_MODULE(Pythia8GeneratorFilter);
1061 
1062 typedef edm::HadronizerFilter<Pythia8Hadronizer, ExternalDecayDriver> Pythia8HadronizerFilter;
1063 DEFINE_FWK_MODULE(Pythia8HadronizerFilter);
1064 
1065 typedef edm::ConcurrentGeneratorFilter<Pythia8Hadronizer, ConcurrentExternalDecayDriver>
1066     Pythia8ConcurrentGeneratorFilter;
1067 DEFINE_FWK_MODULE(Pythia8ConcurrentGeneratorFilter);
1068 
1069 typedef edm::ConcurrentHadronizerFilter<Pythia8Hadronizer, ConcurrentExternalDecayDriver>
1070     Pythia8ConcurrentHadronizerFilter;
1071 DEFINE_FWK_MODULE(Pythia8ConcurrentHadronizerFilter);