CMS 3D CMS Logo

HLTrigReport.cc
Go to the documentation of this file.
1 
10 #include <iomanip>
11 #include <cstring>
12 #include <sstream>
13 
19 #include "HLTrigReportService.h"
20 #include "Math/QuantFuncMathCore.h"
21 
22 #include "HLTrigReport.h"
23 
25  if (value == "never")
26  return NEVER;
27 
28  if (value == "job")
29  return EVERY_JOB;
30 
31  if (value == "run")
32  return EVERY_RUN;
33 
34  if (value == "lumi")
35  return EVERY_LUMI;
36 
37  if (value == "event")
38  return EVERY_EVENT;
39 
40  throw cms::Exception("Configuration") << "Invalid option value \"" << value
41  << "\". Legal values are \"job\", \"run\", \"lumi\", \"event\" and \"never\".";
42 }
43 
44 //
45 // constructors and destructor
46 //
48  : nEvents_(0),
49  nWasRun_(0),
50  nAccept_(0),
51  nErrors_(0),
52  hlWasRun_(0),
53  hltL1s_(0),
54  hltPre_(0),
55  hlAccept_(0),
56  hlAccTot_(0),
57  hlErrors_(0),
58  hlAccTotDS_(0),
59  dsAccTotS_(0) {}
60 
61 hltrigreport::Accumulate::Accumulate(size_t numHLNames,
62  std::vector<std::vector<unsigned int> > const& hlIndex,
63  std::vector<std::vector<unsigned int> > const& dsIndex)
64  : nEvents_(0),
65  nWasRun_(0),
66  nAccept_(0),
67  nErrors_(0),
68  hlWasRun_(numHLNames),
69  hltL1s_(numHLNames),
70  hltPre_(numHLNames),
71  hlAccept_(numHLNames),
72  hlAccTot_(numHLNames),
73  hlErrors_(numHLNames),
74  hlAccTotDS_(hlIndex.size()),
75  hlAllTotDS_(hlIndex.size()),
76  dsAccTotS_(dsIndex.size()),
77  dsAllTotS_(dsIndex.size()) {
78  for (size_t ds = 0; ds < hlIndex.size(); ++ds) {
79  hlAccTotDS_[ds].resize(hlIndex[ds].size());
80  }
81 
82  for (size_t s = 0; s < dsIndex.size(); ++s) {
83  dsAccTotS_[s].resize(dsIndex[s].size());
84  }
85 }
86 
88  nEvents_ += iOther.nEvents_;
89  nWasRun_ += iOther.nWasRun_;
90  nAccept_ += iOther.nAccept_;
91  nErrors_ += iOther.nErrors_;
92 
93  auto vsum = [](auto& to, auto const& from) {
94  for (size_t i = 0; i < from.size(); ++i) {
95  to[i] += from[i];
96  }
97  };
98 
99  assert(hlWasRun_.size() == iOther.hlWasRun_.size());
100  vsum(hlWasRun_, iOther.hlWasRun_);
101  vsum(hltL1s_, iOther.hltL1s_);
102  vsum(hltPre_, iOther.hltPre_);
103  vsum(hlAccept_, iOther.hlAccept_);
104  vsum(hlAccTot_, iOther.hlAccTot_);
105  vsum(hlErrors_, iOther.hlErrors_);
106 
107  assert(hlAllTotDS_.size() == iOther.hlAllTotDS_.size());
108  vsum(hlAllTotDS_, iOther.hlAllTotDS_);
109  vsum(dsAllTotS_, iOther.dsAllTotS_);
110 
111  auto vvsum = [](auto& to, auto const& from) {
112  for (size_t i = 0; i < from.size(); ++i) {
113  assert(from[i].size() == to[i].size());
114  for (size_t j = 0; j < from[i].size(); ++j) {
115  to[i][j] += from[i][j];
116  }
117  }
118  };
119 
120  vvsum(hlAccTotDS_, iOther.hlAccTotDS_);
121  vvsum(dsAccTotS_, iOther.dsAccTotS_);
122 }
123 
125  nEvents_ = 0;
126  nWasRun_ = 0;
127  nAccept_ = 0;
128  nErrors_ = 0;
129 
130  auto vreset = [](auto& to) { std::fill(to.begin(), to.end(), 0); };
131 
132  vreset(hlWasRun_);
133  vreset(hltL1s_);
134  vreset(hltPre_);
135  vreset(hlAccept_);
136  vreset(hlAccTot_);
137  vreset(hlErrors_);
138 
139  vreset(hlAllTotDS_);
140  vreset(dsAllTotS_);
141 
142  auto vvreset = [&vreset](auto& to) {
143  for (auto& e : to) {
144  vreset(e);
145  }
146  };
147 
148  vvreset(hlAccTotDS_);
149  vvreset(dsAccTotS_);
150 }
151 
153  : hlTriggerResults_(iConfig.getParameter<edm::InputTag>("HLTriggerResults")),
154  hlTriggerResultsToken_(consumes<edm::TriggerResults>(hlTriggerResults_)),
155  configured_(false),
156  hlNames_(0),
157  hlIndex_(0),
158  posL1s_(0),
159  posPre_(0),
160  datasetNames_(0),
161  datasetContents_(0),
162  isCustomDatasets_(false),
163  dsIndex_(0),
164  streamNames_(0),
165  streamContents_(0),
166  isCustomStreams_(false),
167  refPath_("HLTriggerFinalPath"),
168  refIndex_(0),
169  refRate_(iConfig.getUntrackedParameter<double>("ReferenceRate", 100.0)),
170  reportBy_(decode(iConfig.getUntrackedParameter<std::string>("reportBy", "job"))),
171  resetBy_(decode(iConfig.getUntrackedParameter<std::string>("resetBy", "never"))),
172  serviceBy_(decode(iConfig.getUntrackedParameter<std::string>("serviceBy", "never"))),
173  hltConfig_() {
174  const edm::ParameterSet customDatasets(
175  iConfig.getUntrackedParameter<edm::ParameterSet>("CustomDatasets", edm::ParameterSet()));
176  isCustomDatasets_ = (customDatasets != edm::ParameterSet());
177  if (isCustomDatasets_) {
178  datasetNames_ = customDatasets.getParameterNamesForType<std::vector<std::string> >();
179  for (std::vector<std::string>::const_iterator name = datasetNames_.begin(); name != datasetNames_.end(); name++) {
180  datasetContents_.push_back(customDatasets.getParameter<std::vector<std::string> >(*name));
181  }
182  }
183 
184  const edm::ParameterSet customStreams(
185  iConfig.getUntrackedParameter<edm::ParameterSet>("CustomStreams", edm::ParameterSet()));
186  isCustomStreams_ = (customStreams != edm::ParameterSet());
187  if (isCustomStreams_) {
188  streamNames_ = customStreams.getParameterNamesForType<std::vector<std::string> >();
189  for (std::vector<std::string>::const_iterator name = streamNames_.begin(); name != streamNames_.end(); name++) {
190  streamContents_.push_back(customStreams.getParameter<std::vector<std::string> >(*name));
191  }
192  }
193 
194  refPath_ = iConfig.getUntrackedParameter<std::string>("ReferencePath", "HLTriggerFinalPath");
195  refIndex_ = 0;
196 
197  LogDebug("HLTrigReport") << "HL TiggerResults: " + hlTriggerResults_.encode()
198  << " using reference path and rate: " + refPath_ + " " << refRate_;
199 
201  edm::Service<HLTrigReportService>()->registerModule(this);
202  }
203 }
204 
205 HLTrigReport::~HLTrigReport() = default;
206 
209  desc.add<edm::InputTag>("HLTriggerResults", edm::InputTag("TriggerResults", "", "HLT"));
210  desc.addUntracked<std::string>("reportBy", "job");
211  desc.addUntracked<std::string>("resetBy", "never");
212  desc.addUntracked<std::string>("serviceBy", "never");
213 
214  edm::ParameterSetDescription customDatasetsParameters;
215  desc.addUntracked<edm::ParameterSetDescription>("CustomDatasets", customDatasetsParameters);
216  edm::ParameterSetDescription customStreamsParameters;
217  desc.addUntracked<edm::ParameterSetDescription>("CustomStreams", customStreamsParameters);
218  desc.addUntracked<std::string>("ReferencePath", "HLTriggerFinalPath");
219  desc.addUntracked<double>("ReferenceRate", 100.0);
220 
221  descriptions.add("hltTrigReport", desc);
222 }
223 
224 //
225 // member functions
226 //
227 
228 const std::vector<std::string>& HLTrigReport::datasetNames() const { return datasetNames_; }
229 const std::vector<std::string>& HLTrigReport::streamNames() const { return streamNames_; }
230 
232  // update trigger names
234 
235  const unsigned int n = hlNames_.size();
236 
237  // find the positions of seeding and prescaler modules
238  posL1s_.resize(n);
239  posPre_.resize(n);
240  for (unsigned int i = 0; i < n; ++i) {
241  posL1s_[i] = -1;
242  posPre_[i] = -1;
243  const std::vector<std::string>& moduleLabels(hltConfig_.moduleLabels(i));
244  for (unsigned int j = 0; j < moduleLabels.size(); ++j) {
245  const std::string& label = hltConfig_.moduleType(moduleLabels[j]);
246  if (label == "HLTLevel1GTSeed")
247  posL1s_[i] = j;
248  else if (label == "HLTPrescaler")
249  posPre_[i] = j;
250  }
251  }
252 
253  // if not overridden, reload the datasets and streams
254  if (not isCustomDatasets_) {
257  }
258  if (not isCustomStreams_) {
261  }
262 
263  // fill the matrices of hlIndex_, hlAccTotDS_
264  hlIndex_.clear();
265  hlIndex_.resize(datasetNames_.size());
266  for (unsigned int ds = 0; ds < datasetNames_.size(); ds++) {
267  unsigned int size = datasetContents_[ds].size();
268  hlIndex_[ds].reserve(size);
269  for (unsigned int p = 0; p < size; ++p) {
270  unsigned int i = hltConfig_.triggerIndex(datasetContents_[ds][p]);
271  if (i < n) {
272  hlIndex_[ds].push_back(i);
273  }
274  }
275  }
276 
277  // fill the matrices of dsIndex_, dsAccTotS_
278  dsIndex_.clear();
279  dsIndex_.resize(streamNames_.size());
280  for (unsigned int s = 0; s < streamNames_.size(); ++s) {
281  unsigned int size = streamContents_[s].size();
282  dsIndex_.reserve(size);
283  for (unsigned int ds = 0; ds < size; ++ds) {
284  unsigned int i = 0;
285  for (; i < datasetNames_.size(); i++)
286  if (datasetNames_[i] == streamContents_[s][ds])
287  break;
288  // report only datasets that have at least one path otherwise crash
289  if (i < datasetNames_.size() and !hlIndex_[i].empty()) {
290  dsIndex_[s].push_back(i);
291  }
292  }
293  }
294 
295  // if needed, update the reference path
297  if (refIndex_ >= n) {
298  refIndex_ = 0;
299  edm::LogWarning("HLTrigReport") << "Requested reference path '" + refPath_ + "' not in HLT menu. "
300  << "Using HLTriggerFinalPath instead.";
301  refPath_ = "HLTriggerFinalPath";
303  if (refIndex_ >= n) {
304  refIndex_ = 0;
305  edm::LogWarning("HLTrigReport") << "Requested reference path '" + refPath_ + "' not in HLT menu. "
306  << "Using first path in table (index=0) instead.";
307  }
308  }
309 
313  }
314 }
315 
317 
319  if (resetBy_ == EVERY_JOB)
320  reset();
321 }
322 
324  if (reportBy_ == EVERY_JOB)
325  dumpReport(accumulate_, "Summary for Job");
326  if (serviceBy_ == EVERY_JOB) {
328  }
329 }
330 
331 void HLTrigReport::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
332  bool changed = true;
333  if (hltConfig_.init(iRun, iSetup, hlTriggerResults_.process(), changed)) {
334  configured_ = true;
335  if (changed) {
336  dumpReport(accumulate_, "Summary for this HLT table");
339  }
340  } else {
341  dumpReport(accumulate_, "Summary for this HLT table");
342  // cannot initialize the HLT menu - reset and clear all counters and tables
343  configured_ = false;
344 
346  }
347 
348  if (resetBy_ == EVERY_RUN)
349  reset();
350 }
351 
353  if (reportBy_ == EVERY_RUN) {
354  std::stringstream stream;
355  stream << "Summary for Run " << run.run();
356  dumpReport(accumulate_, stream.str());
357  }
358  if (serviceBy_ == EVERY_RUN) {
360  }
361 }
362 
363 std::shared_ptr<HLTrigReport::Accumulate> HLTrigReport::globalBeginLuminosityBlock(edm::LuminosityBlock const& lumi,
364  edm::EventSetup const& setup) const {
365  if (useLumiCache()) {
366  if (not configured_) {
367  return std::make_shared<Accumulate>();
368  }
369  return std::make_shared<Accumulate>(hlNames_.size(), hlIndex_, dsIndex_);
370  }
371  return std::shared_ptr<Accumulate>();
372 }
373 
375  if (not useLumiCache()) {
376  return;
377  }
378 
379  if (resetBy_ == EVERY_LUMI and readAfterLumi()) {
380  //we will be reporting the last processed lumi
381  accumulate_ = std::move(*luminosityBlockCache(lumi.index()));
383  //we need to add this lumi's info to the longer lived accumulation
384  accumulate_.accumulate(*luminosityBlockCache(lumi.index()));
385  }
386 
387  if (reportBy_ == EVERY_LUMI) {
388  std::stringstream stream;
389  stream << "Summary for Run " << lumi.run() << ", LumiSection " << lumi.luminosityBlock();
391  dumpReport(*luminosityBlockCache(lumi.index()), stream.str());
392  } else {
393  dumpReport(accumulate_, stream.str());
394  }
395  }
396 
397  if (serviceBy_ == EVERY_LUMI) {
399  updateService(*luminosityBlockCache(lumi.index()));
400  } else {
402  }
403  }
404 }
405 
406 // ------------ method called to produce the data ------------
408  // accumulation of statistics event by event
409 
410  using namespace std;
411  using namespace edm;
412 
413  auto& accumulate = chooseAccumulate(iEvent.getLuminosityBlock().index());
414  if (resetBy_ == EVERY_EVENT) {
415  //NOTE if we have reportBy == lumi/run/job then we will report the last
416  // event in each lumi/run/job
417  accumulate.reset();
418  }
419 
420  accumulate.nEvents_++;
421 
422  // get hold of TriggerResults
424  iEvent.getByToken(hlTriggerResultsToken_, HLTR);
425  if (HLTR.isValid()) {
426  if (HLTR->wasrun())
427  accumulate.nWasRun_++;
428  const bool accept(HLTR->accept());
429  LogDebug("HLTrigReport") << "HLT TriggerResults decision: " << accept;
430  if (accept)
431  ++accumulate.nAccept_;
432  if (HLTR->error())
433  accumulate.nErrors_++;
434  } else {
435  LogDebug("HLTrigReport") << "HLT TriggerResults with label [" + hlTriggerResults_.encode() + "] not found!";
436  accumulate.nErrors_++;
437  return;
438  }
439 
440  // HLTConfigProvider not configured - cannot produce any detailed statistics
441  if (not configured_)
442  return;
443 
444  // decision for each HL algorithm
445  const unsigned int n(hlNames_.size());
446  bool acceptedByPrevoiusPaths = false;
447  for (unsigned int i = 0; i != n; ++i) {
448  if (HLTR->wasrun(i))
449  accumulate.hlWasRun_[i]++;
450  if (HLTR->accept(i)) {
451  acceptedByPrevoiusPaths = true;
452  accumulate.hlAccept_[i]++;
453  }
454  if (acceptedByPrevoiusPaths)
455  accumulate.hlAccTot_[i]++;
456  if (HLTR->error(i))
457  accumulate.hlErrors_[i]++;
458  const int index(static_cast<int>(HLTR->index(i)));
459  if (HLTR->accept(i)) {
460  if (index >= posL1s_[i])
461  accumulate.hltL1s_[i]++;
462  if (index >= posPre_[i])
463  accumulate.hltPre_[i]++;
464  } else {
465  if (index > posL1s_[i])
466  accumulate.hltL1s_[i]++;
467  if (index > posPre_[i])
468  accumulate.hltPre_[i]++;
469  }
470  }
471 
472  // calculate accumulation of accepted events by a path within a dataset
473  std::vector<bool> acceptedByDS(hlIndex_.size(), false);
474  for (size_t ds = 0; ds < hlIndex_.size(); ++ds) {
475  for (size_t p = 0; p < hlIndex_[ds].size(); ++p) {
476  if (acceptedByDS[ds] or HLTR->accept(hlIndex_[ds][p])) {
477  acceptedByDS[ds] = true;
478  accumulate.hlAccTotDS_[ds][p]++;
479  }
480  }
481  if (acceptedByDS[ds])
482  accumulate.hlAllTotDS_[ds]++;
483  }
484 
485  // calculate accumulation of accepted events by a dataset within a stream
486  for (size_t s = 0; s < dsIndex_.size(); ++s) {
487  bool acceptedByS = false;
488  for (size_t ds = 0; ds < dsIndex_[s].size(); ++ds) {
489  if (acceptedByS or acceptedByDS[dsIndex_[s][ds]]) {
490  acceptedByS = true;
491  accumulate.dsAccTotS_[s][ds]++;
492  }
493  }
494  if (acceptedByS)
495  accumulate.dsAllTotS_[s]++;
496  }
497 
498  if (reportBy_ == EVERY_EVENT) {
499  std::stringstream stream;
500  stream << "Summary for Run " << iEvent.run() << ", LumiSection " << iEvent.luminosityBlock() << ", Event "
501  << iEvent.id();
502  dumpReport(accumulate, stream.str());
503  }
504  if (serviceBy_ == EVERY_EVENT) {
505  updateService(accumulate);
506  }
507 }
508 
509 void HLTrigReport::updateService(Accumulate const& accumulate) const {
511  if (s) {
512  s->setDatasetCounts(accumulate.hlAllTotDS_);
513  s->setStreamCounts(accumulate.dsAllTotS_);
514  }
515 }
516 
518  std::string const& header /* = std::string() */) const {
519  // final printout of accumulated statistics
520 
521  using namespace std;
522  using namespace edm;
523  const unsigned int n(hlNames_.size());
524 
525  if ((n == 0) and (accumulate.nEvents_ == 0))
526  return;
527 
528  LogVerbatim("HLTrigReport") << dec << endl;
529  LogVerbatim("HLTrigReport") << "HLT-Report "
530  << "---------- Event Summary ------------" << endl;
531  if (not header.empty())
532  LogVerbatim("HLTrigReport") << "HLT-Report " << header << endl;
533  LogVerbatim("HLTrigReport") << "HLT-Report"
534  << " Events total = " << accumulate.nEvents_ << " wasrun = " << accumulate.nWasRun_
535  << " passed = " << accumulate.nAccept_ << " errors = " << accumulate.nErrors_ << endl;
536 
537  // HLTConfigProvider not configured - cannot produce any detailed statistics
538  if (not configured_)
539  return;
540 
541  double scale = accumulate.hlAccept_[refIndex_] > 0 ? refRate_ / accumulate.hlAccept_[refIndex_] : 0.;
542  double alpha = 1 - (1.0 - .6854) / 2; // for the Clopper-Pearson 68% CI
543 
544  LogVerbatim("HLTrigReport") << endl;
545  LogVerbatim("HLTrigReport") << "HLT-Report "
546  << "---------- HLTrig Summary ------------" << endl;
547  LogVerbatim("HLTrigReport") << "HLT-Report " << right << setw(7) << "HLT #"
548  << " " << right << setw(7) << "WasRun"
549  << " " << right << setw(7) << "L1S"
550  << " " << right << setw(7) << "Pre"
551  << " " << right << setw(7) << "HLT"
552  << " " << right << setw(9) << "%L1sPre"
553  << " " << right << setw(7) << "Rate"
554  << " " << right << setw(7) << "RateHi"
555  << " " << right << setw(7) << "Errors"
556  << " "
557  << "Name" << endl;
558 
559  if (n > 0) {
560  for (unsigned int i = 0; i != n; ++i) {
561  LogVerbatim("HLTrigReport")
562  << "HLT-Report " << right << setw(7) << i << " " << right << setw(7) << accumulate.hlWasRun_[i] << " "
563  << right << setw(7) << accumulate.hltL1s_[i] << " " << right << setw(7) << accumulate.hltPre_[i] << " "
564  << right << setw(7) << accumulate.hlAccept_[i] << " " << right << setw(9) << fixed << setprecision(5)
565  << static_cast<float>(100 * accumulate.hlAccept_[i]) / static_cast<float>(max(accumulate.hltPre_[i], 1u))
566  << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[i] << " " << right
567  << setw(7) << fixed << setprecision(1)
568  << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i] > 0)
569  ? refRate_ * ROOT::Math::beta_quantile(alpha,
570  accumulate.hlAccept_[i] + 1,
571  accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i])
572  : 0)
573  << " " << right << setw(7) << accumulate.hlErrors_[i] << " " << hlNames_[i] << endl;
574  }
575  }
576 
577  LogVerbatim("HLTrigRprtTt") << endl;
578  LogVerbatim("HLTrigRprtTt") << "HLT-Report "
579  << "---------- HLTrig Summary ------------" << endl;
580  LogVerbatim("HLTrigRprtTt") << "HLT-Report " << right << setw(7) << "HLT #"
581  << " " << right << setw(7) << "WasRun"
582  << " " << right << setw(7) << "L1S"
583  << " " << right << setw(7) << "Pre"
584  << " " << right << setw(7) << "HLT"
585  << " " << right << setw(9) << "%L1sPre"
586  << " " << right << setw(7) << "Rate"
587  << " " << right << setw(7) << "RateHi"
588  << " " << right << setw(7) << "HLTtot"
589  << " " << right << setw(7) << "RateTot"
590  << " " << right << setw(7) << "Errors"
591  << " "
592  << "Name" << endl;
593 
594  if (n > 0) {
595  for (unsigned int i = 0; i != n; ++i) {
596  LogVerbatim("HLTrigRprtTt")
597  << "HLT-Report " << right << setw(7) << i << " " << right << setw(7) << accumulate.hlWasRun_[i] << " "
598  << right << setw(7) << accumulate.hltL1s_[i] << " " << right << setw(7) << accumulate.hltPre_[i] << " "
599  << right << setw(7) << accumulate.hlAccept_[i] << " " << right << setw(9) << fixed << setprecision(5)
600  << static_cast<float>(100 * accumulate.hlAccept_[i]) / static_cast<float>(max(accumulate.hltPre_[i], 1u))
601  << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[i] << " " << right
602  << setw(7) << fixed << setprecision(1)
603  << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i] > 0)
604  ? refRate_ * ROOT::Math::beta_quantile(alpha,
605  accumulate.hlAccept_[i] + 1,
606  accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i])
607  : 0)
608  << " " << right << setw(7) << accumulate.hlAccTot_[i] << " " << right << setw(7) << fixed << setprecision(1)
609  << scale * accumulate.hlAccTot_[i] << " " << right << setw(7) << accumulate.hlErrors_[i] << " " << hlNames_[i]
610  << endl;
611  }
612 
613  // now for each dataset
614  for (size_t ds = 0; ds < hlIndex_.size(); ++ds) {
615  LogVerbatim("HLTrigRprtPD") << endl;
616  LogVerbatim("HLTrigRprtPD") << "HLT-Report "
617  << "---------- Dataset Summary: " << datasetNames_[ds] << " ------------"
618  << accumulate.hlAllTotDS_[ds] << endl;
619  LogVerbatim("HLTrigRprtPD") << "HLT-Report " << right << setw(7) << "HLT #"
620  << " " << right << setw(7) << "WasRun"
621  << " " << right << setw(7) << "L1S"
622  << " " << right << setw(7) << "Pre"
623  << " " << right << setw(7) << "HLT"
624  << " " << right << setw(9) << "%L1sPre"
625  << " " << right << setw(7) << "Rate"
626  << " " << right << setw(7) << "RateHi"
627  << " " << right << setw(7) << "HLTtot"
628  << " " << right << setw(7) << "RateTot"
629  << " " << right << setw(7) << "Errors"
630  << " "
631  << "Name" << endl;
632  for (size_t p = 0; p < hlIndex_[ds].size(); ++p) {
633  LogVerbatim("HLTrigRprtPD")
634  << "HLT-Report " << right << setw(7) << p << " " << right << setw(7)
635  << accumulate.hlWasRun_[hlIndex_[ds][p]] << " " << right << setw(7) << accumulate.hltL1s_[hlIndex_[ds][p]]
636  << " " << right << setw(7) << accumulate.hltPre_[hlIndex_[ds][p]] << " " << right << setw(7)
637  << accumulate.hlAccept_[hlIndex_[ds][p]] << " " << right << setw(9) << fixed << setprecision(5)
638  << static_cast<float>(100 * accumulate.hlAccept_[hlIndex_[ds][p]]) /
639  static_cast<float>(max(accumulate.hltPre_[hlIndex_[ds][p]], 1u))
640  << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[hlIndex_[ds][p]]
641  << " " << right << setw(7) << fixed << setprecision(1)
642  << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[hlIndex_[ds][p]] > 0)
643  ? refRate_ * ROOT::Math::beta_quantile(
644  alpha,
645  accumulate.hlAccept_[hlIndex_[ds][p]] + 1,
646  accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[hlIndex_[ds][p]])
647  : 0)
648  << " " << right << setw(7) << accumulate.hlAccTotDS_[ds][p] << " " << right << setw(7) << fixed
649  << setprecision(1) << scale * accumulate.hlAccTotDS_[ds][p] << " " << right << setw(7)
650  << accumulate.hlErrors_[hlIndex_[ds][p]] << " " << hlNames_[hlIndex_[ds][p]] << endl;
651  }
652  }
653 
654  // now for each stream
655  for (size_t s = 0; s < dsIndex_.size(); ++s) {
656  LogVerbatim("HLTrigRprtST") << endl;
657  LogVerbatim("HLTrigRprtST") << "HLT-Report "
658  << "---------- Stream Summary: " << streamNames_[s] << " ------------"
659  << accumulate.dsAllTotS_[s] << endl;
660  LogVerbatim("HLTrigRprtST") << "HLT-Report " << right << setw(10) << "Dataset #"
661  << " " << right << setw(10) << "Individual"
662  << " " << right << setw(10) << "Total"
663  << " " << right << setw(10) << "Rate"
664  << " " << right << setw(10) << "RateHi"
665  << " " << right << setw(10) << "RateTot"
666  << " "
667  << "Name" << endl;
668  for (size_t ds = 0; ds < dsIndex_[s].size(); ++ds) {
669  unsigned int acceptedDS = accumulate.hlAccTotDS_[dsIndex_[s][ds]][hlIndex_[dsIndex_[s][ds]].size() - 1];
670  LogVerbatim("HLTrigRprtST")
671  << "HLT-Report " << right << setw(10) << ds << " " << right << setw(10) << acceptedDS << " " << right
672  << setw(10) << accumulate.dsAccTotS_[s][ds] << " " << right << setw(10) << fixed << setprecision(1)
673  << scale * acceptedDS << " " << right << setw(10) << fixed << setprecision(1)
674  << ((accumulate.hlAccept_[refIndex_] - acceptedDS > 0)
675  ? refRate_ *
676  ROOT::Math::beta_quantile(alpha, acceptedDS + 1, accumulate.hlAccept_[refIndex_] - acceptedDS)
677  : 0)
678  << " " << right << setw(10) << fixed << setprecision(1) << scale * accumulate.dsAccTotS_[s][ds] << " "
679  << datasetNames_[dsIndex_[s][ds]] << endl;
680  }
681  }
682 
683  } else {
684  LogVerbatim("HLTrigReport") << "HLT-Report - No HLT paths found!" << endl;
685  }
686 
687  LogVerbatim("HLTrigReport") << endl;
688  LogVerbatim("HLTrigReport") << "HLT-Report end!" << endl;
689  LogVerbatim("HLTrigReport") << endl;
690 
691  return;
692 }
693 
694 // declare this class as a framework plugin
size
Write out results.
std::vector< std::vector< std::string > > datasetContents_
Definition: HLTrigReport.h:125
unsigned int index(const unsigned int i) const
Get index (slot position) of module giving the decision of the ith path.
bool isCustomStreams_
Definition: HLTrigReport.h:131
bool accept() const
Has at least one path accepted the event?
Log< level::Info, true > LogVerbatim
std::vector< unsigned int > hlAllTotDS_
Definition: HLTrigReport.h:52
bool isCustomDatasets_
Definition: HLTrigReport.h:126
hltrigreport::Accumulate & chooseAccumulate(edm::LuminosityBlockIndex index)
Definition: HLTrigReport.h:96
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< std::string > hlNames_
Definition: HLTrigReport.h:116
const ReportEvery resetBy_
Definition: HLTrigReport.h:137
const ReportEvery reportBy_
Definition: HLTrigReport.h:136
bool error() const
Has any path encountered an error (exception)
std::string encode() const
Definition: InputTag.cc:159
std::vector< unsigned int > hlAccTot_
Definition: HLTrigReport.h:47
std::vector< std::vector< unsigned int > > hlAccTotDS_
Definition: HLTrigReport.h:51
const std::string & moduleType(const std::string &module) const
C++ class name of module.
std::vector< std::string > datasetNames_
Definition: HLTrigReport.h:124
bool useLumiCache() const
Definition: HLTrigReport.h:102
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
assert(be >=bs)
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
const std::vector< std::string > & datasetNames() const
static ReportEvery decode(const std::string &value)
Definition: HLTrigReport.cc:24
void beginRun(edm::Run const &, edm::EventSetup const &) override
void updateService(Accumulate const &accumulate) const
void updateConfigCache()
std::vector< int > posL1s_
Definition: HLTrigReport.h:121
bool wasrun() const
Was at least one path run?
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > posPre_
Definition: HLTrigReport.h:122
char const * label
void endRun(edm::Run const &, edm::EventSetup const &) override
int iEvent
Definition: GenABIO.cc:224
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
const std::vector< std::string > & moduleLabels(unsigned int trigger) const
label(s) of module(s) on a trigger path
std::vector< std::string > getParameterNamesForType(bool trackiness=true) const
Definition: ParameterSet.h:180
std::vector< std::string > streamNames_
Definition: HLTrigReport.h:129
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
Definition: value.py:1
unsigned int refIndex_
Definition: HLTrigReport.h:133
void dumpReport(hltrigreport::Accumulate const &accumulate, std::string const &header=std::string()) const
const std::vector< std::vector< std::string > > & streamContents() const
names of datasets for all streams
const edm::EDGetTokenT< edm::TriggerResults > hlTriggerResultsToken_
Definition: HLTrigReport.h:113
std::string refPath_
Definition: HLTrigReport.h:132
hltrigreport::Accumulate accumulate_
Definition: HLTrigReport.h:141
HLTConfigProvider hltConfig_
Definition: HLTrigReport.h:139
std::vector< std::vector< std::string > > streamContents_
Definition: HLTrigReport.h:130
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
std::shared_ptr< Accumulate > globalBeginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::vector< std::string > & streamNames() const
const std::vector< std::vector< std::string > > & datasetContents() const
names of trigger paths for all datasets
void globalEndLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
std::vector< std::vector< unsigned int > > dsAccTotS_
Definition: HLTrigReport.h:54
const std::vector< std::string > & triggerNames() const
names of trigger paths
bool isValid() const
Definition: HandleBase.h:70
const double refRate_
Definition: HLTrigReport.h:134
HLT enums.
~HLTrigReport() override
void beginJob() override
HLTrigReport(const edm::ParameterSet &)
std::vector< std::vector< unsigned int > > hlIndex_
Definition: HLTrigReport.h:119
void accumulate(Accumulate const &)
Definition: HLTrigReport.cc:87
std::vector< unsigned int > dsAllTotS_
Definition: HLTrigReport.h:55
std::vector< unsigned int > hlWasRun_
Definition: HLTrigReport.h:43
void endJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< unsigned int > hlAccept_
Definition: HLTrigReport.h:46
const edm::InputTag hlTriggerResults_
Definition: HLTrigReport.h:112
Log< level::Warning, false > LogWarning
bool decode(bool &, std::string_view)
Definition: types.cc:72
std::string const & process() const
Definition: InputTag.h:40
const std::vector< std::string > & datasetNames() const
const ReportEvery serviceBy_
Definition: HLTrigReport.h:138
const std::vector< std::string > & streamNames() const
hltrigreport::Accumulate Accumulate
Definition: HLTrigReport.h:69
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
std::vector< std::vector< unsigned int > > dsIndex_
Definition: HLTrigReport.h:128
std::vector< unsigned int > hltPre_
Definition: HLTrigReport.h:45
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
bool readAfterLumi() const
Definition: HLTrigReport.h:103
std::vector< unsigned int > hlErrors_
Definition: HLTrigReport.h:48
std::vector< unsigned int > hltL1s_
Definition: HLTrigReport.h:44
#define LogDebug(id)