Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:23

0001 /** \class HLTrigReport
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Martin Grunewald
0007  *
0008  */
0009 
0010 #include <iomanip>
0011 #include <cstring>
0012 #include <sstream>
0013 
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "FWCore/Common/interface/TriggerNames.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 #include "HLTrigReportService.h"
0020 #include "Math/QuantFuncMathCore.h"
0021 
0022 #include "HLTrigReport.h"
0023 
0024 HLTrigReport::ReportEvery HLTrigReport::decode(const std::string& value) {
0025   if (value == "never")
0026     return NEVER;
0027 
0028   if (value == "job")
0029     return EVERY_JOB;
0030 
0031   if (value == "run")
0032     return EVERY_RUN;
0033 
0034   if (value == "lumi")
0035     return EVERY_LUMI;
0036 
0037   if (value == "event")
0038     return EVERY_EVENT;
0039 
0040   throw cms::Exception("Configuration") << "Invalid option value \"" << value
0041                                         << "\". Legal values are \"job\", \"run\", \"lumi\", \"event\" and \"never\".";
0042 }
0043 
0044 //
0045 // constructors and destructor
0046 //
0047 hltrigreport::Accumulate::Accumulate()
0048     : nEvents_(0),
0049       nWasRun_(0),
0050       nAccept_(0),
0051       nErrors_(0),
0052       hlWasRun_(0),
0053       hltL1s_(0),
0054       hltPre_(0),
0055       hlAccept_(0),
0056       hlAccTot_(0),
0057       hlErrors_(0),
0058       hlAccTotDS_(0),
0059       dsAccTotS_(0) {}
0060 
0061 hltrigreport::Accumulate::Accumulate(size_t numHLNames,
0062                                      std::vector<std::vector<unsigned int> > const& hlIndex,
0063                                      std::vector<std::vector<unsigned int> > const& dsIndex)
0064     : nEvents_(0),
0065       nWasRun_(0),
0066       nAccept_(0),
0067       nErrors_(0),
0068       hlWasRun_(numHLNames),
0069       hltL1s_(numHLNames),
0070       hltPre_(numHLNames),
0071       hlAccept_(numHLNames),
0072       hlAccTot_(numHLNames),
0073       hlErrors_(numHLNames),
0074       hlAccTotDS_(hlIndex.size()),
0075       hlAllTotDS_(hlIndex.size()),
0076       dsAccTotS_(dsIndex.size()),
0077       dsAllTotS_(dsIndex.size()) {
0078   for (size_t ds = 0; ds < hlIndex.size(); ++ds) {
0079     hlAccTotDS_[ds].resize(hlIndex[ds].size());
0080   }
0081 
0082   for (size_t s = 0; s < dsIndex.size(); ++s) {
0083     dsAccTotS_[s].resize(dsIndex[s].size());
0084   }
0085 }
0086 
0087 void hltrigreport::Accumulate::accumulate(hltrigreport::Accumulate const& iOther) {
0088   nEvents_ += iOther.nEvents_;
0089   nWasRun_ += iOther.nWasRun_;
0090   nAccept_ += iOther.nAccept_;
0091   nErrors_ += iOther.nErrors_;
0092 
0093   auto vsum = [](auto& to, auto const& from) {
0094     for (size_t i = 0; i < from.size(); ++i) {
0095       to[i] += from[i];
0096     }
0097   };
0098 
0099   assert(hlWasRun_.size() == iOther.hlWasRun_.size());
0100   vsum(hlWasRun_, iOther.hlWasRun_);
0101   vsum(hltL1s_, iOther.hltL1s_);
0102   vsum(hltPre_, iOther.hltPre_);
0103   vsum(hlAccept_, iOther.hlAccept_);
0104   vsum(hlAccTot_, iOther.hlAccTot_);
0105   vsum(hlErrors_, iOther.hlErrors_);
0106 
0107   assert(hlAllTotDS_.size() == iOther.hlAllTotDS_.size());
0108   vsum(hlAllTotDS_, iOther.hlAllTotDS_);
0109   vsum(dsAllTotS_, iOther.dsAllTotS_);
0110 
0111   auto vvsum = [](auto& to, auto const& from) {
0112     for (size_t i = 0; i < from.size(); ++i) {
0113       assert(from[i].size() == to[i].size());
0114       for (size_t j = 0; j < from[i].size(); ++j) {
0115         to[i][j] += from[i][j];
0116       }
0117     }
0118   };
0119 
0120   vvsum(hlAccTotDS_, iOther.hlAccTotDS_);
0121   vvsum(dsAccTotS_, iOther.dsAccTotS_);
0122 }
0123 
0124 void hltrigreport::Accumulate::reset() {
0125   nEvents_ = 0;
0126   nWasRun_ = 0;
0127   nAccept_ = 0;
0128   nErrors_ = 0;
0129 
0130   auto vreset = [](auto& to) { std::fill(to.begin(), to.end(), 0); };
0131 
0132   vreset(hlWasRun_);
0133   vreset(hltL1s_);
0134   vreset(hltPre_);
0135   vreset(hlAccept_);
0136   vreset(hlAccTot_);
0137   vreset(hlErrors_);
0138 
0139   vreset(hlAllTotDS_);
0140   vreset(dsAllTotS_);
0141 
0142   auto vvreset = [&vreset](auto& to) {
0143     for (auto& e : to) {
0144       vreset(e);
0145     }
0146   };
0147 
0148   vvreset(hlAccTotDS_);
0149   vvreset(dsAccTotS_);
0150 }
0151 
0152 HLTrigReport::HLTrigReport(const edm::ParameterSet& iConfig)
0153     : hlTriggerResults_(iConfig.getParameter<edm::InputTag>("HLTriggerResults")),
0154       hlTriggerResultsToken_(consumes<edm::TriggerResults>(hlTriggerResults_)),
0155       configured_(false),
0156       hlNames_(0),
0157       hlIndex_(0),
0158       posL1s_(0),
0159       posPre_(0),
0160       datasetNames_(0),
0161       datasetContents_(0),
0162       isCustomDatasets_(false),
0163       dsIndex_(0),
0164       streamNames_(0),
0165       streamContents_(0),
0166       isCustomStreams_(false),
0167       refPath_("HLTriggerFinalPath"),
0168       refIndex_(0),
0169       refRate_(iConfig.getUntrackedParameter<double>("ReferenceRate", 100.0)),
0170       reportBy_(decode(iConfig.getUntrackedParameter<std::string>("reportBy", "job"))),
0171       resetBy_(decode(iConfig.getUntrackedParameter<std::string>("resetBy", "never"))),
0172       serviceBy_(decode(iConfig.getUntrackedParameter<std::string>("serviceBy", "never"))),
0173       hltConfig_() {
0174   const edm::ParameterSet customDatasets(
0175       iConfig.getUntrackedParameter<edm::ParameterSet>("CustomDatasets", edm::ParameterSet()));
0176   isCustomDatasets_ = (customDatasets != edm::ParameterSet());
0177   if (isCustomDatasets_) {
0178     datasetNames_ = customDatasets.getParameterNamesForType<std::vector<std::string> >();
0179     for (std::vector<std::string>::const_iterator name = datasetNames_.begin(); name != datasetNames_.end(); name++) {
0180       datasetContents_.push_back(customDatasets.getParameter<std::vector<std::string> >(*name));
0181     }
0182   }
0183 
0184   const edm::ParameterSet customStreams(
0185       iConfig.getUntrackedParameter<edm::ParameterSet>("CustomStreams", edm::ParameterSet()));
0186   isCustomStreams_ = (customStreams != edm::ParameterSet());
0187   if (isCustomStreams_) {
0188     streamNames_ = customStreams.getParameterNamesForType<std::vector<std::string> >();
0189     for (std::vector<std::string>::const_iterator name = streamNames_.begin(); name != streamNames_.end(); name++) {
0190       streamContents_.push_back(customStreams.getParameter<std::vector<std::string> >(*name));
0191     }
0192   }
0193 
0194   refPath_ = iConfig.getUntrackedParameter<std::string>("ReferencePath", "HLTriggerFinalPath");
0195   refIndex_ = 0;
0196 
0197   LogDebug("HLTrigReport") << "HL TiggerResults: " + hlTriggerResults_.encode()
0198                            << " using reference path and rate: " + refPath_ + " " << refRate_;
0199 
0200   if (serviceBy_ != NEVER and edm::Service<HLTrigReportService>()) {
0201     edm::Service<HLTrigReportService>()->registerModule(this);
0202   }
0203 }
0204 
0205 HLTrigReport::~HLTrigReport() = default;
0206 
0207 void HLTrigReport::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0208   edm::ParameterSetDescription desc;
0209   desc.add<edm::InputTag>("HLTriggerResults", edm::InputTag("TriggerResults", "", "HLT"));
0210   desc.addUntracked<std::string>("reportBy", "job");
0211   desc.addUntracked<std::string>("resetBy", "never");
0212   desc.addUntracked<std::string>("serviceBy", "never");
0213 
0214   edm::ParameterSetDescription customDatasetsParameters;
0215   desc.addUntracked<edm::ParameterSetDescription>("CustomDatasets", customDatasetsParameters);
0216   edm::ParameterSetDescription customStreamsParameters;
0217   desc.addUntracked<edm::ParameterSetDescription>("CustomStreams", customStreamsParameters);
0218   desc.addUntracked<std::string>("ReferencePath", "HLTriggerFinalPath");
0219   desc.addUntracked<double>("ReferenceRate", 100.0);
0220 
0221   descriptions.add("hltTrigReport", desc);
0222 }
0223 
0224 //
0225 // member functions
0226 //
0227 
0228 const std::vector<std::string>& HLTrigReport::datasetNames() const { return datasetNames_; }
0229 const std::vector<std::string>& HLTrigReport::streamNames() const { return streamNames_; }
0230 
0231 void HLTrigReport::updateConfigCache() {
0232   // update trigger names
0233   hlNames_ = hltConfig_.triggerNames();
0234 
0235   const unsigned int n = hlNames_.size();
0236 
0237   // find the positions of seeding and prescaler modules
0238   posL1s_.resize(n);
0239   posPre_.resize(n);
0240   for (unsigned int i = 0; i < n; ++i) {
0241     posL1s_[i] = -1;
0242     posPre_[i] = -1;
0243     const std::vector<std::string>& moduleLabels(hltConfig_.moduleLabels(i));
0244     for (unsigned int j = 0; j < moduleLabels.size(); ++j) {
0245       const std::string& label = hltConfig_.moduleType(moduleLabels[j]);
0246       if (label == "HLTLevel1GTSeed")
0247         posL1s_[i] = j;
0248       else if (label == "HLTPrescaler")
0249         posPre_[i] = j;
0250     }
0251   }
0252 
0253   // if not overridden, reload the datasets and streams
0254   if (not isCustomDatasets_) {
0255     datasetNames_ = hltConfig_.datasetNames();
0256     datasetContents_ = hltConfig_.datasetContents();
0257   }
0258   if (not isCustomStreams_) {
0259     streamNames_ = hltConfig_.streamNames();
0260     streamContents_ = hltConfig_.streamContents();
0261   }
0262 
0263   // fill the matrices of hlIndex_, hlAccTotDS_
0264   hlIndex_.clear();
0265   hlIndex_.resize(datasetNames_.size());
0266   for (unsigned int ds = 0; ds < datasetNames_.size(); ds++) {
0267     unsigned int size = datasetContents_[ds].size();
0268     hlIndex_[ds].reserve(size);
0269     for (unsigned int p = 0; p < size; ++p) {
0270       unsigned int i = hltConfig_.triggerIndex(datasetContents_[ds][p]);
0271       if (i < n) {
0272         hlIndex_[ds].push_back(i);
0273       }
0274     }
0275   }
0276 
0277   // fill the matrices of dsIndex_, dsAccTotS_
0278   dsIndex_.clear();
0279   dsIndex_.resize(streamNames_.size());
0280   for (unsigned int s = 0; s < streamNames_.size(); ++s) {
0281     unsigned int size = streamContents_[s].size();
0282     dsIndex_.reserve(size);
0283     for (unsigned int ds = 0; ds < size; ++ds) {
0284       unsigned int i = 0;
0285       for (; i < datasetNames_.size(); i++)
0286         if (datasetNames_[i] == streamContents_[s][ds])
0287           break;
0288       // report only datasets that have at least one path otherwise crash
0289       if (i < datasetNames_.size() and !hlIndex_[i].empty()) {
0290         dsIndex_[s].push_back(i);
0291       }
0292     }
0293   }
0294 
0295   // if needed, update the reference path
0296   refIndex_ = hltConfig_.triggerIndex(refPath_);
0297   if (refIndex_ >= n) {
0298     refIndex_ = 0;
0299     edm::LogWarning("HLTrigReport") << "Requested reference path '" + refPath_ + "' not in HLT menu. "
0300                                     << "Using HLTriggerFinalPath instead.";
0301     refPath_ = "HLTriggerFinalPath";
0302     refIndex_ = hltConfig_.triggerIndex(refPath_);
0303     if (refIndex_ >= n) {
0304       refIndex_ = 0;
0305       edm::LogWarning("HLTrigReport") << "Requested reference path '" + refPath_ + "' not in HLT menu. "
0306                                       << "Using first path in table (index=0) instead.";
0307     }
0308   }
0309 
0310   if (serviceBy_ != NEVER and edm::Service<HLTrigReportService>()) {
0311     edm::Service<HLTrigReportService>()->setDatasetNames(datasetNames_);
0312     edm::Service<HLTrigReportService>()->setStreamNames(streamNames_);
0313   }
0314 }
0315 
0316 void HLTrigReport::reset() { accumulate_.reset(); }
0317 
0318 void HLTrigReport::beginJob() {
0319   if (resetBy_ == EVERY_JOB)
0320     reset();
0321 }
0322 
0323 void HLTrigReport::endJob() {
0324   if (reportBy_ == EVERY_JOB)
0325     dumpReport(accumulate_, "Summary for Job");
0326   if (serviceBy_ == EVERY_JOB) {
0327     updateService(accumulate_);
0328   }
0329 }
0330 
0331 void HLTrigReport::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
0332   bool changed = true;
0333   if (hltConfig_.init(iRun, iSetup, hlTriggerResults_.process(), changed)) {
0334     configured_ = true;
0335     if (changed) {
0336       dumpReport(accumulate_, "Summary for this HLT table");
0337       updateConfigCache();
0338       accumulate_ = Accumulate(hlNames_.size(), hlIndex_, dsIndex_);
0339     }
0340   } else {
0341     dumpReport(accumulate_, "Summary for this HLT table");
0342     // cannot initialize the HLT menu - reset and clear all counters and tables
0343     configured_ = false;
0344 
0345     accumulate_ = hltrigreport::Accumulate();
0346   }
0347 
0348   if (resetBy_ == EVERY_RUN)
0349     reset();
0350 }
0351 
0352 void HLTrigReport::endRun(edm::Run const& run, edm::EventSetup const& setup) {
0353   if (reportBy_ == EVERY_RUN) {
0354     std::stringstream stream;
0355     stream << "Summary for Run " << run.run();
0356     dumpReport(accumulate_, stream.str());
0357   }
0358   if (serviceBy_ == EVERY_RUN) {
0359     updateService(accumulate_);
0360   }
0361 }
0362 
0363 std::shared_ptr<HLTrigReport::Accumulate> HLTrigReport::globalBeginLuminosityBlock(edm::LuminosityBlock const& lumi,
0364                                                                                    edm::EventSetup const& setup) const {
0365   if (useLumiCache()) {
0366     if (not configured_) {
0367       return std::make_shared<Accumulate>();
0368     }
0369     return std::make_shared<Accumulate>(hlNames_.size(), hlIndex_, dsIndex_);
0370   }
0371   return std::shared_ptr<Accumulate>();
0372 }
0373 
0374 void HLTrigReport::globalEndLuminosityBlock(edm::LuminosityBlock const& lumi, edm::EventSetup const& setup) {
0375   if (not useLumiCache()) {
0376     return;
0377   }
0378 
0379   if (resetBy_ == EVERY_LUMI and readAfterLumi()) {
0380     //we will be reporting the last processed lumi
0381     accumulate_ = std::move(*luminosityBlockCache(lumi.index()));
0382   } else if (resetBy_ == NEVER or resetBy_ == EVERY_RUN or resetBy_ == EVERY_JOB or readAfterLumi()) {
0383     //we need to add this lumi's info to the longer lived accumulation
0384     accumulate_.accumulate(*luminosityBlockCache(lumi.index()));
0385   }
0386 
0387   if (reportBy_ == EVERY_LUMI) {
0388     std::stringstream stream;
0389     stream << "Summary for Run " << lumi.run() << ", LumiSection " << lumi.luminosityBlock();
0390     if (resetBy_ == EVERY_LUMI or resetBy_ == EVERY_EVENT) {
0391       dumpReport(*luminosityBlockCache(lumi.index()), stream.str());
0392     } else {
0393       dumpReport(accumulate_, stream.str());
0394     }
0395   }
0396 
0397   if (serviceBy_ == EVERY_LUMI) {
0398     if (resetBy_ == EVERY_LUMI or resetBy_ == EVERY_EVENT) {
0399       updateService(*luminosityBlockCache(lumi.index()));
0400     } else {
0401       updateService(accumulate_);
0402     }
0403   }
0404 }
0405 
0406 // ------------ method called to produce the data  ------------
0407 void HLTrigReport::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0408   // accumulation of statistics event by event
0409 
0410   using namespace std;
0411   using namespace edm;
0412 
0413   auto& accumulate = chooseAccumulate(iEvent.getLuminosityBlock().index());
0414   if (resetBy_ == EVERY_EVENT) {
0415     //NOTE if we have reportBy == lumi/run/job then we will report the last
0416     // event in each lumi/run/job
0417     accumulate.reset();
0418   }
0419 
0420   accumulate.nEvents_++;
0421 
0422   // get hold of TriggerResults
0423   Handle<TriggerResults> HLTR;
0424   iEvent.getByToken(hlTriggerResultsToken_, HLTR);
0425   if (HLTR.isValid()) {
0426     if (HLTR->wasrun())
0427       accumulate.nWasRun_++;
0428     const bool accept(HLTR->accept());
0429     LogDebug("HLTrigReport") << "HLT TriggerResults decision: " << accept;
0430     if (accept)
0431       ++accumulate.nAccept_;
0432     if (HLTR->error())
0433       accumulate.nErrors_++;
0434   } else {
0435     LogDebug("HLTrigReport") << "HLT TriggerResults with label [" + hlTriggerResults_.encode() + "] not found!";
0436     accumulate.nErrors_++;
0437     return;
0438   }
0439 
0440   // HLTConfigProvider not configured - cannot produce any detailed statistics
0441   if (not configured_)
0442     return;
0443 
0444   // decision for each HL algorithm
0445   const unsigned int n(hlNames_.size());
0446   bool acceptedByPrevoiusPaths = false;
0447   for (unsigned int i = 0; i != n; ++i) {
0448     if (HLTR->wasrun(i))
0449       accumulate.hlWasRun_[i]++;
0450     if (HLTR->accept(i)) {
0451       acceptedByPrevoiusPaths = true;
0452       accumulate.hlAccept_[i]++;
0453     }
0454     if (acceptedByPrevoiusPaths)
0455       accumulate.hlAccTot_[i]++;
0456     if (HLTR->error(i))
0457       accumulate.hlErrors_[i]++;
0458     const int index(static_cast<int>(HLTR->index(i)));
0459     if (HLTR->accept(i)) {
0460       if (index >= posL1s_[i])
0461         accumulate.hltL1s_[i]++;
0462       if (index >= posPre_[i])
0463         accumulate.hltPre_[i]++;
0464     } else {
0465       if (index > posL1s_[i])
0466         accumulate.hltL1s_[i]++;
0467       if (index > posPre_[i])
0468         accumulate.hltPre_[i]++;
0469     }
0470   }
0471 
0472   // calculate accumulation of accepted events by a path within a dataset
0473   std::vector<bool> acceptedByDS(hlIndex_.size(), false);
0474   for (size_t ds = 0; ds < hlIndex_.size(); ++ds) {
0475     for (size_t p = 0; p < hlIndex_[ds].size(); ++p) {
0476       if (acceptedByDS[ds] or HLTR->accept(hlIndex_[ds][p])) {
0477         acceptedByDS[ds] = true;
0478         accumulate.hlAccTotDS_[ds][p]++;
0479       }
0480     }
0481     if (acceptedByDS[ds])
0482       accumulate.hlAllTotDS_[ds]++;
0483   }
0484 
0485   // calculate accumulation of accepted events by a dataset within a stream
0486   for (size_t s = 0; s < dsIndex_.size(); ++s) {
0487     bool acceptedByS = false;
0488     for (size_t ds = 0; ds < dsIndex_[s].size(); ++ds) {
0489       if (acceptedByS or acceptedByDS[dsIndex_[s][ds]]) {
0490         acceptedByS = true;
0491         accumulate.dsAccTotS_[s][ds]++;
0492       }
0493     }
0494     if (acceptedByS)
0495       accumulate.dsAllTotS_[s]++;
0496   }
0497 
0498   if (reportBy_ == EVERY_EVENT) {
0499     std::stringstream stream;
0500     stream << "Summary for Run " << iEvent.run() << ", LumiSection " << iEvent.luminosityBlock() << ", Event "
0501            << iEvent.id();
0502     dumpReport(accumulate, stream.str());
0503   }
0504   if (serviceBy_ == EVERY_EVENT) {
0505     updateService(accumulate);
0506   }
0507 }
0508 
0509 void HLTrigReport::updateService(Accumulate const& accumulate) const {
0510   edm::Service<HLTrigReportService> s;
0511   if (s) {
0512     s->setDatasetCounts(accumulate.hlAllTotDS_);
0513     s->setStreamCounts(accumulate.dsAllTotS_);
0514   }
0515 }
0516 
0517 void HLTrigReport::dumpReport(hltrigreport::Accumulate const& accumulate,
0518                               std::string const& header /* = std::string() */) const {
0519   // final printout of accumulated statistics
0520 
0521   using namespace std;
0522   using namespace edm;
0523   const unsigned int n(hlNames_.size());
0524 
0525   if ((n == 0) and (accumulate.nEvents_ == 0))
0526     return;
0527 
0528   LogVerbatim("HLTrigReport") << dec << endl;
0529   LogVerbatim("HLTrigReport") << "HLT-Report "
0530                               << "---------- Event  Summary ------------" << endl;
0531   if (not header.empty())
0532     LogVerbatim("HLTrigReport") << "HLT-Report " << header << endl;
0533   LogVerbatim("HLTrigReport") << "HLT-Report"
0534                               << " Events total = " << accumulate.nEvents_ << " wasrun = " << accumulate.nWasRun_
0535                               << " passed = " << accumulate.nAccept_ << " errors = " << accumulate.nErrors_ << endl;
0536 
0537   // HLTConfigProvider not configured - cannot produce any detailed statistics
0538   if (not configured_)
0539     return;
0540 
0541   double scale = accumulate.hlAccept_[refIndex_] > 0 ? refRate_ / accumulate.hlAccept_[refIndex_] : 0.;
0542   double alpha = 1 - (1.0 - .6854) / 2;  // for the Clopper-Pearson 68% CI
0543 
0544   LogVerbatim("HLTrigReport") << endl;
0545   LogVerbatim("HLTrigReport") << "HLT-Report "
0546                               << "---------- HLTrig Summary ------------" << endl;
0547   LogVerbatim("HLTrigReport") << "HLT-Report " << right << setw(7) << "HLT #"
0548                               << " " << right << setw(7) << "WasRun"
0549                               << " " << right << setw(7) << "L1S"
0550                               << " " << right << setw(7) << "Pre"
0551                               << " " << right << setw(7) << "HLT"
0552                               << " " << right << setw(9) << "%L1sPre"
0553                               << " " << right << setw(7) << "Rate"
0554                               << " " << right << setw(7) << "RateHi"
0555                               << " " << right << setw(7) << "Errors"
0556                               << " "
0557                               << "Name" << endl;
0558 
0559   if (n > 0) {
0560     for (unsigned int i = 0; i != n; ++i) {
0561       LogVerbatim("HLTrigReport")
0562           << "HLT-Report " << right << setw(7) << i << " " << right << setw(7) << accumulate.hlWasRun_[i] << " "
0563           << right << setw(7) << accumulate.hltL1s_[i] << " " << right << setw(7) << accumulate.hltPre_[i] << " "
0564           << right << setw(7) << accumulate.hlAccept_[i] << " " << right << setw(9) << fixed << setprecision(5)
0565           << static_cast<float>(100 * accumulate.hlAccept_[i]) / static_cast<float>(max(accumulate.hltPre_[i], 1u))
0566           << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[i] << " " << right
0567           << setw(7) << fixed << setprecision(1)
0568           << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i] > 0)
0569                   ? refRate_ * ROOT::Math::beta_quantile(alpha,
0570                                                          accumulate.hlAccept_[i] + 1,
0571                                                          accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i])
0572                   : 0)
0573           << " " << right << setw(7) << accumulate.hlErrors_[i] << " " << hlNames_[i] << endl;
0574     }
0575   }
0576 
0577   LogVerbatim("HLTrigRprtTt") << endl;
0578   LogVerbatim("HLTrigRprtTt") << "HLT-Report "
0579                               << "---------- HLTrig Summary ------------" << endl;
0580   LogVerbatim("HLTrigRprtTt") << "HLT-Report " << right << setw(7) << "HLT #"
0581                               << " " << right << setw(7) << "WasRun"
0582                               << " " << right << setw(7) << "L1S"
0583                               << " " << right << setw(7) << "Pre"
0584                               << " " << right << setw(7) << "HLT"
0585                               << " " << right << setw(9) << "%L1sPre"
0586                               << " " << right << setw(7) << "Rate"
0587                               << " " << right << setw(7) << "RateHi"
0588                               << " " << right << setw(7) << "HLTtot"
0589                               << " " << right << setw(7) << "RateTot"
0590                               << " " << right << setw(7) << "Errors"
0591                               << " "
0592                               << "Name" << endl;
0593 
0594   if (n > 0) {
0595     for (unsigned int i = 0; i != n; ++i) {
0596       LogVerbatim("HLTrigRprtTt")
0597           << "HLT-Report " << right << setw(7) << i << " " << right << setw(7) << accumulate.hlWasRun_[i] << " "
0598           << right << setw(7) << accumulate.hltL1s_[i] << " " << right << setw(7) << accumulate.hltPre_[i] << " "
0599           << right << setw(7) << accumulate.hlAccept_[i] << " " << right << setw(9) << fixed << setprecision(5)
0600           << static_cast<float>(100 * accumulate.hlAccept_[i]) / static_cast<float>(max(accumulate.hltPre_[i], 1u))
0601           << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[i] << " " << right
0602           << setw(7) << fixed << setprecision(1)
0603           << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i] > 0)
0604                   ? refRate_ * ROOT::Math::beta_quantile(alpha,
0605                                                          accumulate.hlAccept_[i] + 1,
0606                                                          accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[i])
0607                   : 0)
0608           << " " << right << setw(7) << accumulate.hlAccTot_[i] << " " << right << setw(7) << fixed << setprecision(1)
0609           << scale * accumulate.hlAccTot_[i] << " " << right << setw(7) << accumulate.hlErrors_[i] << " " << hlNames_[i]
0610           << endl;
0611     }
0612 
0613     // now for each dataset
0614     for (size_t ds = 0; ds < hlIndex_.size(); ++ds) {
0615       LogVerbatim("HLTrigRprtPD") << endl;
0616       LogVerbatim("HLTrigRprtPD") << "HLT-Report "
0617                                   << "---------- Dataset Summary: " << datasetNames_[ds] << " ------------"
0618                                   << accumulate.hlAllTotDS_[ds] << endl;
0619       LogVerbatim("HLTrigRprtPD") << "HLT-Report " << right << setw(7) << "HLT #"
0620                                   << " " << right << setw(7) << "WasRun"
0621                                   << " " << right << setw(7) << "L1S"
0622                                   << " " << right << setw(7) << "Pre"
0623                                   << " " << right << setw(7) << "HLT"
0624                                   << " " << right << setw(9) << "%L1sPre"
0625                                   << " " << right << setw(7) << "Rate"
0626                                   << " " << right << setw(7) << "RateHi"
0627                                   << " " << right << setw(7) << "HLTtot"
0628                                   << " " << right << setw(7) << "RateTot"
0629                                   << " " << right << setw(7) << "Errors"
0630                                   << " "
0631                                   << "Name" << endl;
0632       for (size_t p = 0; p < hlIndex_[ds].size(); ++p) {
0633         LogVerbatim("HLTrigRprtPD")
0634             << "HLT-Report " << right << setw(7) << p << " " << right << setw(7)
0635             << accumulate.hlWasRun_[hlIndex_[ds][p]] << " " << right << setw(7) << accumulate.hltL1s_[hlIndex_[ds][p]]
0636             << " " << right << setw(7) << accumulate.hltPre_[hlIndex_[ds][p]] << " " << right << setw(7)
0637             << accumulate.hlAccept_[hlIndex_[ds][p]] << " " << right << setw(9) << fixed << setprecision(5)
0638             << static_cast<float>(100 * accumulate.hlAccept_[hlIndex_[ds][p]]) /
0639                    static_cast<float>(max(accumulate.hltPre_[hlIndex_[ds][p]], 1u))
0640             << " " << right << setw(7) << fixed << setprecision(1) << scale * accumulate.hlAccept_[hlIndex_[ds][p]]
0641             << " " << right << setw(7) << fixed << setprecision(1)
0642             << ((accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[hlIndex_[ds][p]] > 0)
0643                     ? refRate_ * ROOT::Math::beta_quantile(
0644                                      alpha,
0645                                      accumulate.hlAccept_[hlIndex_[ds][p]] + 1,
0646                                      accumulate.hlAccept_[refIndex_] - accumulate.hlAccept_[hlIndex_[ds][p]])
0647                     : 0)
0648             << " " << right << setw(7) << accumulate.hlAccTotDS_[ds][p] << " " << right << setw(7) << fixed
0649             << setprecision(1) << scale * accumulate.hlAccTotDS_[ds][p] << " " << right << setw(7)
0650             << accumulate.hlErrors_[hlIndex_[ds][p]] << " " << hlNames_[hlIndex_[ds][p]] << endl;
0651       }
0652     }
0653 
0654     // now for each stream
0655     for (size_t s = 0; s < dsIndex_.size(); ++s) {
0656       LogVerbatim("HLTrigRprtST") << endl;
0657       LogVerbatim("HLTrigRprtST") << "HLT-Report "
0658                                   << "---------- Stream Summary: " << streamNames_[s] << " ------------"
0659                                   << accumulate.dsAllTotS_[s] << endl;
0660       LogVerbatim("HLTrigRprtST") << "HLT-Report " << right << setw(10) << "Dataset #"
0661                                   << " " << right << setw(10) << "Individual"
0662                                   << " " << right << setw(10) << "Total"
0663                                   << " " << right << setw(10) << "Rate"
0664                                   << " " << right << setw(10) << "RateHi"
0665                                   << " " << right << setw(10) << "RateTot"
0666                                   << " "
0667                                   << "Name" << endl;
0668       for (size_t ds = 0; ds < dsIndex_[s].size(); ++ds) {
0669         unsigned int acceptedDS = accumulate.hlAccTotDS_[dsIndex_[s][ds]][hlIndex_[dsIndex_[s][ds]].size() - 1];
0670         LogVerbatim("HLTrigRprtST")
0671             << "HLT-Report " << right << setw(10) << ds << " " << right << setw(10) << acceptedDS << " " << right
0672             << setw(10) << accumulate.dsAccTotS_[s][ds] << " " << right << setw(10) << fixed << setprecision(1)
0673             << scale * acceptedDS << " " << right << setw(10) << fixed << setprecision(1)
0674             << ((accumulate.hlAccept_[refIndex_] - acceptedDS > 0)
0675                     ? refRate_ *
0676                           ROOT::Math::beta_quantile(alpha, acceptedDS + 1, accumulate.hlAccept_[refIndex_] - acceptedDS)
0677                     : 0)
0678             << " " << right << setw(10) << fixed << setprecision(1) << scale * accumulate.dsAccTotS_[s][ds] << " "
0679             << datasetNames_[dsIndex_[s][ds]] << endl;
0680       }
0681     }
0682 
0683   } else {
0684     LogVerbatim("HLTrigReport") << "HLT-Report - No HLT paths found!" << endl;
0685   }
0686 
0687   LogVerbatim("HLTrigReport") << endl;
0688   LogVerbatim("HLTrigReport") << "HLT-Report end!" << endl;
0689   LogVerbatim("HLTrigReport") << endl;
0690 
0691   return;
0692 }
0693 
0694 // declare this class as a framework plugin
0695 #include "FWCore/Framework/interface/MakerMacros.h"
0696 DEFINE_FWK_MODULE(HLTrigReport);