Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 12:50:10

0001 //
0002 // See header file for description
0003 //
0004 
0005 #include "CommonTools/TriggerUtils/interface/PrescaleWeightProvider.h"
0006 
0007 #include <sstream>
0008 
0009 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtUtils.h"
0010 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0011 #include "DataFormats/Common/interface/TriggerResults.h"
0012 #include "DataFormats/L1GlobalTrigger/interface/L1GtTriggerMenuLite.h"
0013 
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/BranchType.h"
0017 
0018 PrescaleWeightProvider::PrescaleWeightProvider(const edm::ParameterSet& config, edm::ConsumesCollector& iC)
0019     // default values
0020     : init_(false),
0021       verbosity_(0),
0022       triggerResultsTag_("TriggerResults::HLT"),
0023       triggerResultsToken_(iC.mayConsume<edm::TriggerResults>(triggerResultsTag_)),
0024       l1GtTriggerMenuLiteTag_("l1GtTriggerMenuLite"),
0025       l1GtTriggerMenuLiteToken_(iC.mayConsume<L1GtTriggerMenuLite, edm::InRun>(l1GtTriggerMenuLiteTag_)) {
0026   hltPaths_.clear();
0027   if (config.exists("prescaleWeightVerbosityLevel"))
0028     verbosity_ = config.getParameter<unsigned>("prescaleWeightVerbosityLevel");
0029   if (config.exists("prescaleWeightTriggerResults"))
0030     triggerResultsTag_ = config.getParameter<edm::InputTag>("prescaleWeightTriggerResults");
0031   if (config.exists("prescaleWeightL1GtTriggerMenuLite"))
0032     l1GtTriggerMenuLiteTag_ = config.getParameter<edm::InputTag>("prescaleWeightL1GtTriggerMenuLite");
0033   if (config.exists("prescaleWeightHltPaths"))
0034     hltPaths_ = config.getParameter<std::vector<std::string> >("prescaleWeightHltPaths");
0035 
0036   configured_ = true;
0037   if (triggerResultsTag_.process().empty()) {
0038     configured_ = false;
0039     if (verbosity_ > 0)
0040       edm::LogWarning("PrescaleWeightProvider") << "Process name not configured via TriggerResults InputTag";
0041   } else if (triggerResultsTag_.label().empty()) {
0042     configured_ = false;
0043     if (verbosity_ > 0)
0044       edm::LogWarning("PrescaleWeightProvider") << "TriggerResults label not configured";
0045   } else if (l1GtTriggerMenuLiteTag_.label().empty()) {
0046     configured_ = false;
0047     if (verbosity_ > 0)
0048       edm::LogWarning("PrescaleWeightProvider") << "L1GtTriggerMenuLite label not configured";
0049   } else if (hltPaths_.empty()) {
0050     configured_ = false;
0051     if (verbosity_ > 0)
0052       edm::LogError("PrescaleWeightProvider") << "HLT paths of interest not configured";
0053   }
0054   if (configured_) {
0055     triggerResultsToken_ = iC.mayConsume<edm::TriggerResults>(triggerResultsTag_);
0056     l1GtTriggerMenuLiteToken_ = iC.mayConsume<L1GtTriggerMenuLite, edm::InRun>(l1GtTriggerMenuLiteTag_);
0057   }
0058 }
0059 
0060 void PrescaleWeightProvider::initRun(const edm::Run& run, const edm::EventSetup& setup) {
0061   init_ = true;
0062 
0063   if (!configured_) {
0064     init_ = false;
0065     if (verbosity_ > 0)
0066       edm::LogWarning("PrescaleWeightProvider") << "Run initialisation failed due to failing configuration";
0067     return;
0068   }
0069 
0070   HLTConfigProvider const& hltConfig = hltPrescaleProvider_->hltConfigProvider();
0071   bool hltChanged(false);
0072   if (!hltPrescaleProvider_->init(run, setup, triggerResultsTag_.process(), hltChanged)) {
0073     if (verbosity_ > 0)
0074       edm::LogError("PrescaleWeightProvider")
0075           << "HLT config initialization error with process name \"" << triggerResultsTag_.process() << "\"";
0076     init_ = false;
0077   } else if (hltConfig.size() <= 0) {
0078     if (verbosity_ > 0)
0079       edm::LogError("PrescaleWeightProvider") << "HLT config size error";
0080     init_ = false;
0081   } else if (hltChanged) {
0082     if (verbosity_ > 0)
0083       edm::LogInfo("PrescaleWeightProvider") << "HLT configuration changed";
0084   }
0085   if (!init_)
0086     return;
0087 
0088   run.getByToken(l1GtTriggerMenuLiteToken_, triggerMenuLite_);
0089   if (!triggerMenuLite_.isValid()) {
0090     if (verbosity_ > 0)
0091       edm::LogError("PrescaleWeightProvider")
0092           << "L1GtTriggerMenuLite with label \"" << l1GtTriggerMenuLiteTag_.label() << "\" not found";
0093     init_ = false;
0094   }
0095 }
0096 
0097 int PrescaleWeightProvider::prescaleWeight(const edm::Event& event, const edm::EventSetup& setup) {
0098   if (!init_)
0099     return 1;
0100 
0101   // L1
0102   L1GtUtils const& l1GtUtils = hltPrescaleProvider_->l1GtUtils();
0103 
0104   // HLT
0105   HLTConfigProvider const& hltConfig = hltPrescaleProvider_->hltConfigProvider();
0106 
0107   edm::Handle<edm::TriggerResults> triggerResults;
0108   event.getByToken(triggerResultsToken_, triggerResults);
0109   if (!triggerResults.isValid()) {
0110     if (verbosity_ > 0)
0111       edm::LogError("PrescaleWeightProvider::prescaleWeight")
0112           << "TriggerResults product not found for InputTag \"" << triggerResultsTag_.encode() << "\"";
0113     return 1;
0114   }
0115 
0116   const int SENTINEL(-1);
0117   int weight(SENTINEL);
0118 
0119   for (unsigned ui = 0; ui < hltPaths_.size(); ui++) {
0120     const std::string hltPath(hltPaths_.at(ui));
0121     unsigned hltIndex(hltConfig.triggerIndex(hltPath));
0122     if (hltIndex == hltConfig.size()) {
0123       if (verbosity_ > 0)
0124         edm::LogError("PrescaleWeightProvider::prescaleWeight") << "HLT path \"" << hltPath << "\" does not exist";
0125       continue;
0126     }
0127     if (!triggerResults->accept(hltIndex))
0128       continue;
0129 
0130     const std::vector<std::pair<bool, std::string> >& level1Seeds = hltConfig.hltL1GTSeeds(hltPath);
0131     if (level1Seeds.size() != 1) {
0132       if (verbosity_ > 0)
0133         edm::LogError("PrescaleWeightProvider::prescaleWeight")
0134             << "HLT path \"" << hltPath << "\" provides too many L1 seeds";
0135       return 1;
0136     }
0137     parseL1Seeds(level1Seeds.at(0).second);
0138     if (l1SeedPaths_.empty()) {
0139       if (verbosity_ > 0)
0140         edm::LogWarning("PrescaleWeightProvider::prescaleWeight")
0141             << "Failed to parse L1 seeds for HLT path \"" << hltPath << "\"";
0142       continue;
0143     }
0144 
0145     int l1Prescale(SENTINEL);
0146     for (unsigned uj = 0; uj < l1SeedPaths_.size(); uj++) {
0147       int l1TempPrescale(SENTINEL);
0148       int errorCode(0);
0149       if (level1Seeds.at(0).first) {  // technical triggers
0150         unsigned techBit(atoi(l1SeedPaths_.at(uj).c_str()));
0151         const std::string techName(*(triggerMenuLite_->gtTechTrigName(techBit, errorCode)));
0152         if (errorCode != 0)
0153           continue;
0154         if (!l1GtUtils.decision(event, techName, errorCode))
0155           continue;
0156         if (errorCode != 0)
0157           continue;
0158         l1TempPrescale = l1GtUtils.prescaleFactor(event, techName, errorCode);
0159         if (errorCode != 0)
0160           continue;
0161       } else {  // algorithmic triggers
0162         if (!l1GtUtils.decision(event, l1SeedPaths_.at(uj), errorCode))
0163           continue;
0164         if (errorCode != 0)
0165           continue;
0166         l1TempPrescale = l1GtUtils.prescaleFactor(event, l1SeedPaths_.at(uj), errorCode);
0167         if (errorCode != 0)
0168           continue;
0169       }
0170       if (l1TempPrescale > 0) {
0171         if (l1Prescale == SENTINEL || l1Prescale > l1TempPrescale)
0172           l1Prescale = l1TempPrescale;
0173       }
0174     }
0175     if (l1Prescale == SENTINEL) {
0176       if (verbosity_ > 0)
0177         edm::LogError("PrescaleWeightProvider::prescaleWeight")
0178             << "Unable to find the L1 prescale for HLT path \"" << hltPath << "\"";
0179       continue;
0180     }
0181 
0182     int hltPrescale(hltPrescaleProvider_->prescaleValue(event, setup, hltPath));
0183 
0184     if (hltPrescale * l1Prescale > 0) {
0185       if (weight == SENTINEL || weight > hltPrescale * l1Prescale) {
0186         weight = hltPrescale * l1Prescale;
0187       }
0188     }
0189   }
0190 
0191   if (weight == SENTINEL) {
0192     if (verbosity_ > 0)
0193       edm::LogWarning("PrescaleWeightProvider::prescaleWeight")
0194           << "No valid weight for any requested HLT path, returning default weight of 1";
0195     return 1;
0196   }
0197   return weight;
0198 }
0199 
0200 void PrescaleWeightProvider::parseL1Seeds(const std::string& l1Seeds) {
0201   l1SeedPaths_.clear();
0202   std::stringstream ss(l1Seeds);
0203   std::string buf;
0204 
0205   while (ss.good() && !ss.eof()) {
0206     ss >> buf;
0207     if (buf[0] == '(' || buf[buf.size() - 1] == ')' || buf == "AND" || buf == "NOT") {
0208       l1SeedPaths_.clear();
0209       if (verbosity_ > 0)
0210         edm::LogWarning("PrescaleWeightProvider::parseL1Seeds") << "Only supported logical expression is OR";
0211       return;
0212     } else if (buf == "OR")
0213       continue;
0214     else
0215       l1SeedPaths_.push_back(buf);
0216   }
0217 }