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