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 L1 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 
6 // C++ headers
7 #include <string>
8 #include <cstring>
9 
10 // boost headers
11 #include <boost/regex.hpp>
12 #include <boost/format.hpp>
13 
14 // Root headers
15 #include <TH1F.h>
16 
17 // CMSSW headers
40 
41 // helper functions
42 template <typename T>
43 static
44 const T & get(const edm::Event & event, const edm::EDGetTokenT<T> & token) {
46  event.getByToken(token, handle);
47  if (not handle.isValid())
48  throw * handle.whyFailed();
49  return * handle.product();
50 }
51 
52 template <typename R, typename T>
53 static
54 const T & get(const edm::EventSetup & setup) {
56  setup.get<R>().get(handle);
57  return * handle.product();
58 }
59 
60 
62 public:
63  explicit TriggerRatesMonitor(edm::ParameterSet const &);
65 
66  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
67 
68 private:
69  virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
70  virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
71  virtual void analyze(edm::Event const &, edm::EventSetup const &) override;
72 
73  // TCDS trigger types
74  // see https://twiki.cern.ch/twiki/bin/viewauth/CMS/TcdsEventRecord
75  static constexpr const char * const s_tcds_trigger_types[] = {
76  "Empty", // 0 - No trigger
77  "Physics", // 1 - GT trigger
78  "Calibration", // 2 - Sequence trigger (calibration)
79  "Random", // 3 - Random trigger
80  "Auxiliary", // 4 - Auxiliary (CPM front panel NIM input) trigger
81  nullptr, // 5 - reserved
82  nullptr, // 6 - reserved
83  nullptr, // 7 - reserved
84  "Cyclic", // 8 - Cyclic trigger
85  "Bunch-pattern", // 9 - Bunch-pattern trigger
86  "Software", // 10 - Software trigger
87  "TTS", // 11 - TTS-sourced trigger
88  nullptr, // 12 - reserved
89  nullptr, // 13 - reserved
90  nullptr, // 14 - reserved
91  nullptr // 15 - reserved
92  };
93 
94  // module configuration
99 
100  // L1T and HLT configuration
101 
103 
104  struct HLTIndices {
105  unsigned int index_l1_seed;
106  unsigned int index_prescale;
107 
109  index_l1_seed( (unsigned int) -1),
110  index_prescale( (unsigned int) -1)
111  { }
112  };
113 
115  std::vector<HLTIndices> m_hltIndices;
116 
117  std::vector<std::vector<unsigned int>> m_datasets;
118  std::vector<std::vector<unsigned int>> m_streams;
119 
120  // L1T and HLT rate plots
121 
122  struct HLTRatesPlots {
123  TH1F * pass_l1_seed;
125  TH1F * accept;
126  TH1F * reject;
127  TH1F * error;
128  };
129 
130  // overall event count and event types
132  std::vector<TH1F *> m_tcds_counts;
133 
134  // L1T triggers
135  std::vector<TH1F *> m_l1t_counts;
136 
137  // HLT triggers
138  std::vector<std::vector<HLTRatesPlots> > m_hlt_by_dataset_counts;
139 
140  // datasets
141  std::vector<TH1F *> m_dataset_counts;
142 
143  // streams
144  std::vector<TH1F *> m_stream_counts;
145 
146 };
147 
148 // definition
150 
151 
153 {
155  desc.addUntracked<edm::InputTag>( "l1tResults", edm::InputTag("gtStage2Digis"));
156  desc.addUntracked<edm::InputTag>( "hltResults", edm::InputTag("TriggerResults"));
157  desc.addUntracked<std::string>( "dqmPath", "HLT/TriggerRates" );
158  desc.addUntracked<uint32_t>( "lumisectionRange", 2500 ); // ~16 hours
159  descriptions.add("triggerRatesMonitor", desc);
160 }
161 
162 
164  // module configuration
165  m_l1t_results(consumes<GlobalAlgBlkBxCollection>( config.getUntrackedParameter<edm::InputTag>( "l1tResults" ) )),
166  m_hlt_results(consumes<edm::TriggerResults>( config.getUntrackedParameter<edm::InputTag>( "hltResults" ) )),
167  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) ),
168  m_lumisections_range( config.getUntrackedParameter<uint32_t>( "lumisectionRange" ) ),
169  // L1T and HLT configuration
171  m_hltConfig(),
172  m_hltIndices(),
173  m_datasets(),
174  m_streams(),
175  // overall event count and event types
177  m_tcds_counts(),
178  // L1T triggers
179  m_l1t_counts(),
180  // HLT triggers
182  // datasets
184  // streams
186 {
187 }
188 
190 {
191 }
192 
194 {
195  m_events_processed = nullptr;
196  m_tcds_counts.clear();
197  m_tcds_counts.resize(sizeof(s_tcds_trigger_types)/sizeof(const char *), nullptr);
198 
199  // cache the L1 trigger menu
200  m_l1tMenu = & get<L1TUtmTriggerMenuRcd, L1TUtmTriggerMenu>(setup);
201  if (m_l1tMenu) {
202  m_l1t_counts.clear();
204  } else {
205  edm::LogError("TriggerRatesMonitor") << "failed to read the L1 menu from the EventSetup, the L1 trigger rates will not be monitored";
206  }
207 
208  // initialise the HLTConfigProvider
209  bool changed = true;
211  labelsForToken(m_hlt_results, labels);
212  if (m_hltConfig.init(run, setup, labels.process, changed)) {
213  m_hltIndices.resize( m_hltConfig.size(), HLTIndices() );
214 
215  unsigned int datasets = m_hltConfig.datasetNames().size();
216  m_hlt_by_dataset_counts.clear();
217  m_hlt_by_dataset_counts.resize( datasets, {} );
218 
219  m_datasets.clear();
220  m_datasets.resize( datasets, {} );
221  for (unsigned int i = 0; i < datasets; ++i) {
222  auto const & paths = m_hltConfig.datasetContent(i);
223  m_hlt_by_dataset_counts[i].resize( paths.size(), HLTRatesPlots() );
224  m_datasets[i].reserve(paths.size());
225  for (auto const & path: paths) {
227  }
228  }
229  m_dataset_counts.clear();
230  m_dataset_counts.resize( datasets, nullptr );
231 
232  unsigned int streams = m_hltConfig.streamNames().size();
233  m_streams.clear();
234  m_streams.resize( streams, {} );
235  for (unsigned int i = 0; i < streams; ++i) {
236  for (auto const & dataset : m_hltConfig.streamContent(i)) {
237  for (auto const & path : m_hltConfig.datasetContent(dataset))
238  m_streams[i].push_back(m_hltConfig.triggerIndex(path));
239  }
240  std::sort(m_streams[i].begin(), m_streams[i].end());
241  auto unique_end = std::unique(m_streams[i].begin(), m_streams[i].end());
242  m_streams[i].resize(unique_end - m_streams[i].begin());
243  m_streams[i].shrink_to_fit();
244  }
245  m_stream_counts.clear();
246  m_stream_counts.resize( streams, nullptr );
247  } else {
248  // HLTConfigProvider not initialised, skip the the HLT monitoring
249  edm::LogError("TriggerRatesMonitor") << "failed to initialise HLTConfigProvider, the HLT trigger and datasets rates will not be monitored";
250  }
251 }
252 
254 {
255  // book the overall event count and event types histograms
256  booker.setCurrentFolder( m_dqm_path );
257  m_events_processed = booker.book1D("events", "Processed events vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
258  booker.setCurrentFolder( m_dqm_path + "/TCDS" );
259  for (unsigned int i = 0; i < sizeof(s_tcds_trigger_types)/sizeof(const char *); ++i)
260  if (s_tcds_trigger_types[i]) {
261  std::string const & title = (boost::format("%s events vs. lumisection") % s_tcds_trigger_types[i]).str();
263  }
264 
265  if (m_l1tMenu) {
266  // book the rate histograms for the L1 triggers that are included in the L1 menu
267  booker.setCurrentFolder( m_dqm_path + "/L1T" );
268  for (auto const & keyval: m_l1tMenu->getAlgorithmMap()) {
269  unsigned int bit = keyval.second.getIndex();
270  bool masked = false; // FIXME read L1 masks once they will be avaiable in the EventSetup
271  std::string const & name = (boost::format("%s (bit %d)") % keyval.first % bit).str();
272  std::string const & title = (boost::format("%s (bit %d)%s vs. lumisection") % keyval.first % bit % (masked ? " (masked)" : "")).str();
273  m_l1t_counts.at(bit) = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
274  }
275  }
276 
277  if (m_hltConfig.inited()) {
278 
279  auto const & datasets = m_hltConfig.datasetNames();
280 
281  // book the rate histograms for the HLT triggers
282  for (unsigned int d = 0; d < datasets.size(); ++d) {
283  booker.setCurrentFolder( m_dqm_path + "/HLT/" + datasets[d]);
284  for (unsigned int i = 0; i < m_datasets[d].size(); ++i) {
285  unsigned int index = m_datasets[d][i];
286  std::string const & name = m_hltConfig.triggerName(index);
287  m_hlt_by_dataset_counts[d][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)->getTH1F();
288  m_hlt_by_dataset_counts[d][i].pass_prescale = booker.book1D(name + "_pass_prescaler", name + " pass prescaler, vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
289  m_hlt_by_dataset_counts[d][i].accept = booker.book1D(name + "_accept", name + " accept, vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
290  m_hlt_by_dataset_counts[d][i].reject = booker.book1D(name + "_reject", name + " reject, vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
291  m_hlt_by_dataset_counts[d][i].error = booker.book1D(name + "_error", name + " error, vs. lumisection", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
292  }
293 
294  // booker.setCurrentFolder( m_dqm_path + "/HLT/" + datasets[d]);
295  for (unsigned int i: m_datasets[d]) {
296 
297  // look for the index of the (last) L1 seed and prescale module in each path
298  m_hltIndices[i].index_l1_seed = m_hltConfig.size(i);
299  m_hltIndices[i].index_prescale = m_hltConfig.size(i);
300  for (unsigned int j = 0; j < m_hltConfig.size(i); ++j) {
301  std::string const & label = m_hltConfig.moduleLabel(i, j);
302  std::string const & type = m_hltConfig.moduleType(label);
303  if (type == "HLTL1TSeed" or type == "HLTLevel1GTSeed" or type == "HLTLevel1Activity" or type == "HLTLevel1Pattern") {
304  // there might be more L1 seed filters in sequence
305  // keep looking and store the index of the last one
306  m_hltIndices[i].index_l1_seed = j;
307  } else if (type == "HLTPrescaler") {
308  // there should be only one prescaler in a path, and it should follow all L1 seed filters
309  m_hltIndices[i].index_prescale = j;
310  break;
311  }
312  }
313  }
314  }
315 
316  // book the HLT datasets rate histograms
317  booker.setCurrentFolder( m_dqm_path + "/Datasets" );
318  for (unsigned int i = 0; i < datasets.size(); ++i)
320 
321  // book the HLT streams rate histograms
322  booker.setCurrentFolder( m_dqm_path + "/Streams" );
323  auto const & streams = m_hltConfig.streamNames();
324  for (unsigned int i = 0; i < streams.size(); ++i)
325  m_stream_counts[i] = booker.book1D(streams[i], streams[i], m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
326  }
327 
328 }
329 
330 
332 {
333  unsigned int lumisection = event.luminosityBlock();
334 
335  // monitor the overall event count and event types rates
336  m_events_processed->Fill(lumisection);
337  if (m_tcds_counts[event.experimentType()])
338  m_tcds_counts[event.experimentType()]->Fill(lumisection);
339 
340  // monitor the L1 triggers rates
341  if (m_l1tMenu) {
342  auto const & bxvector = get<GlobalAlgBlkBxCollection>(event, m_l1t_results);
343  if (not bxvector.isEmpty(0)) {
344  auto const & results = bxvector.at(0, 0);
345  for (unsigned int i = 0; i < GlobalAlgBlk::maxPhysicsTriggers; ++i)
346  if (results.getAlgoDecisionFinal(i))
347  if (m_l1t_counts[i])
348  m_l1t_counts[i]->Fill(lumisection);
349  }
350  }
351 
352  // monitor the HLT triggers and datsets rates
353  if (m_hltConfig.inited()) {
354  edm::TriggerResults const & hltResults = get<edm::TriggerResults>(event, m_hlt_results);
355  if (hltResults.size() == m_hltIndices.size()) {
356  } else {
357  edm::LogWarning("TriggerRatesMonitor") << "This should never happen: the number of HLT paths has changed since the beginning of the run";
358  }
359 
360  for (unsigned int d = 0; d < m_datasets.size(); ++d) {
361  for (unsigned int i: m_datasets[d])
362  if (hltResults.at(i).accept()) {
363  m_dataset_counts[d]->Fill(lumisection);
364  // ensure each dataset is incremented only once per event
365  break;
366  }
367  for (unsigned int i = 0; i < m_datasets[d].size(); ++i) {
368  unsigned int index = m_datasets[d][i];
369  edm::HLTPathStatus const & path = hltResults.at(index);
370 
371  if (path.index() > m_hltIndices[index].index_l1_seed)
372  m_hlt_by_dataset_counts[d][i].pass_l1_seed->Fill(lumisection);
373  if (path.index() > m_hltIndices[index].index_prescale)
374  m_hlt_by_dataset_counts[d][i].pass_prescale->Fill(lumisection);
375  if (path.accept())
376  m_hlt_by_dataset_counts[d][i].accept->Fill(lumisection);
377  else if (path.error())
378  m_hlt_by_dataset_counts[d][i].error ->Fill(lumisection);
379  else
380  m_hlt_by_dataset_counts[d][i].reject->Fill(lumisection);
381  }
382  }
383 
384  for (unsigned int i = 0; i < m_streams.size(); ++i)
385  for (unsigned int j: m_streams[i])
386  if (hltResults.at(j).accept()) {
387  m_stream_counts[i]->Fill(lumisection);
388  // ensure each stream is incremented only once per event
389  break;
390  }
391  }
392 }
393 
394 
395 //define this as a plug-in
unsigned int size() const
number of trigger paths in trigger table
type
Definition: HCALResponse.h:21
const std::string moduleType(const std::string &module) const
C++ class name of module.
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const edm::EDGetTokenT< GlobalAlgBlkBxCollection > m_l1t_results
const std::string & triggerName(unsigned int triggerIndex) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
std::vector< std::vector< HLTRatesPlots > > m_hlt_by_dataset_counts
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: config.py:1
#define nullptr
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
#define constexpr
static const unsigned int maxPhysicsTriggers
Definition: GlobalAlgBlk.h:52
const std::string & moduleLabel(unsigned int trigger, unsigned int module) const
bool inited() const
Accessors (const methods)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
char const * process
Definition: ProductLabels.h:7
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
static const char *const s_tcds_trigger_types[]
std::vector< std::vector< unsigned int > > m_streams
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
const std::vector< std::string > & streamNames() const
std::vector< TH1F * > m_dataset_counts
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
unsigned int size() const
Get number of paths stored.
const std::vector< std::string > & streamContent(unsigned int stream) const
names of datasets in stream with index i
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
std::vector< std::vector< unsigned int > > m_datasets
virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
#define end
Definition: vmac.h:37
L1TUtmTriggerMenu const * m_l1tMenu
format
Some error handling for the usage.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:74
const HLTPathStatus & at(const unsigned int i) const
HLTConfigProvider m_hltConfig
const std::vector< std::string > & datasetContent(unsigned int dataset) const
names of trigger paths in dataset with index i
bool error() const
has this path encountered an error (exception)?
Definition: HLTPathStatus.h:64
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
T const * product() const
Definition: Handle.h:81
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
TH1F * getTH1F(void) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TH1F * > m_stream_counts
bool accept() const
has this path accepted the event?
Definition: HLTPathStatus.h:62
const std::map< std::string, L1TUtmAlgorithm > & getAlgorithmMap() const
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
#define begin
Definition: vmac.h:30
HLT enums.
TriggerRatesMonitor(edm::ParameterSet const &)
edm::EventAuxiliary::ExperimentType experimentType() const
Definition: EventBase.h:63
std::vector< TH1F * > m_l1t_counts
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:106
long double T
std::vector< HLTIndices > m_hltIndices
std::vector< TH1F * > m_tcds_counts
T const * product() const
Definition: ESHandle.h:86
Definition: event.py:1
Definition: Run.h:42
const std::vector< std::string > & datasetNames() const
unsigned int index() const
Definition: HLTPathStatus.h:55