Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-15 22:37:05

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