Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-02 05:09:45

0001 #include <memory>
0002 #include <string>
0003 
0004 #include "GeneratorInterface/Pythia8Interface/interface/Py8InterfaceBase.h"
0005 #include "GeneratorInterface/Pythia8Interface/interface/SLHAReaderBase.h"
0006 
0007 #include "FWCore/Utilities/interface/Exception.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include <filesystem>
0010 
0011 // EvtGen plugin
0012 //
0013 //#include "Pythia8Plugins/EvtGen.h"
0014 
0015 using namespace Pythia8;
0016 
0017 namespace gen {
0018 
0019   Py8InterfaceBase::Py8InterfaceBase(edm::ParameterSet const& ps)
0020       : BaseHadronizer(ps), useEvtGen(false), evtgenDecays(nullptr) {
0021     fParameters = ps;
0022 
0023     pythiaPylistVerbosity = ps.getUntrackedParameter<int>("pythiaPylistVerbosity", 0);
0024     pythiaHepMCVerbosity = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false);
0025     pythiaHepMCVerbosityParticles = ps.getUntrackedParameter<bool>("pythiaHepMCVerbosityParticles", false);
0026     maxEventsToPrint = ps.getUntrackedParameter<int>("maxEventsToPrint", 0);
0027     p8RndmEngine_ = std::make_shared<P8RndmEngine>();
0028 
0029     if (pythiaHepMCVerbosityParticles)
0030       ascii_io = new HepMC::IO_AsciiParticles("cout", std::ios::out);
0031 
0032     if (ps.exists("useEvtGenPlugin")) {
0033       useEvtGen = true;
0034       string evtgenpath(std::getenv("EVTGENDATA"));
0035       evtgenDecFile = evtgenpath + string("/DECAY_2010.DEC");
0036       evtgenPdlFile = evtgenpath + string("/evt.pdl");
0037 
0038       if (ps.exists("evtgenDecFile")) {
0039         edm::FileInPath decay_table(ps.getParameter<std::string>("evtgenDecFile"));
0040         evtgenDecFile = decay_table.fullPath();
0041       }
0042 
0043       if (ps.exists("evtgenPdlFile")) {
0044         edm::FileInPath pdt(ps.getParameter<std::string>("evtgenPdlFile"));
0045         evtgenPdlFile = pdt.fullPath();
0046       }
0047 
0048       if (ps.exists("evtgenUserFile")) {
0049         std::vector<std::string> user_decays = ps.getParameter<std::vector<std::string> >("evtgenUserFile");
0050         for (unsigned int i = 0; i < user_decays.size(); i++) {
0051           edm::FileInPath user_decay(user_decays.at(i));
0052           evtgenUserFiles.push_back(user_decay.fullPath());
0053         }
0054         //evtgenUserFiles = ps.getParameter< std::vector<std::string> >("evtgenUserFile");
0055       }
0056 
0057       if (ps.exists("evtgenUserFileEmbedded")) {
0058         std::vector<std::string> user_decay_lines =
0059             ps.getParameter<std::vector<std::string> >("evtgenUserFileEmbedded");
0060         char tempslhaname[] = "pythia8evtgenXXXXXX";
0061         int fd = mkstemp(tempslhaname);
0062 
0063         for (unsigned int i = 0; i < user_decay_lines.size(); i++) {
0064           user_decay_lines.at(i) += "\n";
0065           write(fd, user_decay_lines.at(i).c_str(), user_decay_lines.at(i).size());
0066         }
0067         close(fd);
0068         evtgenUserFiles.push_back(std::string(tempslhaname));
0069       }
0070     }
0071   }
0072 
0073   bool Py8InterfaceBase::readSettings(int) {
0074     //Pythia 8's default value for first argument to constructor
0075     const string xmlDir = "../share/Pythia8/xmldoc";
0076     bool printBanner = true;
0077     if (fParameters.exists("printBanner")) {
0078       printBanner = fParameters.getUntrackedParameter<bool>("printBanner");
0079     }
0080     if (!fMasterGen.get())
0081       fMasterGen = std::make_unique<Pythia>(xmlDir, printBanner);
0082     fDecayer = std::make_unique<Pythia>(xmlDir, printBanner);
0083 
0084     //add settings for resonance decay filter
0085     fMasterGen->settings.addFlag("BiasedTauDecayer:filter", false);
0086     fMasterGen->settings.addFlag("BiasedTauDecayer:eDecays", true);
0087     fMasterGen->settings.addFlag("BiasedTauDecayer:muDecays", true);
0088 
0089     //add settings for resonance decay filter
0090     fMasterGen->settings.addFlag("ResonanceDecayFilter:filter", false);
0091     fMasterGen->settings.addFlag("ResonanceDecayFilter:exclusive", false);
0092     fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuAsEquivalent", false);
0093     fMasterGen->settings.addFlag("ResonanceDecayFilter:eMuTauAsEquivalent", false);
0094     fMasterGen->settings.addFlag("ResonanceDecayFilter:allNuAsEquivalent", false);
0095     fMasterGen->settings.addFlag("ResonanceDecayFilter:udscAsEquivalent", false);
0096     fMasterGen->settings.addFlag("ResonanceDecayFilter:udscbAsEquivalent", false);
0097     fMasterGen->settings.addFlag("ResonanceDecayFilter:wzAsEquivalent", false);
0098     fMasterGen->settings.addMVec("ResonanceDecayFilter:mothers", std::vector<int>(), false, false, 0, 0);
0099     fMasterGen->settings.addMVec("ResonanceDecayFilter:daughters", std::vector<int>(), false, false, 0, 0);
0100 
0101     //add settings for PT filter
0102     fMasterGen->settings.addFlag("PTFilter:filter", false);
0103     fMasterGen->settings.addMode("PTFilter:quarkToFilter", 5, true, true, 3, 6);
0104     fMasterGen->settings.addParm("PTFilter:scaleToFilter", 0.4, true, true, 0.0, 10.);
0105     fMasterGen->settings.addParm("PTFilter:quarkRapidity", 10.0, true, true, 0.0, 10.);
0106     fMasterGen->settings.addParm("PTFilter:quarkPt", -.1, true, true, -.1, 100.);
0107 
0108     //add settings for RecoilToTop tool
0109     fMasterGen->settings.addFlag("TopRecoilHook:doTopRecoilIn", false);
0110     fMasterGen->settings.addFlag("TopRecoilHook:useOldDipoleIn", false);
0111     fMasterGen->settings.addFlag("TopRecoilHook:doListIn", false);
0112 
0113     //add settings for powheg resonance scale calculation
0114     fMasterGen->settings.addFlag("POWHEGres:calcScales", false);
0115     fMasterGen->settings.addFlag("POWHEG:bb4l", false);
0116     fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:veto", false);
0117     fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoQED", false);
0118     fMasterGen->settings.addFlag("POWHEG:bb4l:DEBUG", false);
0119     fMasterGen->settings.addFlag("POWHEG:bb4l:ScaleResonance:veto", false);
0120     fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:vetoDipoleFrame", false);
0121     fMasterGen->settings.addFlag("POWHEG:bb4l:FSREmission:pTpythiaVeto", false);
0122     fMasterGen->settings.addParm("POWHEG:bb4l:pTminVeto", 10.0, true, true, 0.0, 10.);
0123     fMasterGen->settings.addFlag("POWHEG:bb4l:vetoAllRadtypes", false);
0124 
0125     fMasterGen->setRndmEnginePtr(p8RndmEngine_);
0126     fDecayer->setRndmEnginePtr(p8RndmEngine_);
0127 
0128     fMasterGen->readString("Next:numberShowEvent = 0");
0129     fDecayer->readString("Next:numberShowEvent = 0");
0130 
0131     edm::ParameterSet currentParameters;
0132     if (randomIndex() >= 0) {
0133       std::vector<edm::ParameterSet> randomizedParameters =
0134           fParameters.getParameter<std::vector<edm::ParameterSet> >("RandomizedParameters");
0135       currentParameters = randomizedParameters[randomIndex()];
0136     } else {
0137       currentParameters = fParameters;
0138     }
0139 
0140     ParameterCollector pCollector = currentParameters.getParameter<edm::ParameterSet>("PythiaParameters");
0141 
0142     for (ParameterCollector::const_iterator line = pCollector.begin(); line != pCollector.end(); ++line) {
0143       if (line->find("Random:") != std::string::npos)
0144         throw cms::Exception("PythiaError") << "Attempted to set random number "
0145                                                "using Pythia commands. Please use "
0146                                                "the RandomNumberGeneratorService."
0147                                             << std::endl;
0148 
0149       if (!fMasterGen->readString(*line))
0150         throw cms::Exception("PythiaError") << "Pythia 8 did not accept \"" << *line << "\"." << std::endl;
0151 
0152       if (line->find("ParticleDecays:") != std::string::npos) {
0153         if (!fDecayer->readString(*line))
0154           throw cms::Exception("PythiaError") << "Pythia 8 Decayer did not accept \"" << *line << "\"." << std::endl;
0155       }
0156     }
0157 
0158     slhafile_.clear();
0159 
0160     if (currentParameters.exists("SLHAFileForPythia8")) {
0161       std::string slhafilenameshort = currentParameters.getParameter<std::string>("SLHAFileForPythia8");
0162       edm::FileInPath f1(slhafilenameshort);
0163 
0164       fMasterGen->settings.mode("SLHA:readFrom", 2);
0165       fMasterGen->settings.word("SLHA:file", f1.fullPath());
0166     } else if (currentParameters.exists("SLHATableForPythia8")) {
0167       std::string slhatable = currentParameters.getParameter<std::string>("SLHATableForPythia8");
0168 
0169       makeTmpSLHA(slhatable);
0170     } else if (currentParameters.exists("SLHATreeForPythia8")) {
0171       auto slhaReaderParams = currentParameters.getParameter<edm::ParameterSet>("SLHATreeForPythia8");
0172       std::unique_ptr<SLHAReaderBase> reader =
0173           SLHAReaderFactory::get()->create(slhaReaderParams.getParameter<std::string>("name"), slhaReaderParams);
0174 
0175       makeTmpSLHA(reader->getSLHA(currentParameters.getParameter<std::string>("ConfigDescription")));
0176     }
0177 
0178     return true;
0179   }
0180 
0181   void Py8InterfaceBase::makeTmpSLHA(const std::string& slhatable) {
0182     char tempslhaname[] = "pythia8SLHAtableXXXXXX";
0183     int fd = mkstemp(tempslhaname);
0184     write(fd, slhatable.c_str(), slhatable.size());
0185     close(fd);
0186 
0187     slhafile_ = tempslhaname;
0188 
0189     fMasterGen->settings.mode("SLHA:readFrom", 2);
0190     fMasterGen->settings.word("SLHA:file", slhafile_);
0191   }
0192 
0193   bool Py8InterfaceBase::declareStableParticles(const std::vector<int>& pdgIds) {
0194     for (size_t i = 0; i < pdgIds.size(); i++) {
0195       // FIXME: need to double check if PID's are the same in Py6 & Py8,
0196       //        because the HepPDT translation tool is actually for **Py6**
0197       //
0198       // well, actually it looks like Py8 operates in PDT id's rather than Py6's
0199       //
0200       //    int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
0201       int PyID = pdgIds[i];
0202       std::ostringstream pyCard;
0203       pyCard << PyID << ":mayDecay=false";
0204 
0205       if (fMasterGen->particleData.isParticle(PyID)) {
0206         fMasterGen->readString(pyCard.str());
0207       } else {
0208         edm::LogWarning("DataNotUnderstood") << "Pythia8 does not "
0209                                              << "recognize particle id = " << PyID << std::endl;
0210       }
0211       // alternative:
0212       // set the 2nd input argument warn=false
0213       // - this way Py8 will NOT print warnings about unknown particle code(s)
0214       // fMasterPtr->readString( pyCard.str(), false )
0215     }
0216 
0217     return true;
0218   }
0219 
0220   bool Py8InterfaceBase::declareSpecialSettings(const std::vector<std::string>& settings) {
0221     for (unsigned int iss = 0; iss < settings.size(); iss++) {
0222       if (settings[iss].find("QED-brem-off") != std::string::npos) {
0223         fMasterGen->readString("TimeShower:QEDshowerByL=off");
0224       } else {
0225         size_t fnd1 = settings[iss].find("Pythia8:");
0226         if (fnd1 != std::string::npos) {
0227           std::string value = settings[iss].substr(fnd1 + 8);
0228           fDecayer->readString(value);
0229         }
0230       }
0231     }
0232     return true;
0233   }
0234 
0235   void Py8InterfaceBase::statistics() {
0236     fMasterGen->stat();
0237     return;
0238   }
0239 
0240 }  // namespace gen