Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:58:33

0001 #include "DQMOffline/Trigger/interface/EgHLTOfflineSummaryClient.h"
0002 
0003 #include "FWCore/Framework/interface/Run.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ServiceRegistry/interface/Service.h"
0007 
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 
0010 #include "DQMServices/Core/interface/DQMStore.h"
0011 
0012 #include "DQMOffline/Trigger/interface/EgHLTTrigTools.h"
0013 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0014 #include <boost/algorithm/string.hpp>
0015 #include <fnmatch.h>
0016 
0017 EgHLTOfflineSummaryClient::EgHLTOfflineSummaryClient(const edm::ParameterSet& iConfig)
0018     : egHLTSumHistName_("egHLTTrigSum"), isSetup_(false) {
0019   usesResource("DQMStore");
0020   dirName_ = iConfig.getParameter<std::string>(
0021       "DQMDirName");  //only one chance to get this, if we every have another shot, remember to check isSetup is okay
0022   dbe_ = edm::Service<DQMStore>().operator->();
0023   if (!dbe_) {
0024     edm::LogError("EgHLTOfflineSummaryClient")
0025         << "unable to get DQMStore service, no summary histograms will be produced";
0026   } else {
0027     dbe_->setCurrentFolder(dirName_);
0028   }
0029 
0030   eleHLTFilterNames_ = iConfig.getParameter<std::vector<std::string> >("eleHLTFilterNames");
0031   phoHLTFilterNames_ = iConfig.getParameter<std::vector<std::string> >("phoHLTFilterNames");
0032   eleHLTFilterNamesForSumBit_ = iConfig.getParameter<std::vector<std::string> >("eleHLTFilterNamesForSumBit");
0033   phoHLTFilterNamesForSumBit_ = iConfig.getParameter<std::vector<std::string> >("phoHLTFilterNamesForSumBit");
0034 
0035   filterInactiveTriggers_ = iConfig.getParameter<bool>("filterInactiveTriggers");
0036   hltTag_ = iConfig.getParameter<std::string>("hltTag");
0037 
0038   usePathNames_ = iConfig.getParameter<bool>("usePathNames");
0039 
0040   //std::vector<std::string> egHLTSumQTests = iConfig.getParameter<std::vector<std::string> >("egHLTSumQTests");
0041   // splitStringsToPairs_(egHLTSumQTests,egHLTSumHistXBins_);
0042 
0043   fillQTestData_(iConfig, egHLTSumHistXBins_, "egHLTSumQTests");
0044   fillQTestData_(iConfig, eleQTestsForSumBit_, "egHLTEleQTestsForSumBit");
0045   fillQTestData_(iConfig, phoQTestsForSumBit_, "egHLTPhoQTestsForSumBit");
0046 
0047   runClientEndLumiBlock_ = iConfig.getParameter<bool>("runClientEndLumiBlock");
0048   runClientEndRun_ = iConfig.getParameter<bool>("runClientEndRun");
0049   runClientEndJob_ = iConfig.getParameter<bool>("runClientEndJob");
0050 
0051   //egHLTSumHistXBins_.push_back(std::make_pair("Ele Rel Trig Eff",&EgHLTOfflineSummaryClient::eleTrigRelEffQTestResult_));
0052   //egHLTSumHistXBins_.push_back(std::make_pair("Pho Rel Trig Eff",&EgHLTOfflineSummaryClient::phoTrigRelEffQTestResult_));
0053   //egHLTSumHistXBins_.push_back(std::make_pair("Ele T&P Trig Eff",&EgHLTOfflineSummaryClient::eleTrigTPEffQTestResult_));
0054   //egHLTSumHistXBins_.push_back(std::make_pair("Triggered Ele",&EgHLTOfflineSummaryClient::trigEleQTestResult_));
0055   //egHLTSumHistXBins_.push_back(std::make_pair("Triggered Pho",&EgHLTOfflineSummaryClient::trigPhoQTestResult_));
0056 }
0057 
0058 EgHLTOfflineSummaryClient::~EgHLTOfflineSummaryClient() = default;
0059 
0060 void EgHLTOfflineSummaryClient::beginJob() {}
0061 
0062 void EgHLTOfflineSummaryClient::endJob() {
0063   if (runClientEndJob_)
0064     runClient_();
0065 }
0066 
0067 void EgHLTOfflineSummaryClient::beginRun(const edm::Run& run, const edm::EventSetup& c) {
0068   if (!isSetup_) {
0069     bool changed;
0070     HLTConfigProvider hltConfig;
0071     hltConfig.init(run, c, hltTag_, changed);
0072     if (filterInactiveTriggers_) {
0073       std::vector<std::string> activeFilters;
0074       std::vector<std::string> activeEleFilters;
0075       std::vector<std::string> activeEle2LegFilters;
0076       std::vector<std::string> activePhoFilters;
0077       std::vector<std::string> activePho2LegFilters;
0078 
0079       egHLT::trigTools::getActiveFilters(
0080           hltConfig, activeFilters, activeEleFilters, activeEle2LegFilters, activePhoFilters, activePho2LegFilters);
0081 
0082       egHLT::trigTools::filterInactiveTriggers(eleHLTFilterNames_, activeFilters);
0083       egHLT::trigTools::filterInactiveTriggers(phoHLTFilterNames_, activePhoFilters);
0084       egHLT::trigTools::filterInactiveTriggers(eleHLTFilterNamesForSumBit_, activeEleFilters);
0085       egHLT::trigTools::filterInactiveTriggers(phoHLTFilterNamesForSumBit_, activePhoFilters);
0086     }
0087     getEgHLTFiltersToMon_(egHLTFiltersToMon_);
0088 
0089     if (usePathNames_)
0090       egHLT::trigTools::translateFiltersToPathNames(hltConfig, egHLTFiltersToMon_, egHLTFiltersToMonPaths_);
0091     isSetup_ = true;
0092   }
0093 }
0094 
0095 void EgHLTOfflineSummaryClient::endRun(const edm::Run& run, const edm::EventSetup& c) {
0096   if (runClientEndRun_)
0097     runClient_();
0098 }
0099 
0100 //dummy analysis function
0101 void EgHLTOfflineSummaryClient::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {}
0102 
0103 void EgHLTOfflineSummaryClient::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& c) {
0104   if (runClientEndLumiBlock_)
0105     runClient_();
0106 }
0107 
0108 void EgHLTOfflineSummaryClient::runClient_() {
0109   MonitorElement* egHLTSumME = getEgHLTSumHist_();
0110 
0111   for (size_t filterNr = 0; filterNr < egHLTFiltersToMon_.size(); filterNr++) {
0112     for (size_t xBinNr = 0; xBinNr < egHLTSumHistXBins_.size(); xBinNr++) {
0113       //egHLTSumHist->setBinContent(xBinNr+1,filterNr+1,(*egHLTSumHistXBins_[xBinNr].second)(egHLTFiltersToMon_[filterNr].c_str()));
0114       egHLTSumME->setBinContent(
0115           xBinNr + 1,
0116           filterNr + 1,
0117           getQTestResults_(egHLTFiltersToMon_[filterNr], egHLTSumHistXBins_[xBinNr].qTestPatterns));
0118     }
0119   }
0120 
0121   MonitorElement* hltEleSumBit = dbe_->get("HLT/EventInfo/reportSummaryContents/HLT_Electron");
0122   MonitorElement* hltPhoSumBit = dbe_->get("HLT/EventInfo/reportSummaryContents/HLT_Photon");
0123   dbe_->setCurrentFolder("HLT/EventInfo/reportSummaryContents/");
0124   if (hltEleSumBit == nullptr)
0125     hltEleSumBit = dbe_->bookFloat("HLT_Electron");
0126   if (hltPhoSumBit == nullptr)
0127     hltPhoSumBit = dbe_->bookFloat("HLT_Photon");
0128 
0129   float eleSumBit = 1.;
0130   for (size_t filterNr = 0; filterNr < eleHLTFilterNamesForSumBit_.size() && eleSumBit == 1;
0131        filterNr++) {  //breaks as soon as a test fails
0132     for (size_t testNr = 0; testNr < eleQTestsForSumBit_.size() && eleSumBit == 1; testNr++) {
0133       if (getQTestResults_(eleHLTFilterNamesForSumBit_[filterNr], eleQTestsForSumBit_[testNr].qTestPatterns) == 0)
0134         eleSumBit = 0;
0135     }
0136   }
0137   hltEleSumBit->Fill(eleSumBit);
0138 
0139   float phoSumBit = 1.;
0140   for (size_t filterNr = 0; filterNr < phoHLTFilterNamesForSumBit_.size() && phoSumBit == 1;
0141        filterNr++) {  //breaks as soon as a test fails
0142     for (size_t testNr = 0; testNr < phoQTestsForSumBit_.size() && phoSumBit == 1; testNr++) {
0143       if (getQTestResults_(phoHLTFilterNamesForSumBit_[filterNr], phoQTestsForSumBit_[testNr].qTestPatterns) == 0)
0144         phoSumBit = 0;
0145     }
0146   }
0147   hltPhoSumBit->Fill(phoSumBit);
0148 }
0149 void EgHLTOfflineSummaryClient::splitStringsToPairs_(const std::vector<std::string>& stringsToSplit,
0150                                                      std::vector<std::pair<std::string, std::string> >& splitStrings) {
0151   splitStrings.clear();
0152   splitStrings.reserve(stringsToSplit.size());
0153   for (auto const& stringNr : stringsToSplit) {
0154     std::vector<std::string> tempSplitStrings;
0155     boost::split(tempSplitStrings, stringNr, boost::is_any_of(std::string(":")));
0156     if (tempSplitStrings.size() == 2) {
0157       splitStrings.push_back(std::make_pair(tempSplitStrings[0], tempSplitStrings[1]));
0158     } else {
0159       edm::LogWarning("EgHLTOfflineSummaryClient")
0160           << " Error : entry " << stringNr
0161           << " is not of form A:B, ignoring (ie this quailty test isnt being included in the sumamry hist) ";
0162     }
0163   }
0164 }
0165 
0166 EgHLTOfflineSummaryClient::MonitorElement* EgHLTOfflineSummaryClient::getEgHLTSumHist_() {
0167   MonitorElement* egHLTSumHist = dbe_->get(dirName_ + "/" + egHLTSumHistName_);
0168   if (egHLTSumHist == nullptr) {
0169     auto* hist = new TH2F(egHLTSumHistName_.c_str(),
0170                           "E/g HLT Offline Summary",
0171                           egHLTSumHistXBins_.size(),
0172                           0.,
0173                           1.,
0174                           egHLTFiltersToMon_.size(),
0175                           0.,
0176                           1.);
0177     for (size_t xBinNr = 0; xBinNr < egHLTSumHistXBins_.size(); xBinNr++) {
0178       hist->GetXaxis()->SetBinLabel(xBinNr + 1, egHLTSumHistXBins_[xBinNr].name.c_str());
0179     }
0180 
0181     for (size_t yBinNr = 0; yBinNr < egHLTFiltersToMon_.size(); yBinNr++) {
0182       if (usePathNames_ && egHLTFiltersToMonPaths_.size() == egHLTFiltersToMon_.size()) {
0183         hist->GetYaxis()->SetBinLabel(yBinNr + 1, egHLTFiltersToMonPaths_[yBinNr].c_str());
0184       } else {
0185         hist->GetYaxis()->SetBinLabel(yBinNr + 1, egHLTFiltersToMon_[yBinNr].c_str());
0186       }
0187     }
0188     for (size_t xBinNr = 0; xBinNr < egHLTSumHistXBins_.size(); xBinNr++) {
0189       for (size_t yBinNr = 0; yBinNr < egHLTFiltersToMon_.size(); yBinNr++) {
0190         hist->SetBinContent(xBinNr + 1, yBinNr + 1, -2);
0191       }
0192     }
0193 
0194     dbe_->setCurrentFolder(dirName_);
0195     egHLTSumHist = dbe_->book2D(egHLTSumHistName_, hist);
0196   }
0197   return egHLTSumHist;
0198 }
0199 
0200 //this function puts every e/g trigger monitored in a std::vector
0201 //this is *very* similar to EgHLTOfflineSource::getHLTFilterNamesUsed but
0202 //it differs in the fact it only gets the E/g primary triggers not the backups
0203 //due to the design, to ensure we get every filter, filters will be inserted multiple times
0204 //eg electron filters will contain photon triggers which are also in the photon filters
0205 //but only want one copy in the vector
0206 //this function is intended to be called once per job so some inefficiency can can be tolerated
0207 //therefore we will use a std::set to ensure ensure that each filtername is only inserted once
0208 //and then convert to a std::vector
0209 void EgHLTOfflineSummaryClient::getEgHLTFiltersToMon_(std::vector<std::string>& filterNames) const {
0210   std::set<std::string> filterNameSet;
0211   for (auto const& eleHLTFilterName : eleHLTFilterNames_)
0212     filterNameSet.insert(eleHLTFilterName);
0213   for (auto const& phoHLTFilterName : phoHLTFilterNames_)
0214     filterNameSet.insert(phoHLTFilterName);
0215 
0216   //right all the triggers are inserted once and only once in the set, convert to vector
0217   //very lazy, create a new vector so can use the constructor and then use swap to transfer
0218   std::vector<std::string>(filterNameSet.begin(), filterNameSet.end()).swap(filterNames);
0219 }
0220 
0221 //only returns 0 or 1, 0 is bad, one is good and if the test is not found defaults to good
0222 //(this is because its a dumb algorithm, photon tests are run for electron triggers which unsurprisingly are not found)
0223 int EgHLTOfflineSummaryClient::getQTestResults_(const std::string& filterName,
0224                                                 const std::vector<std::string>& patterns) const {
0225   int nrFail = 0;
0226   int nrQTests = 0;
0227   for (auto const& pattern : patterns) {
0228     auto filterpattern = filterName + pattern;
0229     std::vector<MonitorElement*> monElems = dbe_->getAllContents(dirName_);
0230     // std::cout <<"mon elem "<<dirName_+"/"+filterName+patterns[patternNr]<<"nr monElems "<<monElems.size()<<std::endl;
0231     for (auto& monElem : monElems) {
0232       const auto& name = monElem->getName();
0233       int match = fnmatch(filterpattern.c_str(), name.c_str(), 0);
0234       if (match == FNM_NOMATCH)
0235         continue;
0236 
0237       std::vector<QReport*> qTests = monElem->getQReports();
0238       nrQTests += qTests.size();
0239       //  std::cout <<monElems[monElemNr]->getName()<<" "<<monElems[monElemNr]->hasError()<<" nr test "<<qTests.size()<<std::endl;
0240       if (monElem->hasError())
0241         nrFail++;
0242     }
0243   }
0244   if (nrQTests == 0)
0245     return -1;
0246   else if (nrFail == 0)
0247     return 1;
0248   else
0249     return 0;
0250 }
0251 
0252 void EgHLTOfflineSummaryClient::fillQTestData_(const edm::ParameterSet& iConfig,
0253                                                std::vector<SumHistBinData>& qTests,
0254                                                const std::string& label) {
0255   std::vector<edm::ParameterSet> qTestPara = iConfig.getParameter<std::vector<edm::ParameterSet> >(label);
0256   qTests.resize(qTestPara.size());
0257   for (size_t testNr = 0; testNr < qTestPara.size(); testNr++) {
0258     qTests[testNr].name = qTestPara[testNr].getParameter<std::string>("name");
0259     qTests[testNr].qTestPatterns = qTestPara[testNr].getParameter<std::vector<std::string> >("qTestsToCheck");
0260   }
0261 }
0262 
0263 // int EgHLTOfflineSummaryClient::eleTrigRelEffQTestResult_(const std::string& filterName)const
0264 // {
0265 
0266 // }
0267 
0268 // int EgHLTOfflineSummaryClient::phoTrigRelEffQTestResult_(const std::string& filterName)const
0269 // {
0270 
0271 // }
0272 
0273 // int EgHLTOfflineSummaryClient::eleTrigTPEffQTestResult_(const std::string& filterName)const
0274 // {
0275 
0276 // }
0277 
0278 // int EgHLTOfflineSummaryClient::trigEleQTestResult_(const std::string& filterName)const
0279 // {
0280 
0281 // }
0282 
0283 // int EgHLTOfflineSummaryClient::trigPhoQTestResult_(const std::string& filterName)const
0284 // {
0285 
0286 // }
0287 
0288 #include "FWCore/Framework/interface/MakerMacros.h"
0289 DEFINE_FWK_MODULE(EgHLTOfflineSummaryClient);