CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
42 
43 // helper functions
44 template <typename T>
45 static
46 const T * get(const edm::Event & event, const edm::EDGetTokenT<T> & token) {
48  event.getByToken(token, handle);
49  if (not handle.isValid())
50  throw * handle.whyFailed();
51  return handle.product();
52 }
53 
54 template <typename R, typename T>
55 static
56 const T * get(const edm::EventSetup & setup) {
58  setup.get<R>().get(handle);
59  return handle.product();
60 }
61 
62 
64 public:
65  explicit TriggerRatesMonitor(edm::ParameterSet const &);
67 
68  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
69 
70 private:
71  virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
72  virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
73  virtual void analyze(edm::Event const &, edm::EventSetup const &) override;
74 
75  // module configuration
80 
81  // L1T and HLT configuration
82  struct HLTIndices {
83  unsigned int index_l1_seed;
84  unsigned int index_prescale;
85 
87  index_l1_seed( (unsigned int) -1),
88  index_prescale( (unsigned int) -1)
89  { }
90  };
91 
96  std::vector<HLTIndices> m_hltIndices;
97 
98  std::vector<std::vector<unsigned int>> m_datasets;
99  std::vector<std::vector<unsigned int>> m_streams;
100 
101 
102  struct HLTRatesPlots {
103  TH1F * wasrun;
104  TH1F * pass_l1_seed;
106  TH1F * accept;
107  TH1F * reject;
108  TH1F * error;
109  };
110 
111  // overall event count and event types
116 
117  // L1T triggers
118  std::vector<TH1F *> m_l1t_algo_counts;
119  std::vector<TH1F *> m_l1t_tech_counts;
120 
121  // HLT triggers
122  std::vector<HLTRatesPlots> m_hlt_counts;
123 
124  // datasets
125  std::vector<TH1F *> m_dataset_counts;
126 
127  // streams
128  std::vector<TH1F *> m_stream_counts;
129 
130 };
131 
132 
133 
135 {
137  desc.addUntracked<edm::InputTag>( "l1tResults", edm::InputTag("gtDigis"));
138  desc.addUntracked<edm::InputTag>( "hltResults", edm::InputTag("TriggerResults"));
139  desc.addUntracked<std::string>( "dqmPath", "HLT/TriggerRates" );
140  desc.addUntracked<uint32_t>( "lumisectionRange", 2500 ); // ~16 hours
141  descriptions.add("triggerRatesMonitor", desc);
142 }
143 
144 
146  // module configuration
147  m_l1t_results( consumes<L1GlobalTriggerReadoutRecord>( config.getUntrackedParameter<edm::InputTag>( "l1tResults" ) ) ),
148  m_hlt_results( consumes<edm::TriggerResults>( config.getUntrackedParameter<edm::InputTag>( "hltResults" ) ) ),
149  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) ),
150  m_lumisections_range( config.getUntrackedParameter<uint32_t>( "lumisectionRange" ) ),
151  // L1T and HLT configuration
152  m_l1tMenu( nullptr ),
153  m_l1tAlgoMask( nullptr ),
154  m_l1tTechMask( nullptr),
155  m_hltConfig(),
156  m_hltIndices(),
157  m_datasets(),
158  m_streams(),
159  // overall event count and event types
160  m_events_processed( nullptr ),
161  m_events_physics( nullptr ),
162  m_events_calibration( nullptr ),
163  m_events_random( nullptr ),
164  // L1T triggers
165  m_l1t_algo_counts(),
166  m_l1t_tech_counts(),
167  // HLT triggers
168  m_hlt_counts(),
169  // datasets
170  m_dataset_counts(),
171  // streams
172  m_stream_counts()
173 {
174 }
175 
177 {
178 }
179 
181 {
182  m_events_processed = nullptr;
183  m_events_physics = nullptr;
184  m_events_calibration = nullptr;
185  m_events_random = nullptr;
186 
187  // cache the L1 trigger menu
188  m_l1tMenu = get<L1GtTriggerMenuRcd, L1GtTriggerMenu>(setup);
189  m_l1tAlgoMask = get<L1GtTriggerMaskAlgoTrigRcd, L1GtTriggerMask>(setup);
190  m_l1tTechMask = get<L1GtTriggerMaskTechTrigRcd, L1GtTriggerMask>(setup);
191  if (m_l1tMenu and m_l1tAlgoMask and m_l1tTechMask) {
192  m_l1t_algo_counts.clear();
193  m_l1t_algo_counts.resize( m_l1tAlgoMask->gtTriggerMask().size(), nullptr );
194  m_l1t_tech_counts.clear();
195  m_l1t_tech_counts.resize( m_l1tTechMask->gtTriggerMask().size(), nullptr );
196  } else {
197  // L1GtUtils not initialised, skip the the L1T monitoring
198  edm::LogError("TriggerRatesMonitor") << "failed to read the L1 menu or masks from the EventSetup, the L1 trigger rates will not be monitored";
199  }
200 
201  // initialise the HLTConfigProvider
202  bool changed = true;
204  labelsForToken(m_hlt_results, labels);
205  if (m_hltConfig.init(run, setup, labels.process, changed)) {
206  m_hlt_counts.clear();
208  m_hltIndices.resize( m_hltConfig.size(), HLTIndices() );
209 
210  unsigned int datasets = m_hltConfig.datasetNames().size();
211  m_datasets.clear();
212  m_datasets.resize( datasets, {} );
213  for (unsigned int i = 0; i < datasets; ++i) {
214  auto const & paths = m_hltConfig.datasetContent(i);
215  m_datasets[i].reserve(paths.size());
216  for (auto const & path: paths)
218  }
219  m_dataset_counts.clear();
220  m_dataset_counts.resize( datasets, nullptr );
221 
222  unsigned int streams = m_hltConfig.streamNames().size();
223  m_streams.clear();
224  m_streams.resize( streams, {} );
225  for (unsigned int i = 0; i < streams; ++i) {
226  for (auto const & dataset : m_hltConfig.streamContent(i)) {
227  for (auto const & path : m_hltConfig.datasetContent(dataset))
228  m_streams[i].push_back(m_hltConfig.triggerIndex(path));
229  }
231  auto unique_end = std::unique(m_streams[i].begin(), m_streams[i].end());
232  m_streams[i].resize(unique_end - m_streams[i].begin());
233  m_streams[i].shrink_to_fit();
234  }
235  m_stream_counts.clear();
236  m_stream_counts.resize( streams, nullptr );
237  } else {
238  // HLTConfigProvider not initialised, skip the the HLT monitoring
239  edm::LogError("TriggerRatesMonitor") << "failed to initialise HLTConfigProvider, the HLT trigger and datasets rates will not be monitored";
240  }
241 }
242 
244 {
245  // book the overall event count and event types histograms
246  booker.setCurrentFolder( m_dqm_path );
247  m_events_processed = booker.book1D("processed", "Processed events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
248  m_events_physics = booker.book1D("physics", "Physics evenst", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
249  m_events_calibration = booker.book1D("calibration", "Calibration events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
250  m_events_random = booker.book1D("random", "Random events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
251 
252  if (m_l1tMenu and m_l1tAlgoMask) {
253  // book the rate histograms for the L1 Algorithm triggers
254  booker.setCurrentFolder( m_dqm_path + "/L1 Algo" );
255 
256  // book the histograms for L1 algo triggers that are included in the L1 menu
257  for (auto const & keyval: m_l1tMenu->gtAlgorithmAliasMap()) {
258  int bit = keyval.second.algoBitNumber();
259  // check if the trigger is unmasked in *any* partition
260  bool masked = ((m_l1tAlgoMask->gtTriggerMask().at(bit) & 0xff) == 0xff);
261  std::string const & name = (boost::format("%s (bit %d)") % keyval.first.substr(0, keyval.first.find_first_of(".")) % bit).str();
262  std::string const & title = (boost::format("%s (bit %d)%s") % keyval.first % bit % (masked ? " (masked)" : "")).str();
263  m_l1t_algo_counts.at(bit) = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
264  }
265  // book the histograms for L1 algo triggers that are not included in the L1 menu
266  for (unsigned int bit = 0; bit < m_l1tAlgoMask->gtTriggerMask().size(); ++bit) if (not m_l1t_algo_counts.at(bit)) {
267  // check if the trigger is unmasked in *any* partition
268  bool masked = ((m_l1tAlgoMask->gtTriggerMask().at(bit) & 0xff) == 0xff);
269  std::string const & name = (boost::format("L1 Algo (bit %d)") % bit).str();
270  std::string const & title = (boost::format("L1 Algo (bit %d)%s") % bit % (masked ? " (masked)" : "")).str();
271  m_l1t_algo_counts.at(bit) = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
272  }
273  }
274 
275  if (m_l1tMenu and m_l1tTechMask) {
276  // book the rate histograms for the L1 Technical triggers
277  booker.setCurrentFolder( m_dqm_path + "/L1 Tech" );
278 
279  // book the histograms for L1 tech triggers that are included in the L1 menu
280  for (auto const & keyval: m_l1tMenu->gtTechnicalTriggerMap()) {
281  int bit = keyval.second.algoBitNumber();
282  // check if the trigger is unmasked in *any* partition
283  bool masked = ((m_l1tTechMask->gtTriggerMask().at(bit) & 0xff) == 0xff);
284  std::string const & name = (boost::format("%s (bit %d)") % keyval.first.substr(0, keyval.first.find_first_of(".")) % bit).str();
285  std::string const & title = (boost::format("%s (bit %d)%s") % keyval.first % bit % (masked ? " (masked)" : "")).str();
286  m_l1t_tech_counts.at(bit) = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
287  }
288  // book the histograms for L1 tech triggers that are not included in the L1 menu
289  for (unsigned int bit = 0; bit < m_l1tTechMask->gtTriggerMask().size(); ++bit) if (not m_l1t_tech_counts.at(bit)) {
290  // check if the trigger is unmasked in *any* partition
291  bool masked = ((m_l1tTechMask->gtTriggerMask().at(bit) & 0xff) == 0xff);
292  std::string const & name = (boost::format("L1 Tech (bit %d)") % bit).str();
293  std::string const & title = (boost::format("L1 Tech (bit %d)%s") % bit % (masked ? " (masked)" : "")).str();
294  m_l1t_tech_counts.at(bit) = booker.book1D(name, title, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
295  }
296 
297  }
298 
299  if (m_hltConfig.inited()) {
300  // book the HLT triggers rate histograms
301  booker.setCurrentFolder( m_dqm_path + "/HLT" );
302  for (unsigned int i = 0; i < m_hltConfig.size(); ++i) {
304  m_hlt_counts[i].wasrun = booker.book1D(name + " counts", name + " counts", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
305  m_hlt_counts[i].pass_l1_seed = booker.book1D(name + " pass L1 seed", name + " pass L1 seed", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
306  m_hlt_counts[i].pass_prescale = booker.book1D(name + " pass prescaler", name + " pass prescaler", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
307  m_hlt_counts[i].accept = booker.book1D(name + " accept", name + " accept", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
308  m_hlt_counts[i].reject = booker.book1D(name + " reject", name + " reject", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
309  m_hlt_counts[i].error = booker.book1D(name + " error", name + " error", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
310  // look for the index of the (last) L1 seed and prescale module in each path
311  m_hltIndices[i].index_l1_seed = m_hltConfig.size(i);
312  m_hltIndices[i].index_prescale = m_hltConfig.size(i);
313  for (unsigned int j = 0; j < m_hltConfig.size(i); ++j) {
315  std::string const & type = m_hltConfig.moduleType(label);
316  if (type == "HLTLevel1GTSeed" or type == "HLTLevel1Activity" or type == "HLTLevel1Pattern") {
317  // there might be more L1 seed filters in sequence
318  // keep looking and store the index of the last one
319  m_hltIndices[i].index_l1_seed = j;
320  } else if (type == "HLTPrescaler") {
321  // there should be only one prescaler in a path, and it should follow all L1 seed filters
322  m_hltIndices[i].index_prescale = j;
323  break;
324  }
325 
326  }
327  }
328 
329  // book the HLT datasets rate histograms
330  booker.setCurrentFolder( m_dqm_path + "/Datasets" );
331  auto const & datasets = m_hltConfig.datasetNames();
332  for (unsigned int i = 0; i < datasets.size(); ++i)
334 
335  // book the HLT streams rate histograms
336  booker.setCurrentFolder( m_dqm_path + "/Streams" );
337  auto const & streams = m_hltConfig.streamNames();
338  for (unsigned int i = 0; i < streams.size(); ++i)
339  m_stream_counts[i] = booker.book1D(streams[i], streams[i], m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
340  }
341 }
342 
343 
345 {
346  unsigned int lumisection = event.luminosityBlock();
347 
348  // book the overall event count and event types rates
349  m_events_processed->Fill(lumisection);
350  switch (event.experimentType()) {
352  m_events_physics->Fill(lumisection);
353  break;
355  m_events_calibration->Fill(lumisection);
356  break;
358  m_events_random->Fill(lumisection);
359  break;
365  // ignore these event types
366  break;
367  }
368 
369  // monitor the L1 triggers rates
370  if (m_l1tMenu and m_l1tAlgoMask and m_l1tTechMask) {
371  L1GlobalTriggerReadoutRecord const & l1tResults = * get<L1GlobalTriggerReadoutRecord>(event, m_l1t_results);
372 
373  const std::vector<bool> & algoword = l1tResults.decisionWord();
374  if (algoword.size() == m_l1t_algo_counts.size()) {
375  for (unsigned int i = 0; i < m_l1t_algo_counts.size(); ++i)
376  if (algoword[i])
377  m_l1t_algo_counts[i]->Fill(lumisection);
378  } else {
379  edm::LogWarning("TriggerRatesMonitor") << "This should never happen: the size of the L1 Algo Trigger mask does not match the number of L1 Algo Triggers";
380  }
381 
382  const std::vector<bool> & techword = l1tResults.technicalTriggerWord();
383  if (techword.size() == m_l1t_tech_counts.size()) {
384  for (unsigned int i = 0; i < m_l1t_tech_counts.size(); ++i)
385  if (techword[i])
386  m_l1t_tech_counts[i]->Fill(lumisection);
387  } else {
388  edm::LogWarning("TriggerRatesMonitor") << "This should never happen: the size of the L1 Tech Trigger mask does not match the number of L1 Tech Triggers";
389  }
390  }
391 
392  // monitor the HLT triggers and datsets rates
393  if (m_hltConfig.inited()) {
394  edm::TriggerResults const & hltResults = * get<edm::TriggerResults>(event, m_hlt_results);
395  if (hltResults.size() == m_hlt_counts.size()) {
396  for (unsigned int i = 0; i < m_hlt_counts.size(); ++i) {
397  edm::HLTPathStatus const & path = hltResults.at(i);
398  if (path.wasrun())
399  m_hlt_counts[i].wasrun->Fill(lumisection);
400  if (path.index() > m_hltIndices[i].index_l1_seed)
401  m_hlt_counts[i].pass_l1_seed->Fill(lumisection);
402  if (path.index() > m_hltIndices[i].index_prescale)
403  m_hlt_counts[i].pass_prescale->Fill(lumisection);
404  if (path.accept())
405  m_hlt_counts[i].accept->Fill(lumisection);
406  else if (path.error())
407  m_hlt_counts[i].error ->Fill(lumisection);
408  else
409  m_hlt_counts[i].reject->Fill(lumisection);
410  }
411  } else {
412  edm::LogWarning("TriggerRatesMonitor") << "This should never happen: the number of HLT paths has changed since the beginning of the run";
413  }
414 
415  for (unsigned int i = 0; i < m_datasets.size(); ++i)
416  for (unsigned int j: m_datasets[i])
417  if (hltResults.at(j).accept()) {
418  m_dataset_counts[i]->Fill(lumisection);
419  // ensure each dataset is incremented only once per event
420  break;
421  }
422 
423  for (unsigned int i = 0; i < m_streams.size(); ++i)
424  for (unsigned int j: m_streams[i])
425  if (hltResults.at(j).accept()) {
426  m_stream_counts[i]->Fill(lumisection);
427  // ensure each stream is incremented only once per event
428  break;
429  }
430  }
431 }
432 
433 
434 //define this as a plug-in
unsigned int size() const
number of trigger paths in trigger table
type
Definition: HCALResponse.h:21
int i
Definition: DBlmapReader.cc:9
const std::string moduleType(const std::string &module) const
C++ class name of module.
const TechnicalTriggerWord & technicalTriggerWord(int bxInEventValue) const
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const std::string & triggerName(unsigned int triggerIndex) const
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
bool wasrun() const
was this path run?
Definition: HLTPathStatus.h:60
edm::EDGetTokenT< L1GlobalTriggerReadoutRecord > m_l1t_results
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< TH1F * > m_l1t_algo_counts
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
string format
Some error handling for the usage.
#define nullptr
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:25
L1GtTriggerMask const * m_l1tAlgoMask
const std::string & moduleLabel(unsigned int trigger, unsigned int module) const
bool inited() const
Accessors (const methods)
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)
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
std::vector< std::vector< unsigned int > > m_datasets
std::vector< TH1F * > m_l1t_tech_counts
virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
tuple handle
Definition: patZpeak.py:22
L1GtTriggerMask const * m_l1tTechMask
int j
Definition: DBlmapReader.cc:9
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
#define end
Definition: vmac.h:37
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:75
const DecisionWord & decisionWord(int bxInEventValue) const
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:274
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
tuple dataset
Definition: dataset.py:855
T const * product() const
Definition: ESHandle.h:86
TH1F * getTH1F(void) const
std::vector< HLTRatesPlots > m_hlt_counts
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
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
L1GtTriggerMenu const * m_l1tMenu
edm::EDGetTokenT< edm::TriggerResults > m_hlt_results
#define begin
Definition: vmac.h:30
TriggerRatesMonitor(edm::ParameterSet const &)
edm::EventAuxiliary::ExperimentType experimentType() const
Definition: EventBase.h:65
const AlgorithmMap & gtTechnicalTriggerMap() const
get / set the technical trigger map
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:110
const AlgorithmMap & gtAlgorithmAliasMap() const
get / set the algorithm map (by alias)
long double T
std::vector< HLTIndices > m_hltIndices
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
Definition: Run.h:41
const std::vector< std::string > & datasetNames() const
unsigned int index() const
Definition: HLTPathStatus.h:55