Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-11 02:22:01

0001 // Note to self: the implementation uses TH1F's to store the L1T and HLT rates.
0002 // Assuming a maximum rate of 100 kHz times a period of 23.31 s, one needs to store counts up to ~2.3e6.
0003 // A "float" has 24 bits of precision, so it can store up to 2**24 ~ 16.7e6 without loss of precision.
0004 
0005 // C++ headers
0006 #include <algorithm>
0007 #include <string>
0008 #include <vector>
0009 
0010 #include <fmt/printf.h>
0011 
0012 // boost headers
0013 #include <boost/regex.hpp>
0014 
0015 // CMSSW headers
0016 #include "CondFormats/DataRecord/interface/L1TUtmTriggerMenuRcd.h"
0017 #include "CondFormats/L1TObjects/interface/L1TUtmTriggerMenu.h"
0018 #include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
0019 #include "DataFormats/Common/interface/TriggerResults.h"
0020 #include "DataFormats/L1TGlobal/interface/GlobalAlgBlk.h"
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/Framework/interface/EventSetup.h"
0023 #include "FWCore/Framework/interface/LuminosityBlock.h"
0024 #include "FWCore/Framework/interface/Run.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0029 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0030 
0031 namespace {
0032 
0033   using dqm::reco::MonitorElement;
0034 
0035   struct RunBasedHistograms {
0036     // HLT configuration
0037     struct HLTIndices {
0038       unsigned int index_l1_seed;
0039       unsigned int index_prescale;
0040 
0041       HLTIndices() : index_l1_seed((unsigned int)-1), index_prescale((unsigned int)-1) {}
0042     };
0043 
0044     HLTConfigProvider hltConfig;
0045     std::vector<HLTIndices> hltIndices;
0046 
0047     std::vector<std::vector<unsigned int>> datasets;
0048     std::vector<std::vector<unsigned int>> streams;
0049 
0050     // L1T and HLT rate plots
0051 
0052     // per-path HLT plots
0053     struct HLTRatesPlots {
0054       MonitorElement *pass_l1_seed;
0055       MonitorElement *pass_prescale;
0056       MonitorElement *accept;
0057       MonitorElement *reject;
0058       MonitorElement *error;
0059     };
0060 
0061     // overall event count and event types
0062     MonitorElement *events_processed;
0063     std::vector<MonitorElement *> tcds_counts;
0064 
0065     // L1T triggers
0066     std::vector<MonitorElement *> l1t_counts;
0067 
0068     // HLT triggers
0069     std::vector<std::vector<HLTRatesPlots>> hlt_by_dataset_counts;
0070 
0071     // datasets
0072     std::vector<MonitorElement *> dataset_counts;
0073 
0074     // streams
0075     std::vector<MonitorElement *> stream_counts;
0076 
0077     RunBasedHistograms()
0078         :  // L1T and HLT configuration
0079           hltConfig(),
0080           hltIndices(),
0081           datasets(),
0082           streams(),
0083           // overall event count and event types
0084           events_processed(),
0085           tcds_counts(),
0086           // L1T triggers
0087           l1t_counts(),
0088           // HLT triggers
0089           hlt_by_dataset_counts(),
0090           // datasets
0091           dataset_counts(),
0092           // streams
0093           stream_counts() {}
0094   };
0095 }  // namespace
0096 
0097 class TriggerRatesMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
0098 public:
0099   explicit TriggerRatesMonitor(edm::ParameterSet const &);
0100   ~TriggerRatesMonitor() override = default;
0101 
0102   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0103 
0104 private:
0105   void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override;
0106   void bookHistograms(DQMStore::IBooker &,
0107                       edm::Run const &,
0108                       edm::EventSetup const &,
0109                       RunBasedHistograms &) const override;
0110   void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override;
0111 
0112   // TCDS trigger types
0113   // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
0114   static constexpr const char *const s_tcds_trigger_types[] = {
0115       "Empty",          //  0 - No trigger
0116       "Physics",        //  1 - GT trigger
0117       "Calibration",    //  2 - Sequence trigger (calibration)
0118       "Random",         //  3 - Random trigger
0119       "Auxiliary",      //  4 - Auxiliary (CPM front panel NIM input) trigger
0120       nullptr,          //  5 - reserved
0121       nullptr,          //  6 - reserved
0122       nullptr,          //  7 - reserved
0123       "Cyclic",         //  8 - Cyclic trigger
0124       "Bunch-pattern",  //  9 - Bunch-pattern trigger
0125       "Software",       // 10 - Software trigger
0126       "TTS",            // 11 - TTS-sourced trigger
0127       nullptr,          // 12 - reserved
0128       nullptr,          // 13 - reserved
0129       nullptr,          // 14 - reserved
0130       nullptr           // 15 - reserved
0131   };
0132 
0133   // module configuration
0134   const edm::ESGetToken<L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd> m_l1tMenu_token;
0135   const edm::InputTag m_l1t_results_inputTag;
0136   const edm::EDGetTokenT<GlobalAlgBlkBxCollection> m_l1t_results_token;
0137   const edm::EDGetTokenT<edm::TriggerResults> m_hlt_results_token;
0138   const std::string m_dqm_path;
0139   const uint32_t m_lumisections_range;
0140 };
0141 
0142 // definition
0143 constexpr const char *const TriggerRatesMonitor::s_tcds_trigger_types[];
0144 
0145 void TriggerRatesMonitor::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0146   edm::ParameterSetDescription desc;
0147   desc.addUntracked<edm::InputTag>("l1tResults", edm::InputTag("gtStage2Digis"));
0148   desc.addUntracked<edm::InputTag>("hltResults", edm::InputTag("TriggerResults"));
0149   desc.addUntracked<std::string>("dqmPath", "HLT/TriggerRates");
0150   desc.addUntracked<uint32_t>("lumisectionRange", 2500);  // ~16 hours
0151   descriptions.add("triggerRatesMonitor", desc);
0152 }
0153 
0154 TriggerRatesMonitor::TriggerRatesMonitor(edm::ParameterSet const &config)
0155     :  // module configuration
0156       m_l1tMenu_token{esConsumes<edm::Transition::BeginRun>()},
0157       m_l1t_results_inputTag{config.getUntrackedParameter<edm::InputTag>("l1tResults")},
0158       m_l1t_results_token{consumes(m_l1t_results_inputTag)},
0159       m_hlt_results_token{consumes(config.getUntrackedParameter<edm::InputTag>("hltResults"))},
0160       m_dqm_path{config.getUntrackedParameter<std::string>("dqmPath")},
0161       m_lumisections_range{config.getUntrackedParameter<uint32_t>("lumisectionRange")} {}
0162 
0163 void TriggerRatesMonitor::dqmBeginRun(edm::Run const &run,
0164                                       edm::EventSetup const &setup,
0165                                       RunBasedHistograms &histograms) const {
0166   histograms.tcds_counts.clear();
0167   histograms.tcds_counts.resize(sizeof(s_tcds_trigger_types) / sizeof(const char *));
0168 
0169   // cache the L1T trigger menu
0170   histograms.l1t_counts.clear();
0171   histograms.l1t_counts.resize(GlobalAlgBlk::maxPhysicsTriggers);
0172 
0173   // initialise the HLTConfigProvider
0174   bool changed = true;
0175   edm::EDConsumerBase::Labels labels;
0176   labelsForToken(m_hlt_results_token, labels);
0177   if (histograms.hltConfig.init(run, setup, labels.process, changed)) {
0178     // number of trigger paths in labels.process
0179     auto const nTriggers = histograms.hltConfig.size();
0180     histograms.hltIndices.resize(nTriggers);
0181 
0182     unsigned int const nDatasets = histograms.hltConfig.datasetNames().size();
0183     histograms.hlt_by_dataset_counts.clear();
0184     histograms.hlt_by_dataset_counts.resize(nDatasets);
0185 
0186     histograms.datasets.clear();
0187     histograms.datasets.resize(nDatasets);
0188     for (unsigned int i = 0; i < nDatasets; ++i) {
0189       auto const &paths = histograms.hltConfig.datasetContent(i);
0190       histograms.hlt_by_dataset_counts[i].resize(paths.size());
0191       histograms.datasets[i].reserve(paths.size());
0192       for (auto const &path : paths) {
0193         auto const triggerIdx = histograms.hltConfig.triggerIndex(path);
0194         if (triggerIdx < nTriggers)
0195           histograms.datasets[i].push_back(triggerIdx);
0196         else
0197           LogDebug("TriggerRatesMonitor")
0198               << "The rates of the HLT path \"" << path << "\" (dataset: \"" << histograms.hltConfig.datasetName(i)
0199               << "\") will not be monitored for this run.\nThis HLT path is not available in the process \""
0200               << labels.process << "\", but it is listed in its \"datasets\" PSet.";
0201       }
0202     }
0203     histograms.dataset_counts.clear();
0204     histograms.dataset_counts.resize(nDatasets);
0205 
0206     unsigned int const nStreams = histograms.hltConfig.streamNames().size();
0207     histograms.streams.clear();
0208     histograms.streams.resize(nStreams);
0209     for (unsigned int i = 0; i < nStreams; ++i) {
0210       for (auto const &dataset : histograms.hltConfig.streamContent(i)) {
0211         for (auto const &path : histograms.hltConfig.datasetContent(dataset)) {
0212           auto const triggerIdx = histograms.hltConfig.triggerIndex(path);
0213           if (triggerIdx < nTriggers)
0214             histograms.streams[i].push_back(triggerIdx);
0215           else
0216             LogDebug("TriggerRatesMonitor")
0217                 << "The rates of the HLT path \"" << path << "\" (stream: \"" << histograms.hltConfig.streamName(i)
0218                 << "\", dataset: \"" << dataset << "\") will not be monitored for this run.\n"
0219                 << "This HLT path is not available in the process \"" << labels.process
0220                 << "\", but it is listed in its \"datasets\" PSet.";
0221         }
0222       }
0223       std::sort(histograms.streams[i].begin(), histograms.streams[i].end());
0224       auto unique_end = std::unique(histograms.streams[i].begin(), histograms.streams[i].end());
0225       histograms.streams[i].resize(unique_end - histograms.streams[i].begin());
0226       histograms.streams[i].shrink_to_fit();
0227     }
0228     histograms.stream_counts.clear();
0229     histograms.stream_counts.resize(nStreams);
0230   } else {
0231     // HLTConfigProvider not initialised, skip the the HLT monitoring
0232     edm::LogError("TriggerRatesMonitor") << "Failed to initialise HLTConfigProvider: the rates of HLT triggers, "
0233                                             "datasets and streams will not be monitored for this run.";
0234   }
0235 }
0236 
0237 void TriggerRatesMonitor::bookHistograms(DQMStore::IBooker &booker,
0238                                          edm::Run const &run,
0239                                          edm::EventSetup const &setup,
0240                                          RunBasedHistograms &histograms) const {
0241   // book histograms for the overall event count, and trigger types
0242   booker.setCurrentFolder(m_dqm_path);
0243   histograms.events_processed = booker.book1D(
0244       "events", "Processed events vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0245   booker.setCurrentFolder(m_dqm_path + "/TCDS");
0246   unsigned int const sizeof_tcds_trigger_types = sizeof(s_tcds_trigger_types) / sizeof(const char *);
0247   if (sizeof_tcds_trigger_types == histograms.tcds_counts.size()) {
0248     for (unsigned int i = 0; i < sizeof_tcds_trigger_types; ++i)
0249       if (s_tcds_trigger_types[i]) {
0250         std::string const &title = fmt::sprintf("%s events vs. lumisection", s_tcds_trigger_types[i]);
0251         histograms.tcds_counts[i] =
0252             booker.book1D(s_tcds_trigger_types[i], title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0253       }
0254   } else
0255     edm::LogError("TriggerRatesMonitor")
0256         << "This should never happen: size of \"s_tcds_trigger_types\" array (" << sizeof_tcds_trigger_types
0257         << ") differs from size of \"histograms.tcds_counts\" vector (size=" << histograms.tcds_counts.size()
0258         << ").\nRate histograms of TCDS trigger types will not be booked for this run.";
0259 
0260   // book the rate histograms for the L1T triggers that are included in the L1T menu
0261   booker.setCurrentFolder(m_dqm_path + "/L1T");
0262   auto const &l1tMenu = setup.getData(m_l1tMenu_token);
0263   for (auto const &keyval : l1tMenu.getAlgorithmMap()) {
0264     unsigned int const bit = keyval.second.getIndex();
0265     if (bit >= histograms.l1t_counts.size()) {
0266       edm::LogError("TriggerRatesMonitor")
0267           << "This should never happen: bit of L1T algorithm (bit=" << bit << ", name=\"" << keyval.first
0268           << "\") is not smaller than size of \"histograms.l1t_counts\" vector (size=" << histograms.l1t_counts.size()
0269           << ").\nRate histogram of this L1T algorithm will not be booked for this run.";
0270       continue;
0271     }
0272     bool masked = false;  // FIXME read L1T masks once they will be avaiable in the EventSetup
0273     std::string const &name = fmt::sprintf("%s (bit %d)", keyval.first, bit);
0274     std::string const &title =
0275         fmt::sprintf("%s (bit %d)%s vs. lumisection", keyval.first, bit, (masked ? " (masked)" : ""));
0276     histograms.l1t_counts[bit] = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0277   }
0278 
0279   if (histograms.hltConfig.inited()) {
0280     auto const &datasetNames = histograms.hltConfig.datasetNames();
0281 
0282     // book the rate histograms for the HLT triggers
0283     for (unsigned int d = 0; d < datasetNames.size(); ++d) {
0284       booker.setCurrentFolder(m_dqm_path + "/HLT/" + datasetNames[d]);
0285       for (unsigned int i = 0; i < histograms.datasets[d].size(); ++i) {
0286         unsigned int index = histograms.datasets[d][i];
0287         std::string const &name = histograms.hltConfig.triggerName(index);
0288         histograms.hlt_by_dataset_counts[d][i].pass_l1_seed = booker.book1D(name + "_pass_L1_seed",
0289                                                                             name + " pass L1 seed, vs. lumisection",
0290                                                                             m_lumisections_range + 1,
0291                                                                             -0.5,
0292                                                                             m_lumisections_range + 0.5);
0293         histograms.hlt_by_dataset_counts[d][i].pass_prescale = booker.book1D(name + "_pass_prescaler",
0294                                                                              name + " pass prescaler, vs. lumisection",
0295                                                                              m_lumisections_range + 1,
0296                                                                              -0.5,
0297                                                                              m_lumisections_range + 0.5);
0298         histograms.hlt_by_dataset_counts[d][i].accept = booker.book1D(name + "_accept",
0299                                                                       name + " accept, vs. lumisection",
0300                                                                       m_lumisections_range + 1,
0301                                                                       -0.5,
0302                                                                       m_lumisections_range + 0.5);
0303         histograms.hlt_by_dataset_counts[d][i].reject = booker.book1D(name + "_reject",
0304                                                                       name + " reject, vs. lumisection",
0305                                                                       m_lumisections_range + 1,
0306                                                                       -0.5,
0307                                                                       m_lumisections_range + 0.5);
0308         histograms.hlt_by_dataset_counts[d][i].error = booker.book1D(name + "_error",
0309                                                                      name + " error, vs. lumisection",
0310                                                                      m_lumisections_range + 1,
0311                                                                      -0.5,
0312                                                                      m_lumisections_range + 0.5);
0313       }
0314 
0315       for (unsigned int i : histograms.datasets[d]) {
0316         // look for the index of the (last) L1T seed and prescale module in each path
0317         histograms.hltIndices[i].index_l1_seed = histograms.hltConfig.size(i);
0318         histograms.hltIndices[i].index_prescale = histograms.hltConfig.size(i);
0319         for (unsigned int j = 0; j < histograms.hltConfig.size(i); ++j) {
0320           std::string const &label = histograms.hltConfig.moduleLabel(i, j);
0321           std::string const &type = histograms.hltConfig.moduleType(label);
0322           if (type == "HLTL1TSeed" or type == "HLTLevel1GTSeed" or type == "HLTLevel1Activity" or
0323               type == "HLTLevel1Pattern") {
0324             // there might be more L1T seed filters in sequence
0325             // keep looking and store the index of the last one
0326             histograms.hltIndices[i].index_l1_seed = j;
0327           } else if (type == "HLTPrescaler") {
0328             // there should be only one prescaler in a path, and it should follow all L1T seed filters
0329             histograms.hltIndices[i].index_prescale = j;
0330             break;
0331           }
0332         }
0333       }
0334     }
0335 
0336     // book the rate histograms for the HLT datasets
0337     booker.setCurrentFolder(m_dqm_path + "/Datasets");
0338     for (unsigned int i = 0; i < datasetNames.size(); ++i)
0339       histograms.dataset_counts[i] =
0340           booker.book1D(datasetNames[i], datasetNames[i], m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0341 
0342     // book the rate histograms for the HLT streams
0343     booker.setCurrentFolder(m_dqm_path + "/Streams");
0344     auto const &streamNames = histograms.hltConfig.streamNames();
0345     for (unsigned int i = 0; i < streamNames.size(); ++i)
0346       histograms.stream_counts[i] =
0347           booker.book1D(streamNames[i], streamNames[i], m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
0348   }
0349 }
0350 
0351 void TriggerRatesMonitor::dqmAnalyze(edm::Event const &event,
0352                                      edm::EventSetup const &setup,
0353                                      RunBasedHistograms const &histograms) const {
0354   unsigned int lumisection = event.luminosityBlock();
0355 
0356   // monitor the overall event count and event types rates
0357   histograms.events_processed->Fill(lumisection);
0358   if (histograms.tcds_counts[event.experimentType()])
0359     histograms.tcds_counts[event.experimentType()]->Fill(lumisection);
0360 
0361   // monitor the rates of L1T triggers
0362   auto const &algBlkBxVecHandle = event.getHandle(m_l1t_results_token);
0363   if (not algBlkBxVecHandle.isValid()) {
0364     edm::LogError("TriggerRatesMonitor")
0365         << "L1 trigger results with label [" << m_l1t_results_inputTag.encode()
0366         << "] not present or invalid. MonitorElements of L1T results not filled for this event.";
0367   } else if (algBlkBxVecHandle->isEmpty(0)) {
0368     edm::LogError("TriggerRatesMonitor")
0369         << "L1 trigger results with label [" << m_l1t_results_inputTag.encode()
0370         << "] empty for BX=0. MonitorElements of L1T results not filled for this event.";
0371   } else {
0372     auto const &results = algBlkBxVecHandle->at(0, 0);
0373     for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i)
0374       if (results.getAlgoDecisionFinal(i))
0375         if (histograms.l1t_counts[i])
0376           histograms.l1t_counts[i]->Fill(lumisection);
0377   }
0378 
0379   // monitor the rates of HLT triggers, datasets and streams
0380   if (histograms.hltConfig.inited()) {
0381     auto const &hltResults = event.get(m_hlt_results_token);
0382     if (hltResults.size() != histograms.hltIndices.size()) {
0383       edm::LogError("TriggerRatesMonitor")
0384           << "This should never happen: the number of HLT paths has changed since the beginning of the run"
0385           << " (from " << histograms.hltIndices.size() << " to " << hltResults.size() << ").\n"
0386           << "Histograms for rates of HLT paths, datasets and streams will not be filled for this event.";
0387       return;
0388     }
0389 
0390     for (unsigned int d = 0; d < histograms.datasets.size(); ++d) {
0391       for (unsigned int i : histograms.datasets[d])
0392         if (hltResults[i].accept()) {
0393           histograms.dataset_counts[d]->Fill(lumisection);
0394           // ensure each dataset is incremented only once per event
0395           break;
0396         }
0397       for (unsigned int i = 0; i < histograms.datasets[d].size(); ++i) {
0398         unsigned int const index = histograms.datasets[d][i];
0399         edm::HLTPathStatus const &path = hltResults[index];
0400 
0401         if (path.index() > histograms.hltIndices[index].index_l1_seed)
0402           histograms.hlt_by_dataset_counts[d][i].pass_l1_seed->Fill(lumisection);
0403         if (path.index() > histograms.hltIndices[index].index_prescale)
0404           histograms.hlt_by_dataset_counts[d][i].pass_prescale->Fill(lumisection);
0405         if (path.accept())
0406           histograms.hlt_by_dataset_counts[d][i].accept->Fill(lumisection);
0407         else if (path.error())
0408           histograms.hlt_by_dataset_counts[d][i].error->Fill(lumisection);
0409         else
0410           histograms.hlt_by_dataset_counts[d][i].reject->Fill(lumisection);
0411       }
0412     }
0413 
0414     for (unsigned int i = 0; i < histograms.streams.size(); ++i)
0415       for (unsigned int j : histograms.streams[i])
0416         if (hltResults[j].accept()) {
0417           histograms.stream_counts[i]->Fill(lumisection);
0418           // ensure each stream is incremented only once per event
0419           break;
0420         }
0421   }
0422 }
0423 
0424 // define this as a plugin
0425 #include "FWCore/Framework/interface/MakerMacros.h"
0426 DEFINE_FWK_MODULE(TriggerRatesMonitor);