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 
13 // Root headers
14 #include <TH1F.h>
15 
16 // CMSSW headers
41 
42 // length of a lumisections, corresponding to 2**18 LHC orbits, or 23.31 seconds
43 static const double SECS_PER_LUMI = 23.31040958083832;
44 
45 // helper functions
46 template <typename T>
47 static
48 const T * get(const edm::Event & event, const edm::EDGetTokenT<T> & token) {
50  event.getByToken(token, handle);
51  if (not handle.isValid())
52  throw * handle.whyFailed();
53  return handle.product();
54 }
55 
56 template <typename R, typename T>
57 static
58 const T * get(const edm::EventSetup & setup) {
60  setup.get<R>().get(handle);
61  return handle.product();
62 }
63 
64 
66 public:
67  explicit TriggerRatesMonitor(edm::ParameterSet const &);
69 
70  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
71 
72 private:
73  virtual void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override;
74  virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
75  virtual void analyze(edm::Event const &, edm::EventSetup const &) override;
76 
77  // module configuration
82 
83  // L1T and HLT configuration
84  struct HLTIndices {
85  unsigned int index_l1_seed;
86  unsigned int index_prescale;
87 
89  index_l1_seed( (unsigned int) -1),
90  index_prescale( (unsigned int) -1)
91  { }
92  };
93 
98  std::vector<HLTIndices> m_hltIndices;
99 
100  std::vector<std::vector<unsigned int>> m_datasets;
101 
102 
103  struct HLTRatesPlots {
104  TH1F * wasrun;
105  TH1F * pass_l1_seed;
107  TH1F * accept;
108  TH1F * reject;
109  TH1F * error;
110  };
111 
112  // overall event count and event types
117 
118  // L1T triggers
119  std::vector<TH1F *> m_l1t_algo_counts;
120  std::vector<TH1F *> m_l1t_tech_counts;
121 
122  // HLT triggers
123  std::vector<HLTRatesPlots> m_hlt_counts;
124 
125  // datasets
126  std::vector<TH1F *> m_dataset_counts;
127 
128 };
129 
130 
131 
133 {
135  desc.addUntracked<edm::InputTag>( "l1tResults", edm::InputTag("gtDigis"));
136  desc.addUntracked<edm::InputTag>( "hltResults", edm::InputTag("TriggerResults"));
137  desc.addUntracked<std::string>( "dqmPath", "HLT/TriggerRates" );
138  desc.addUntracked<uint32_t>( "lumisectionRange", 2500 ); // ~16 hours
139  descriptions.add("triggerRatesMonitor", desc);
140 }
141 
142 
144  // module configuration
145  m_l1t_results( consumes<L1GlobalTriggerReadoutRecord>( config.getUntrackedParameter<edm::InputTag>( "l1tResults" ) ) ),
146  m_hlt_results( consumes<edm::TriggerResults>( config.getUntrackedParameter<edm::InputTag>( "hltResults" ) ) ),
147  m_dqm_path( config.getUntrackedParameter<std::string>( "dqmPath" ) ),
148  m_lumisections_range( config.getUntrackedParameter<uint32_t>( "lumisectionRange" ) ),
149  // L1T and HLT configuration
150  m_l1tMenu( nullptr ),
151  m_l1tAlgoMask( nullptr ),
152  m_l1tTechMask( nullptr),
153  m_hltConfig(),
154  m_hltIndices(),
155  m_datasets(),
156  // overall event count and event types
157  m_events_processed( nullptr ),
158  m_events_physics( nullptr ),
159  m_events_calibration( nullptr ),
160  m_events_random( nullptr ),
161  // L1T triggers
162  m_l1t_algo_counts(),
163  m_l1t_tech_counts(),
164  // HLT triggers
165  m_hlt_counts(),
166  // datasets
167  m_dataset_counts()
168 {
169 }
170 
172 {
173 }
174 
176 {
177  m_events_processed = nullptr;
178  m_events_physics = nullptr;
179  m_events_calibration = nullptr;
180  m_events_random = nullptr;
181 
182  // cache the L1 trigger menu
183  m_l1tMenu = get<L1GtTriggerMenuRcd, L1GtTriggerMenu>(setup);
184  // FIXME - do we really need this ?
185  //(const_cast<L1GtTriggerMenu *>(m_l1tMenu))->buildGtConditionMap();
186  m_l1tAlgoMask = get<L1GtTriggerMaskAlgoTrigRcd, L1GtTriggerMask>(setup);
187  m_l1tTechMask = get<L1GtTriggerMaskTechTrigRcd, L1GtTriggerMask>(setup);
188  if (m_l1tMenu and m_l1tAlgoMask and m_l1tTechMask) {
189  m_l1t_algo_counts.clear();
190  m_l1t_algo_counts.resize( m_l1tAlgoMask->gtTriggerMask().size(), nullptr );
191  m_l1t_tech_counts.clear();
192  m_l1t_tech_counts.resize( m_l1tTechMask->gtTriggerMask().size(), nullptr );
193  } else {
194  // L1GtUtils not initialised, skip the the L1T monitoring
195  edm::LogError("TriggerRatesMonitor") << "failed to read the L1 menu or masks from the EventSetup, the L1 trigger rates will not be monitored";
196  }
197 
198  // initialise the HLTConfigProvider
199  bool changed = true;
201  labelsForToken(m_hlt_results, labels);
202  if (m_hltConfig.init(run, setup, labels.process, changed)) {
203  m_hlt_counts.clear();
205  m_hltIndices.resize( m_hltConfig.size(), HLTIndices() );
206 
207  unsigned int datasets = m_hltConfig.datasetNames().size();
208  m_datasets.clear();
209  m_datasets.resize( datasets, {} );
210  for (unsigned int i = 0; i < datasets; ++i) {
211  auto const & paths = m_hltConfig.datasetContent(i);
212  m_datasets[i].reserve(paths.size());
213  for (auto const & path: paths)
215  }
216  m_dataset_counts.clear();
217  m_dataset_counts.resize( datasets, nullptr );
218  } else {
219  // HLTConfigProvider not initialised, skip the the HLT monitoring
220  edm::LogError("TriggerRatesMonitor") << "failed to initialise HLTConfigProvider, the HLT trigger and datasets rates will not be monitored";
221  }
222 }
223 
225 {
226  // book the overall event count and event types histograms
227  booker.setCurrentFolder( m_dqm_path );
228  m_events_processed = booker.book1D("processed", "Processed events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
229  m_events_physics = booker.book1D("physics", "Physics evenst", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
230  m_events_calibration = booker.book1D("calibration", "Calibration events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
231  m_events_random = booker.book1D("random", "Random events", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
232 
233  if (m_l1tMenu and m_l1tAlgoMask and m_l1tTechMask) {
234  // book the L1T rate histograms
235  booker.setCurrentFolder( m_dqm_path + "/L1" );
236 
237  for (auto const & keyval: m_l1tMenu->gtAlgorithmAliasMap()) {
238  int bit = keyval.second.algoBitNumber();
239  std::string name = keyval.first.substr( 0, keyval.first.find_first_of(".") );
240  // check if the trigger is unmasked in *any* partition
241  if ((m_l1tAlgoMask->gtTriggerMask()[bit] & 0xff) != 0xff)
242  m_l1t_algo_counts[bit] = booker.book1D(name, keyval.first, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
243  else
244  m_l1t_algo_counts[bit] = booker.book1D(name, keyval.first + " (masked)", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
245  }
246 
247  for (auto const & keyval: m_l1tMenu->gtTechnicalTriggerMap()) {
248  int bit = keyval.second.algoBitNumber();
249  std::string name = keyval.first.substr( 0, keyval.first.find_first_of(".") );
250  // check if the trigger is unmasked in *any* partition
251  if ((m_l1tTechMask->gtTriggerMask()[bit] & 0xff) != 0xff)
252  m_l1t_tech_counts[bit] = booker.book1D(name, keyval.first, m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
253  else
254  m_l1t_tech_counts[bit] = booker.book1D(name, keyval.first + " (masked)", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
255  }
256  }
257 
258  if (m_hltConfig.inited()) {
259  // book the HLT triggers rate histograms
260  booker.setCurrentFolder( m_dqm_path + "/HLT" );
261  for (unsigned int i = 0; i < m_hltConfig.size(); ++i) {
263  m_hlt_counts[i].wasrun = booker.book1D(name + "_wasrun", name + " counts", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
264  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();
265  m_hlt_counts[i].pass_prescale = booker.book1D(name + "_pass_prescale", name + " pass prescaler", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
266  m_hlt_counts[i].accept = booker.book1D(name + "_accept", name + " trigger", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
267  m_hlt_counts[i].reject = booker.book1D(name + "_reject", name + " reject", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
268  m_hlt_counts[i].error = booker.book1D(name + "_error", name + " error count", m_lumisections_range + 1, -0.5, m_lumisections_range + 0.5)->getTH1F();
269  // look for the index of the (last) L1 seed and prescale module in each path
270  m_hltIndices[i].index_l1_seed = m_hltConfig.size(i);
271  m_hltIndices[i].index_prescale = m_hltConfig.size(i);
272  for (unsigned int j = 0; j < m_hltConfig.size(i); ++j) {
274  std::string const & type = m_hltConfig.moduleType(label);
275  if (type == "HLTLevel1GTSeed" or type == "HLTLevel1Activity" or type == "HLTLevel1Pattern") {
276  // there might be more L1 seed filters in sequence
277  // keep looking and store the index of the last one
278  m_hltIndices[i].index_l1_seed = j;
279  } else if (type == "HLTPrescaler") {
280  // there should be only one prescaler in a path, and it should follow all L1 seed filters
281  m_hltIndices[i].index_prescale = j;
282  break;
283  }
284 
285  }
286  }
287 
288  // book the HLT datasets rate histograms
289  booker.setCurrentFolder( m_dqm_path + "/Datasets" );
290  auto const & datasets = m_hltConfig.datasetNames();
291  for (unsigned int i = 0; i < datasets.size(); ++i)
293  }
294 }
295 
296 
298 {
299  unsigned int lumisection = event.luminosityBlock();
300 
301  // book the overall event count and event types rates
302  m_events_processed->Fill(lumisection);
303  switch (event.experimentType()) {
305  m_events_physics->Fill(lumisection);
306  break;
308  m_events_calibration->Fill(lumisection);
309  break;
311  m_events_random->Fill(lumisection);
312  break;
318  // ignore these event types
319  break;
320  }
321 
322  // monitor the L1 triggers rates
323  if (m_l1tMenu and m_l1tAlgoMask and m_l1tTechMask) {
324  L1GlobalTriggerReadoutRecord const & l1tResults = * get<L1GlobalTriggerReadoutRecord>(event, m_l1t_results);
325 
326  const std::vector<bool> & algoword = l1tResults.decisionWord();
327  assert(algoword.size() == m_l1t_algo_counts.size());
328  for (unsigned int i = 0; i < m_l1t_algo_counts.size(); ++i)
329  if (algoword[i])
330  m_l1t_algo_counts[i]->Fill(lumisection);
331 
332  const std::vector<bool> & techword = l1tResults.technicalTriggerWord();
333  assert(techword.size() == m_l1t_tech_counts.size());
334  for (unsigned int i = 0; i < m_l1t_tech_counts.size(); ++i)
335  if (techword[i])
336  m_l1t_tech_counts[i]->Fill(lumisection);
337  }
338 
339  // monitor the HLT triggers and datsets rates
340  if (m_hltConfig.inited()) {
341  edm::TriggerResults const & hltResults = * get<edm::TriggerResults>(event, m_hlt_results);
342  assert(hltResults.size() == m_hlt_counts.size());
343  for (unsigned int i = 0; i < m_hlt_counts.size(); ++i) {
344  edm::HLTPathStatus const & path = hltResults.at(i);
345  if (path.wasrun())
346  m_hlt_counts[i].wasrun->Fill(lumisection);
347  if (path.index() > m_hltIndices[i].index_l1_seed)
348  m_hlt_counts[i].pass_l1_seed->Fill(lumisection);
349  if (path.index() > m_hltIndices[i].index_prescale)
350  m_hlt_counts[i].pass_prescale->Fill(lumisection);
351  if (path.accept())
352  m_hlt_counts[i].accept->Fill(lumisection);
353  else if (path.error())
354  m_hlt_counts[i].error ->Fill(lumisection);
355  else
356  m_hlt_counts[i].reject->Fill(lumisection);
357  }
358 
359  for (unsigned int i = 0; i < m_datasets.size(); ++i)
360  for (unsigned int j: m_datasets[i])
361  if (hltResults.at(j).accept()) {
362  m_dataset_counts[i]->Fill(lumisection);
363  // ensure each dataset is incremented only once per event
364  break;
365  }
366  }
367 }
368 
369 
370 //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
#define nullptr
std::vector< TH1F * > m_l1t_algo_counts
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:26
L1GtTriggerMask const * m_l1tAlgoMask
const std::string & moduleLabel(unsigned int trigger, unsigned int module) const
bool inited() const
Accessors (const methods)
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
tuple TriggerResults
Definition: old-fu_pass.py:28
tuple path
else: Piece not in the list, fine.
std::vector< TH1F * > m_dataset_counts
const std::vector< unsigned int > & gtTriggerMask() const
get the trigger mask
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
unsigned int size() const
Get number of paths stored.
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:113
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
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:76
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
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
T const * product() const
Definition: ESHandle.h:62
TH1F * getTH1F(void) const
std::vector< HLTRatesPlots > m_hlt_counts
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * product() const
Definition: Handle.h:81
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
TriggerRatesMonitor(edm::ParameterSet const &)
edm::EventAuxiliary::ExperimentType experimentType() const
Definition: EventBase.h:61
static const double SECS_PER_LUMI
const AlgorithmMap & gtTechnicalTriggerMap() const
get / set the technical trigger map
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:111
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