CMS 3D CMS Logo

TriggerRatesMonitor.cc
Go to the documentation of this file.
1 // Note to self: the implementation uses TH1F's to store the L1T and HLT rates.
2 // Assuming a maximum rate of 100 kHz times a period of 23.31 s, one needs to store counts up to ~2.3e6.
3 // A "float" has 24 bits of precision, so it can store up to 2**24 ~ 16.7e6 without loss of precision.
4 
5 // C++ headers
6 #include <algorithm>
7 #include <string>
8 #include <vector>
9 
10 #include <fmt/printf.h>
11 
12 // boost headers
13 #include <boost/regex.hpp>
14 
15 // CMSSW headers
30 
31 namespace {
32 
34 
35  struct RunBasedHistograms {
36  // HLT configuration
37  struct HLTIndices {
38  unsigned int index_l1_seed;
39  unsigned int index_prescale;
40 
41  HLTIndices() : index_l1_seed((unsigned int)-1), index_prescale((unsigned int)-1) {}
42  };
43 
45  std::vector<HLTIndices> hltIndices;
46 
47  std::vector<std::vector<unsigned int>> datasets;
48  std::vector<std::vector<unsigned int>> streams;
49 
50  // L1T and HLT rate plots
51 
52  // per-path HLT plots
53  struct HLTRatesPlots {
54  MonitorElement *pass_l1_seed;
55  MonitorElement *pass_prescale;
59  };
60 
61  // overall event count and event types
62  MonitorElement *events_processed;
63  std::vector<MonitorElement *> tcds_counts;
64 
65  // L1T triggers
66  std::vector<MonitorElement *> l1t_counts;
67 
68  // HLT triggers
69  std::vector<std::vector<HLTRatesPlots>> hlt_by_dataset_counts;
70 
71  // datasets
72  std::vector<MonitorElement *> dataset_counts;
73 
74  // streams
75  std::vector<MonitorElement *> stream_counts;
76 
77  RunBasedHistograms()
78  : // L1T and HLT configuration
79  hltConfig(),
80  hltIndices(),
81  datasets(),
82  streams(),
83  // overall event count and event types
84  events_processed(),
85  tcds_counts(),
86  // L1T triggers
87  l1t_counts(),
88  // HLT triggers
89  hlt_by_dataset_counts(),
90  // datasets
91  dataset_counts(),
92  // streams
93  stream_counts() {}
94  };
95 } // namespace
96 
97 class TriggerRatesMonitor : public DQMGlobalEDAnalyzer<RunBasedHistograms> {
98 public:
99  explicit TriggerRatesMonitor(edm::ParameterSet const &);
100  ~TriggerRatesMonitor() override = default;
101 
102  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
103 
104 private:
105  void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override;
107  edm::Run const &,
108  edm::EventSetup const &,
109  RunBasedHistograms &) const override;
110  void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override;
111 
112  // TCDS trigger types
113  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
114  static constexpr const char *const s_tcds_trigger_types[] = {
115  "Empty", // 0 - No trigger
116  "Physics", // 1 - GT trigger
117  "Calibration", // 2 - Sequence trigger (calibration)
118  "Random", // 3 - Random trigger
119  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
120  nullptr, // 5 - reserved
121  nullptr, // 6 - reserved
122  nullptr, // 7 - reserved
123  "Cyclic", // 8 - Cyclic trigger
124  "Bunch-pattern", // 9 - Bunch-pattern trigger
125  "Software", // 10 - Software trigger
126  "TTS", // 11 - TTS-sourced trigger
127  nullptr, // 12 - reserved
128  nullptr, // 13 - reserved
129  nullptr, // 14 - reserved
130  nullptr // 15 - reserved
131  };
132 
133  // module configuration
139  const uint32_t m_lumisections_range;
140 };
141 
142 // definition
143 constexpr const char *const TriggerRatesMonitor::s_tcds_trigger_types[];
144 
147  desc.addUntracked<edm::InputTag>("l1tResults", edm::InputTag("gtStage2Digis"));
148  desc.addUntracked<edm::InputTag>("hltResults", edm::InputTag("TriggerResults"));
149  desc.addUntracked<std::string>("dqmPath", "HLT/TriggerRates");
150  desc.addUntracked<uint32_t>("lumisectionRange", 2500); // ~16 hours
151  descriptions.add("triggerRatesMonitor", desc);
152 }
153 
155  : // module configuration
156  m_l1tMenu_token{esConsumes<edm::Transition::BeginRun>()},
157  m_l1t_results_inputTag{config.getUntrackedParameter<edm::InputTag>("l1tResults")},
158  m_l1t_results_token{consumes(m_l1t_results_inputTag)},
159  m_hlt_results_token{consumes(config.getUntrackedParameter<edm::InputTag>("hltResults"))},
160  m_dqm_path{config.getUntrackedParameter<std::string>("dqmPath")},
161  m_lumisections_range{config.getUntrackedParameter<uint32_t>("lumisectionRange")} {}
162 
164  edm::EventSetup const &setup,
165  RunBasedHistograms &histograms) const {
166  histograms.tcds_counts.clear();
167  histograms.tcds_counts.resize(sizeof(s_tcds_trigger_types) / sizeof(const char *));
168 
169  // cache the L1T trigger menu
170  histograms.l1t_counts.clear();
171  histograms.l1t_counts.resize(GlobalAlgBlk::maxPhysicsTriggers);
172 
173  // initialise the HLTConfigProvider
174  bool changed = true;
177  if (histograms.hltConfig.init(run, setup, labels.process, changed)) {
178  // number of trigger paths in labels.process
179  auto const nTriggers = histograms.hltConfig.size();
180  histograms.hltIndices.resize(nTriggers);
181 
182  unsigned int const nDatasets = histograms.hltConfig.datasetNames().size();
183  histograms.hlt_by_dataset_counts.clear();
184  histograms.hlt_by_dataset_counts.resize(nDatasets);
185 
186  histograms.datasets.clear();
187  histograms.datasets.resize(nDatasets);
188  for (unsigned int i = 0; i < nDatasets; ++i) {
189  auto const &paths = histograms.hltConfig.datasetContent(i);
190  histograms.hlt_by_dataset_counts[i].resize(paths.size());
191  histograms.datasets[i].reserve(paths.size());
192  for (auto const &path : paths) {
193  auto const triggerIdx = histograms.hltConfig.triggerIndex(path);
194  if (triggerIdx < nTriggers)
195  histograms.datasets[i].push_back(triggerIdx);
196  else
197  LogDebug("TriggerRatesMonitor")
198  << "The rates of the HLT path \"" << path << "\" (dataset: \"" << histograms.hltConfig.datasetName(i)
199  << "\") will not be monitored for this run.\nThis HLT path is not available in the process \""
200  << labels.process << "\", but it is listed in its \"datasets\" PSet.";
201  }
202  }
203  histograms.dataset_counts.clear();
204  histograms.dataset_counts.resize(nDatasets);
205 
206  unsigned int const nStreams = histograms.hltConfig.streamNames().size();
207  histograms.streams.clear();
208  histograms.streams.resize(nStreams);
209  for (unsigned int i = 0; i < nStreams; ++i) {
210  for (auto const &dataset : histograms.hltConfig.streamContent(i)) {
211  for (auto const &path : histograms.hltConfig.datasetContent(dataset)) {
212  auto const triggerIdx = histograms.hltConfig.triggerIndex(path);
213  if (triggerIdx < nTriggers)
214  histograms.streams[i].push_back(triggerIdx);
215  else
216  LogDebug("TriggerRatesMonitor")
217  << "The rates of the HLT path \"" << path << "\" (stream: \"" << histograms.hltConfig.streamName(i)
218  << "\", dataset: \"" << dataset << "\") will not be monitored for this run.\n"
219  << "This HLT path is not available in the process \"" << labels.process
220  << "\", but it is listed in its \"datasets\" PSet.";
221  }
222  }
223  std::sort(histograms.streams[i].begin(), histograms.streams[i].end());
224  auto unique_end = std::unique(histograms.streams[i].begin(), histograms.streams[i].end());
225  histograms.streams[i].resize(unique_end - histograms.streams[i].begin());
226  histograms.streams[i].shrink_to_fit();
227  }
228  histograms.stream_counts.clear();
229  histograms.stream_counts.resize(nStreams);
230  } else {
231  // HLTConfigProvider not initialised, skip the the HLT monitoring
232  edm::LogError("TriggerRatesMonitor") << "Failed to initialise HLTConfigProvider: the rates of HLT triggers, "
233  "datasets and streams will not be monitored for this run.";
234  }
235 }
236 
238  edm::Run const &run,
239  edm::EventSetup const &setup,
240  RunBasedHistograms &histograms) const {
241  // book histograms for the overall event count, and trigger types
243  histograms.events_processed = booker.book1D(
244  "events", "Processed events vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
245  booker.setCurrentFolder(m_dqm_path + "/TCDS");
246  unsigned int const sizeof_tcds_trigger_types = sizeof(s_tcds_trigger_types) / sizeof(const char *);
247  if (sizeof_tcds_trigger_types == histograms.tcds_counts.size()) {
248  for (unsigned int i = 0; i < sizeof_tcds_trigger_types; ++i)
249  if (s_tcds_trigger_types[i]) {
250  std::string const &title = fmt::sprintf("%s events vs. lumisection", s_tcds_trigger_types[i]);
251  histograms.tcds_counts[i] =
253  }
254  } else
255  edm::LogError("TriggerRatesMonitor")
256  << "This should never happen: size of \"s_tcds_trigger_types\" array (" << sizeof_tcds_trigger_types
257  << ") differs from size of \"histograms.tcds_counts\" vector (size=" << histograms.tcds_counts.size()
258  << ").\nRate histograms of TCDS trigger types will not be booked for this run.";
259 
260  // book the rate histograms for the L1T triggers that are included in the L1T menu
261  booker.setCurrentFolder(m_dqm_path + "/L1T");
262  auto const &l1tMenu = setup.getData(m_l1tMenu_token);
263  for (auto const &keyval : l1tMenu.getAlgorithmMap()) {
264  unsigned int const bit = keyval.second.getIndex();
265  if (bit >= histograms.l1t_counts.size()) {
266  edm::LogError("TriggerRatesMonitor")
267  << "This should never happen: bit of L1T algorithm (bit=" << bit << ", name=\"" << keyval.first
268  << "\") is not smaller than size of \"histograms.l1t_counts\" vector (size=" << histograms.l1t_counts.size()
269  << ").\nRate histogram of this L1T algorithm will not be booked for this run.";
270  continue;
271  }
272  bool masked = false; // FIXME read L1T masks once they will be avaiable in the EventSetup
273  std::string const &name = fmt::sprintf("%s (bit %d)", keyval.first, bit);
274  std::string const &title =
275  fmt::sprintf("%s (bit %d)%s vs. lumisection", keyval.first, bit, (masked ? " (masked)" : ""));
276  histograms.l1t_counts[bit] = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
277  }
278 
279  if (histograms.hltConfig.inited()) {
280  auto const &datasetNames = histograms.hltConfig.datasetNames();
281 
282  // book the rate histograms for the HLT triggers
283  for (unsigned int d = 0; d < datasetNames.size(); ++d) {
284  booker.setCurrentFolder(m_dqm_path + "/HLT/" + datasetNames[d]);
285  for (unsigned int i = 0; i < histograms.datasets[d].size(); ++i) {
286  unsigned int index = histograms.datasets[d][i];
287  std::string const &name = histograms.hltConfig.triggerName(index);
288  histograms.hlt_by_dataset_counts[d][i].pass_l1_seed = booker.book1D(name + "_pass_L1_seed",
289  name + " pass L1 seed, vs. lumisection",
291  -0.5,
292  m_lumisections_range + 0.5);
293  histograms.hlt_by_dataset_counts[d][i].pass_prescale = booker.book1D(name + "_pass_prescaler",
294  name + " pass prescaler, vs. lumisection",
296  -0.5,
297  m_lumisections_range + 0.5);
298  histograms.hlt_by_dataset_counts[d][i].accept = booker.book1D(name + "_accept",
299  name + " accept, vs. lumisection",
301  -0.5,
302  m_lumisections_range + 0.5);
303  histograms.hlt_by_dataset_counts[d][i].reject = booker.book1D(name + "_reject",
304  name + " reject, vs. lumisection",
306  -0.5,
307  m_lumisections_range + 0.5);
308  histograms.hlt_by_dataset_counts[d][i].error = booker.book1D(name + "_error",
309  name + " error, vs. lumisection",
311  -0.5,
312  m_lumisections_range + 0.5);
313  }
314 
315  for (unsigned int i : histograms.datasets[d]) {
316  // look for the index of the (last) L1T seed and prescale module in each path
317  histograms.hltIndices[i].index_l1_seed = histograms.hltConfig.size(i);
318  histograms.hltIndices[i].index_prescale = histograms.hltConfig.size(i);
319  for (unsigned int j = 0; j < histograms.hltConfig.size(i); ++j) {
320  std::string const &label = histograms.hltConfig.moduleLabel(i, j);
321  std::string const &type = histograms.hltConfig.moduleType(label);
322  if (type == "HLTL1TSeed" or type == "HLTLevel1GTSeed" or type == "HLTLevel1Activity" or
323  type == "HLTLevel1Pattern") {
324  // there might be more L1T seed filters in sequence
325  // keep looking and store the index of the last one
326  histograms.hltIndices[i].index_l1_seed = j;
327  } else if (type == "HLTPrescaler") {
328  // there should be only one prescaler in a path, and it should follow all L1T seed filters
329  histograms.hltIndices[i].index_prescale = j;
330  break;
331  }
332  }
333  }
334  }
335 
336  // book the rate histograms for the HLT datasets
337  booker.setCurrentFolder(m_dqm_path + "/Datasets");
338  for (unsigned int i = 0; i < datasetNames.size(); ++i)
339  histograms.dataset_counts[i] =
341 
342  // book the rate histograms for the HLT streams
343  booker.setCurrentFolder(m_dqm_path + "/Streams");
344  auto const &streamNames = histograms.hltConfig.streamNames();
345  for (unsigned int i = 0; i < streamNames.size(); ++i)
346  histograms.stream_counts[i] =
347  booker.book1D(streamNames[i], streamNames[i], m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5);
348  }
349 }
350 
352  edm::EventSetup const &setup,
353  RunBasedHistograms const &histograms) const {
354  unsigned int lumisection = event.luminosityBlock();
355 
356  // monitor the overall event count and event types rates
357  histograms.events_processed->Fill(lumisection);
358  if (histograms.tcds_counts[event.experimentType()])
359  histograms.tcds_counts[event.experimentType()]->Fill(lumisection);
360 
361  // monitor the rates of L1T triggers
362  auto const &algBlkBxVecHandle = event.getHandle(m_l1t_results_token);
363  if (not algBlkBxVecHandle.isValid()) {
364  edm::LogError("TriggerRatesMonitor")
365  << "L1 trigger results with label [" << m_l1t_results_inputTag.encode()
366  << "] not present or invalid. MonitorElements of L1T results not filled for this event.";
367  } else if (algBlkBxVecHandle->isEmpty(0)) {
368  edm::LogError("TriggerRatesMonitor")
369  << "L1 trigger results with label [" << m_l1t_results_inputTag.encode()
370  << "] empty for BX=0. MonitorElements of L1T results not filled for this event.";
371  } else {
372  auto const &results = algBlkBxVecHandle->at(0, 0);
373  for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i)
374  if (results.getAlgoDecisionFinal(i))
375  if (histograms.l1t_counts[i])
376  histograms.l1t_counts[i]->Fill(lumisection);
377  }
378 
379  // monitor the rates of HLT triggers, datasets and streams
380  if (histograms.hltConfig.inited()) {
381  auto const &hltResults = event.get(m_hlt_results_token);
382  if (hltResults.size() != histograms.hltIndices.size()) {
383  edm::LogError("TriggerRatesMonitor")
384  << "This should never happen: the number of HLT paths has changed since the beginning of the run"
385  << " (from " << histograms.hltIndices.size() << " to " << hltResults.size() << ").\n"
386  << "Histograms for rates of HLT paths, datasets and streams will not be filled for this event.";
387  return;
388  }
389 
390  for (unsigned int d = 0; d < histograms.datasets.size(); ++d) {
391  for (unsigned int i : histograms.datasets[d])
392  if (hltResults[i].accept()) {
393  histograms.dataset_counts[d]->Fill(lumisection);
394  // ensure each dataset is incremented only once per event
395  break;
396  }
397  for (unsigned int i = 0; i < histograms.datasets[d].size(); ++i) {
398  unsigned int const index = histograms.datasets[d][i];
400 
401  if (path.index() > histograms.hltIndices[index].index_l1_seed)
402  histograms.hlt_by_dataset_counts[d][i].pass_l1_seed->Fill(lumisection);
403  if (path.index() > histograms.hltIndices[index].index_prescale)
404  histograms.hlt_by_dataset_counts[d][i].pass_prescale->Fill(lumisection);
405  if (path.accept())
406  histograms.hlt_by_dataset_counts[d][i].accept->Fill(lumisection);
407  else if (path.error())
408  histograms.hlt_by_dataset_counts[d][i].error->Fill(lumisection);
409  else
410  histograms.hlt_by_dataset_counts[d][i].reject->Fill(lumisection);
411  }
412  }
413 
414  for (unsigned int i = 0; i < histograms.streams.size(); ++i)
415  for (unsigned int j : histograms.streams[i])
416  if (hltResults[j].accept()) {
417  histograms.stream_counts[i]->Fill(lumisection);
418  // ensure each stream is incremented only once per event
419  break;
420  }
421  }
422 }
423 
424 // define this as a plugin
void dqmBeginRun(edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
const std::string m_dqm_path
dqm::impl::MonitorElement MonitorElement
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::string encode() const
Definition: InputTag.cc:159
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results_token
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: config.py:1
Log< level::Error, false > LogError
const uint32_t m_lumisections_range
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, RunBasedHistograms &) const override
static constexpr const char *const s_tcds_trigger_types[]
char const * label
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
def unique(seq, keepstr=True)
Definition: tier0.py:24
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results_token
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
d
Definition: ztail.py:151
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, RunBasedHistograms const &) const override
const edm::ESGetToken< L1TUtmTriggerMenu, L1TUtmTriggerMenuRcd > m_l1tMenu_token
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TriggerRatesMonitor(edm::ParameterSet const &)
results
Definition: mysort.py:8
static constexpr unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:52
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
~TriggerRatesMonitor() override=default
const edm::InputTag m_l1t_results_inputTag
Definition: event.py:1
Definition: Run.h:45
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
#define LogDebug(id)