Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-09-30 03:09:08

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/ServiceRegistry/interface/Service.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0010 
0011 #include "DataFormats/Common/interface/TriggerResults.h"
0012 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0013 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
0014 
0015 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0016 
0017 #include <string>
0018 #include <vector>
0019 #include <map>
0020 #include <algorithm>
0021 
0022 class HLTFiltersDQMonitor : public DQMEDAnalyzer {
0023 public:
0024   explicit HLTFiltersDQMonitor(const edm::ParameterSet&);
0025   ~HLTFiltersDQMonitor() override = default;
0026   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0027 
0028 private:
0029   void dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0030   void bookHistograms(DQMStore::IBooker&, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0031   void analyze(const edm::Event&, const edm::EventSetup&) override;
0032 
0033   bool skipStreamByName(const std::string& streamName) const;
0034   bool skipPathMonitorElement(const std::string& pathName) const;
0035   bool skipModuleByEDMType(const std::string& moduleEDMType) const;
0036   bool skipModuleByType(const std::string& moduleType) const;
0037 
0038   const std::string folderName_;
0039   const std::string efficPlotNamePrefix_;
0040   std::string processName_;
0041   bool initFailed_;
0042   bool skipRun_;
0043 
0044   MonitorElement* meMenu_;
0045   std::map<std::string, MonitorElement*> meDatasetMap_;
0046   std::map<std::string, MonitorElement*> mePathMap_;
0047 
0048   edm::EDGetTokenT<edm::TriggerResults> triggerResultsToken_;
0049   edm::EDGetTokenT<trigger::TriggerEvent> triggerSummaryTokenAOD_;
0050   edm::EDGetTokenT<trigger::TriggerEventWithRefs> triggerSummaryTokenRAW_;
0051 
0052   HLTConfigProvider hltConfigProvider_;
0053 };
0054 
0055 HLTFiltersDQMonitor::HLTFiltersDQMonitor(const edm::ParameterSet& iConfig)
0056     : folderName_(iConfig.getParameter<std::string>("folderName")),
0057       efficPlotNamePrefix_(iConfig.getParameter<std::string>("efficPlotNamePrefix")),
0058       processName_(""),
0059       initFailed_(false),
0060       skipRun_(false),
0061       meMenu_(nullptr) {
0062   auto const triggerResultsInputTag(iConfig.getParameter<edm::InputTag>("triggerResults"));
0063 
0064   if (triggerResultsInputTag.process().empty()) {
0065     edm::LogError("Input") << "process not specified in HLT TriggerResults InputTag \""
0066                            << triggerResultsInputTag.encode() << "\" -> plugin will not produce DQM outputs";
0067     initFailed_ = true;
0068     return;
0069   } else {
0070     processName_ = triggerResultsInputTag.process();
0071 
0072     triggerResultsToken_ = consumes<edm::TriggerResults>(triggerResultsInputTag);
0073 
0074     auto triggerSummaryAODInputTag(iConfig.getParameter<edm::InputTag>("triggerSummaryAOD"));
0075     if (triggerSummaryAODInputTag.process().empty()) {
0076       triggerSummaryAODInputTag =
0077           edm::InputTag(triggerSummaryAODInputTag.label(), triggerSummaryAODInputTag.instance(), processName_);
0078     } else if (triggerSummaryAODInputTag.process() != processName_) {
0079       edm::LogWarning("Input") << "edm::TriggerResults process name '" << processName_
0080                                << "' differs from trigger::TriggerEvent process name '"
0081                                << triggerSummaryAODInputTag.process() << "' -> plugin will not produce DQM outputs";
0082       initFailed_ = true;
0083       return;
0084     }
0085     triggerSummaryTokenAOD_ = consumes<trigger::TriggerEvent>(triggerSummaryAODInputTag);
0086 
0087     auto triggerSummaryRAWInputTag(iConfig.getParameter<edm::InputTag>("triggerSummaryRAW"));
0088     if (triggerSummaryRAWInputTag.process().empty()) {
0089       triggerSummaryRAWInputTag =
0090           edm::InputTag(triggerSummaryRAWInputTag.label(), triggerSummaryRAWInputTag.instance(), processName_);
0091     } else if (triggerSummaryRAWInputTag.process() != processName_) {
0092       edm::LogWarning("Input") << "edm::TriggerResults process name '" << processName_
0093                                << "' differs from trigger::TriggerEventWithRefs process name '"
0094                                << triggerSummaryRAWInputTag.process() << "' -> plugin will not produce DQM outputs";
0095       initFailed_ = true;
0096       return;
0097     }
0098     triggerSummaryTokenRAW_ = mayConsume<trigger::TriggerEventWithRefs>(triggerSummaryRAWInputTag);
0099   }
0100 }
0101 
0102 void HLTFiltersDQMonitor::dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
0103   if (initFailed_) {
0104     return;
0105   }
0106 
0107   LogTrace("")
0108       << "[HLTFiltersDQMonitor] "
0109       << "----------------------------------------------------------------------------------------------------";
0110   LogTrace("") << "[HLTFiltersDQMonitor::dqmBeginRun] Run = " << iRun.id();
0111 
0112   // reset data members holding information from the previous run
0113   skipRun_ = false;
0114 
0115   bool hltChanged(true);
0116   if (hltConfigProvider_.init(iRun, iSetup, processName_, hltChanged)) {
0117     LogTrace("") << "[HLTFiltersDQMonitor::dqmBeginRun] HLTConfigProvider initialized [processName() = "
0118                  << hltConfigProvider_.processName() << ", tableName() = " << hltConfigProvider_.tableName()
0119                  << ", size() = " << hltConfigProvider_.size() << "]";
0120   } else {
0121     edm::LogError("Input") << "initialization of HLTConfigProvider failed for Run=" << iRun.id() << " (process=\""
0122                            << processName_ << "\") -> plugin will not produce DQM outputs for this run";
0123     skipRun_ = true;
0124     return;
0125   }
0126 }
0127 
0128 void HLTFiltersDQMonitor::bookHistograms(DQMStore::IBooker& iBooker,
0129                                          edm::Run const& iRun,
0130                                          edm::EventSetup const& iSetup) {
0131   if (skipRun_ or initFailed_) {
0132     return;
0133   }
0134 
0135   iBooker.setCurrentFolder(folderName_);
0136 
0137   iBooker.bookString("HLTMenu", hltConfigProvider_.tableName().c_str());
0138 
0139   auto hltMenuName(hltConfigProvider_.tableName());
0140   std::replace(hltMenuName.begin(), hltMenuName.end(), '/', '_');
0141   std::replace(hltMenuName.begin(), hltMenuName.end(), '.', '_');
0142   while (hltMenuName.front() == '_') {
0143     hltMenuName.erase(0, 1);
0144   }
0145 
0146   auto const& triggerNames(hltConfigProvider_.triggerNames());
0147 
0148   meMenu_ = iBooker.bookProfile(efficPlotNamePrefix_ + hltMenuName,
0149                                 "Path Efficiency",
0150                                 triggerNames.size(),
0151                                 0.,
0152                                 triggerNames.size(),
0153                                 -0.1,
0154                                 1.1,
0155                                 "");
0156   if (meMenu_ and meMenu_->getTProfile() and meMenu_->getTProfile()->GetXaxis()) {
0157     for (size_t idx = 0; idx < triggerNames.size(); ++idx) {
0158       meMenu_->getTProfile()->GetXaxis()->SetBinLabel(idx + 1, triggerNames.at(idx).c_str());
0159     }
0160   }
0161 
0162   LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms] HLTConfigProvider::size() = " << hltConfigProvider_.size()
0163                << ", HLTConfigProvider::triggerNames().size() = " << triggerNames.size();
0164 
0165   for (auto const& istream : hltConfigProvider_.streamNames()) {
0166     LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms] Stream = \"" << istream << "\"";
0167 
0168     if (not this->skipStreamByName(istream)) {
0169       auto const& dsets(hltConfigProvider_.streamContent(istream));
0170       for (auto const& idset : dsets) {
0171         const std::vector<std::string>& dsetPathNames = hltConfigProvider_.datasetContent(idset);
0172         iBooker.setCurrentFolder(folderName_ + "/" + idset);
0173         LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms]   Dataset = \"" << idset << "\"";
0174         const std::string meDatasetName(efficPlotNamePrefix_ + idset);
0175         meDatasetMap_[meDatasetName] = iBooker.bookProfile(
0176             meDatasetName.c_str(), meDatasetName.c_str(), dsetPathNames.size(), 0., dsetPathNames.size(), -0.1, 1.1, "");
0177         TProfile* meDatasetTProf(nullptr);
0178         if (meDatasetMap_.at(meDatasetName)) {
0179           meDatasetTProf = meDatasetMap_.at(meDatasetName)->getTProfile();
0180         }
0181         for (size_t idxPath = 0; idxPath < dsetPathNames.size(); ++idxPath) {
0182           auto const& iPathName(dsetPathNames.at(idxPath));
0183           if (std::find(triggerNames.begin(), triggerNames.end(), iPathName) == triggerNames.end()) {
0184             continue;
0185           }
0186           if (meDatasetTProf and meDatasetTProf->GetXaxis()) {
0187             meDatasetTProf->GetXaxis()->SetBinLabel(idxPath + 1, iPathName.c_str());
0188           }
0189           if (this->skipPathMonitorElement(iPathName)) {
0190             continue;
0191           }
0192 
0193           LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms]     Path = \"" << iPathName << "\"";
0194 
0195           auto const& moduleLabels(hltConfigProvider_.moduleLabels(iPathName));
0196           std::vector<std::string> mePath_binLabels;
0197           mePath_binLabels.reserve(moduleLabels.size());
0198           for (size_t iMod = 0; iMod < moduleLabels.size(); ++iMod) {
0199             auto const& moduleLabel(moduleLabels.at(iMod));
0200             if (this->skipModuleByEDMType(hltConfigProvider_.moduleEDMType(moduleLabel)) or
0201                 this->skipModuleByType(hltConfigProvider_.moduleType(moduleLabel))) {
0202               LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms]       [-] Module = \"" << moduleLabel << "\"";
0203               continue;
0204             }
0205             LogTrace("") << "[HLTFiltersDQMonitor::bookHistograms]       [bin=" << mePath_binLabels.size() + 1
0206                          << "] Module = \"" << moduleLabel << "\"";
0207 
0208             if (std::find(mePath_binLabels.begin(), mePath_binLabels.end(), moduleLabel) != mePath_binLabels.end()) {
0209               edm::LogInfo("Input") << "module \"" << moduleLabel << "\" included multiple times in path \""
0210                                     << iPathName << "\""
0211                                     << "-> only 1 bin labelled \"" << moduleLabel
0212                                     << "\" will be created in the DQM MonitorElement of the path";
0213               continue;
0214             }
0215 
0216             mePath_binLabels.emplace_back(moduleLabel);
0217           }
0218 
0219           if (mePath_binLabels.empty()) {
0220             continue;
0221           }
0222 
0223           const std::string mePathName(efficPlotNamePrefix_ + idset + "_" + iPathName);
0224           mePathMap_[mePathName] = iBooker.bookProfile(mePathName.c_str(),
0225                                                        iPathName.c_str(),
0226                                                        mePath_binLabels.size(),
0227                                                        0.,
0228                                                        mePath_binLabels.size(),
0229                                                        -0.1,
0230                                                        1.1,
0231                                                        "");
0232 
0233           if (mePathMap_.at(mePathName)) {
0234             auto* const mePathTProf(mePathMap_.at(mePathName)->getTProfile());
0235             if (mePathTProf and mePathTProf->GetXaxis()) {
0236               for (size_t iMod = 0; iMod < mePath_binLabels.size(); ++iMod) {
0237                 mePathTProf->GetXaxis()->SetBinLabel(iMod + 1, mePath_binLabels.at(iMod).c_str());
0238               }
0239             }
0240           }
0241         }
0242       }
0243     }
0244   }
0245 }
0246 
0247 void HLTFiltersDQMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0248   if (skipRun_ or initFailed_) {
0249     return;
0250   }
0251 
0252   LogTrace("") << "[HLTFiltersDQMonitor::analyze] --------------------------------------------------------";
0253   LogTrace("") << "[HLTFiltersDQMonitor::analyze] Run = " << iEvent.id().run()
0254                << ", LuminosityBlock = " << iEvent.id().luminosityBlock() << ", Event = " << iEvent.id().event();
0255 
0256   auto const& triggerResults(iEvent.getHandle(triggerResultsToken_));
0257 
0258   if (not triggerResults.isValid()) {
0259     edm::LogWarning("Input") << "invalid handle to edm::TriggerResults (InputTag: \"triggerResults\")"
0260                              << " -> plugin will not fill DQM outputs for this event";
0261     return;
0262   }
0263 
0264   // fill MonitorElement: HLT-Menu (bin: path)
0265   if (meMenu_ and meMenu_->getTProfile() and meMenu_->getTProfile()->GetXaxis()) {
0266     auto const& triggerNames(hltConfigProvider_.triggerNames());
0267     for (auto const& iPathName : triggerNames) {
0268       const uint pathIndex(hltConfigProvider_.triggerIndex(iPathName));
0269       if (pathIndex >= triggerResults->size()) {
0270         edm::LogError("Logic") << "[HLTFiltersDQMonitor::analyze]       "
0271                                << "index associated to path \"" << iPathName << "\" (" << pathIndex
0272                                << ") is inconsistent with triggerResults::size() (" << triggerResults->size()
0273                                << ") -> plugin will not fill bin associated to this path in HLT-Menu MonitorElement";
0274         continue;
0275       }
0276       auto const pathAccept(triggerResults->accept(pathIndex));
0277       auto* const axis(meMenu_->getTProfile()->GetXaxis());
0278       auto const ibin(axis->FindBin(iPathName.c_str()));
0279       if ((0 < ibin) and (ibin <= axis->GetNbins())) {
0280         meMenu_->Fill(ibin - 0.5, pathAccept);
0281       }
0282     }
0283   }
0284 
0285   auto const& triggerEventAOD(iEvent.getHandle(triggerSummaryTokenAOD_));
0286   edm::Handle<trigger::TriggerEventWithRefs> triggerEventRAW;
0287 
0288   bool useTriggerEventAOD(true);
0289   if (not triggerEventAOD.isValid()) {
0290     useTriggerEventAOD = false;
0291     edm::LogInfo("Input") << "invalid handle to trigger::TriggerEvent (InputTag: \"triggerSummaryAOD\"),"
0292                           << " will attempt to access trigger::TriggerEventWithRefs (InputTag: \"triggerSummaryRAW\")";
0293 
0294     triggerEventRAW = iEvent.getHandle(triggerSummaryTokenRAW_);
0295     if (not triggerEventRAW.isValid()) {
0296       edm::LogWarning("Input") << "invalid handle to trigger::TriggerEventWithRefs (InputTag: \"triggerSummaryRAW\")"
0297                                << " -> plugin will not fill DQM outputs for this event";
0298       return;
0299     }
0300   }
0301 
0302   auto const triggerEventSize(useTriggerEventAOD ? triggerEventAOD->sizeFilters() : triggerEventRAW->size());
0303   LogTrace("") << "[HLTFiltersDQMonitor::analyze] useTriggerEventAOD = " << useTriggerEventAOD
0304                << ", triggerEventSize = " << triggerEventSize;
0305 
0306   // fill MonitorElements for PrimaryDatasets and Paths
0307   // loop over Streams
0308   for (auto const& istream : hltConfigProvider_.streamNames()) {
0309     LogTrace("") << "[HLTFiltersDQMonitor::analyze]   Stream = \"" << istream << "\"";
0310     auto const& dsets(hltConfigProvider_.streamContent(istream));
0311     // loop over PrimaryDatasets in Stream
0312     for (auto const& idset : dsets) {
0313       LogTrace("") << "[HLTFiltersDQMonitor::analyze]     Dataset = \"" << idset << "\"";
0314       TProfile* meDatasetTProf(nullptr);
0315       const std::string meDatasetName(efficPlotNamePrefix_ + idset);
0316       if (meDatasetMap_.find(meDatasetName) != meDatasetMap_.end()) {
0317         meDatasetTProf = meDatasetMap_.at(meDatasetName)->getTProfile();
0318       }
0319       auto const& dsetPathNames(hltConfigProvider_.datasetContent(idset));
0320       // loop over Paths in PrimaryDataset
0321       for (auto const& iPathName : dsetPathNames) {
0322         const uint pathIndex(hltConfigProvider_.triggerIndex(iPathName));
0323         if (pathIndex >= triggerResults->size()) {
0324           edm::LogError("Logic") << "[HLTFiltersDQMonitor::analyze]       "
0325                                  << "index associated to path \"" << iPathName << "\" (" << pathIndex
0326                                  << ") is inconsistent with triggerResults::size() (" << triggerResults->size()
0327                                  << ") -> plugin will not fill DQM info related to this path";
0328           continue;
0329         }
0330         auto const pathAccept(triggerResults->accept(pathIndex));
0331         LogTrace("") << "[HLTFiltersDQMonitor::analyze]       "
0332                      << "Path = \"" << iPathName << "\", HLTConfigProvider::triggerIndex(\"" << iPathName
0333                      << "\") = " << pathIndex << ", Accept = " << pathAccept;
0334         // fill MonitorElement: PrimaryDataset (bin: path)
0335         if (meDatasetTProf and meDatasetTProf->GetXaxis()) {
0336           auto* const axis(meDatasetTProf->GetXaxis());
0337           auto const ibin(axis->FindBin(iPathName.c_str()));
0338           if ((0 < ibin) and (ibin <= axis->GetNbins())) {
0339             meDatasetTProf->Fill(ibin - 0.5, pathAccept);
0340           }
0341         }
0342         // fill MonitorElement: Path (bin: filter)
0343         auto const mePathName(efficPlotNamePrefix_ + idset + "_" + iPathName);
0344         if (mePathMap_.find(mePathName) != mePathMap_.end()) {
0345           auto* const mePathTProf(mePathMap_.at(mePathName)->getTProfile());
0346           if (mePathTProf) {
0347             auto* const axis(mePathTProf->GetXaxis());
0348             if (axis) {
0349               unsigned indexLastFilterPathModules(triggerResults->index(pathIndex) + 1);
0350               LogTrace("") << "[HLTFiltersDQMonitor::analyze]         "
0351                            << "indexLastFilterPathModules = " << indexLastFilterPathModules;
0352               // identify module corresponding to last filter executed in the path
0353               while (indexLastFilterPathModules > 0) {
0354                 --indexLastFilterPathModules;
0355                 const std::string& labelLastFilterPathModules(
0356                     hltConfigProvider_.moduleLabel(pathIndex, indexLastFilterPathModules));
0357                 const uint indexLastFilterFilters =
0358                     useTriggerEventAOD
0359                         ? triggerEventAOD->filterIndex(edm::InputTag(labelLastFilterPathModules, "", processName_))
0360                         : triggerEventRAW->filterIndex(edm::InputTag(labelLastFilterPathModules, "", processName_));
0361                 LogTrace("") << "[HLTFiltersDQMonitor::analyze]           "
0362                              << "indexLastFilterPathModules = " << indexLastFilterPathModules
0363                              << ", labelLastFilterPathModules = " << labelLastFilterPathModules
0364                              << ", indexLastFilterFilters = " << indexLastFilterFilters
0365                              << " (triggerEventSize = " << triggerEventSize << ")";
0366                 if (indexLastFilterFilters < triggerEventSize) {
0367                   if (this->skipModuleByType(hltConfigProvider_.moduleType(labelLastFilterPathModules))) {
0368                     continue;
0369                   }
0370                   break;
0371                 }
0372               }
0373               // number of modules in the path
0374               const unsigned sizeModulesPath(hltConfigProvider_.size(pathIndex));
0375               LogTrace("") << "[HLTFiltersDQMonitor::analyze]         "
0376                            << "-> selected indexLastFilterPathModules = " << indexLastFilterPathModules
0377                            << " (HLTConfigProvider::size(" << pathIndex << ") = " << sizeModulesPath << ")";
0378               if (indexLastFilterPathModules >= sizeModulesPath) {
0379                 edm::LogError("Logic") << " selected index (" << indexLastFilterPathModules
0380                                        << ") for last filter of path \"" << iPathName
0381                                        << "\" is inconsistent with number of modules in the path (" << sizeModulesPath
0382                                        << ")";
0383                 continue;
0384               }
0385               // store decision of previous filter
0386               bool previousFilterAccept(true);
0387               for (size_t modIdx = 0; modIdx < sizeModulesPath; ++modIdx) {
0388                 // each filter-bin is filled, with a 0 or 1, only when all previous filters in the path have passed
0389                 if (not previousFilterAccept) {
0390                   break;
0391                 }
0392                 // consider only selected EDFilter modules
0393                 auto const& moduleLabel(hltConfigProvider_.moduleLabel(pathIndex, modIdx));
0394                 if (this->skipModuleByEDMType(hltConfigProvider_.moduleEDMType(moduleLabel)) or
0395                     this->skipModuleByType(hltConfigProvider_.moduleType(moduleLabel))) {
0396                   continue;
0397                 }
0398                 // index of the module in the path [0,sizeModulesPath)
0399                 const unsigned slotModule(hltConfigProvider_.moduleIndex(pathIndex, moduleLabel));
0400                 bool filterAccept(false);
0401                 if (slotModule < indexLastFilterPathModules) {
0402                   filterAccept = true;
0403                 } else if (slotModule == indexLastFilterPathModules) {
0404                   filterAccept = pathAccept;
0405                 }
0406                 LogTrace("") << "[HLTFiltersDQMonitor::analyze]         "
0407                              << "HLTConfigProvider::moduleLabel(" << pathIndex << ", " << modIdx << ") = \""
0408                              << moduleLabel << "\", HLTConfigProvider::moduleIndex(" << pathIndex << ", \""
0409                              << moduleLabel << "\") = " << slotModule << ", filterAccept = " << filterAccept
0410                              << ", previousFilterAccept = " << previousFilterAccept;
0411                 auto const ibin(axis->FindBin(moduleLabel.c_str()));
0412                 if ((0 < ibin) and (ibin <= axis->GetNbins())) {
0413                   mePathTProf->Fill(ibin - 0.5, filterAccept);
0414                 }
0415                 previousFilterAccept = filterAccept;
0416               }
0417             }
0418           }
0419         }
0420       }
0421     }
0422   }
0423 }
0424 
0425 bool HLTFiltersDQMonitor::skipStreamByName(const std::string& streamName) const {
0426   if ((streamName.find("Physics") != std::string::npos) || (streamName.find("Scouting") != std::string::npos) ||
0427       (streamName.find("Parking") != std::string::npos) || (streamName == "A")) {
0428     return false;
0429   }
0430   return true;
0431 }
0432 
0433 bool HLTFiltersDQMonitor::skipPathMonitorElement(const std::string& pathName) const {
0434   if ((pathName.find("HLT_") == std::string::npos) || (pathName.find("HLT_Physics") != std::string::npos) ||
0435       (pathName.find("HLT_Random") != std::string::npos)) {
0436     return true;
0437   }
0438   return false;
0439 }
0440 
0441 bool HLTFiltersDQMonitor::skipModuleByEDMType(const std::string& moduleEDMType) const {
0442   return (moduleEDMType != "EDFilter");
0443 }
0444 
0445 bool HLTFiltersDQMonitor::skipModuleByType(const std::string& moduleType) const { return (moduleType == "HLTBool"); }
0446 
0447 void HLTFiltersDQMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0448   edm::ParameterSetDescription desc;
0449   desc.add<std::string>("folderName", "HLT/Filters");
0450   desc.add<std::string>("efficPlotNamePrefix", "effic_");
0451   desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults::HLT"));
0452   desc.add<edm::InputTag>("triggerSummaryAOD", edm::InputTag("hltTriggerSummaryAOD::HLT"));
0453   desc.add<edm::InputTag>("triggerSummaryRAW", edm::InputTag("hltTriggerSummaryRAW::HLT"));
0454   descriptions.add("hltFiltersDQMonitor", desc);
0455 }
0456 
0457 DEFINE_FWK_MODULE(HLTFiltersDQMonitor);