CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes

HLTrigReport Class Reference

#include <HLTrigReport.h>

Inheritance diagram for HLTrigReport:
edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void beginLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
virtual void beginRun (edm::Run const &, edm::EventSetup const &)
const std::vector< unsigned int > & datasetCounts () const
const std::vector< std::string > & datasetNames () const
virtual void endJob ()
virtual void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
virtual void endRun (edm::Run const &, edm::EventSetup const &)
 HLTrigReport (const edm::ParameterSet &)
void reset (bool changed=false)
const std::vector< unsigned int > & streamCounts () const
const std::vector< std::string > & streamNames () const
 ~HLTrigReport ()

Static Public Member Functions

static ReportEvery decode (const std::string &value)

Private Types

enum  ReportEvery {
  NEVER = 0, EVERY_EVENT = 1, EVERY_LUMI = 2, EVERY_RUN = 3,
  EVERY_JOB = 4
}

Private Member Functions

void dumpReport (std::string const &header=std::string())

Private Attributes

bool configured_
std::vector< std::vector
< std::string > > 
datasetContents_
std::vector< std::string > datasetNames_
std::vector< std::vector
< unsigned int > > 
dsAccTotS_
std::vector< unsigned int > dsAllTotS_
std::vector< std::vector
< unsigned int > > 
dsIndex_
std::vector< unsigned int > hlAccept_
std::vector< unsigned int > hlAccTot_
std::vector< std::vector
< unsigned int > > 
hlAccTotDS_
std::vector< unsigned int > hlAllTotDS_
std::vector< unsigned int > hlErrors_
std::vector< std::vector
< unsigned int > > 
hlIndex_
std::vector< std::string > hlNames_
HLTConfigProvider hltConfig_
std::vector< unsigned int > hltL1s_
std::vector< unsigned int > hltPre_
edm::InputTag hlTriggerResults_
std::vector< unsigned int > hlWasRun_
bool isCustomDatasets_
bool isCustomStreams_
unsigned int nAccept_
unsigned int nErrors_
unsigned int nEvents_
unsigned int nWasRun_
std::vector< int > posL1s_
std::vector< int > posPre_
unsigned int refIndex_
std::string refPath_
double refRate_
ReportEvery reportBy_
ReportEvery resetBy_
ReportEvery serviceBy_
std::vector< std::vector
< std::string > > 
streamContents_
std::vector< std::string > streamNames_

Detailed Description

This class is an EDAnalyzer implementing TrigReport (statistics printed to log file) for HL triggers

Date:
2011/03/30 16:17:33
Revision:
1.17
Author:
Martin Grunewald

See header file for documentation

Date:
2011/03/30 16:17:33
Revision:
1.33
Author:
Martin Grunewald

Definition at line 29 of file HLTrigReport.h.


Member Enumeration Documentation

enum HLTrigReport::ReportEvery [private]
Enumerator:
NEVER 
EVERY_EVENT 
EVERY_LUMI 
EVERY_RUN 
EVERY_JOB 

Definition at line 31 of file HLTrigReport.h.

                       {
        NEVER       = 0,
        EVERY_EVENT = 1,
        EVERY_LUMI  = 2,
        EVERY_RUN   = 3,
        EVERY_JOB   = 4
      };

Constructor & Destructor Documentation

HLTrigReport::HLTrigReport ( const edm::ParameterSet iConfig) [explicit]

Definition at line 53 of file HLTrigReport.cc.

References datasetContents_, datasetNames_, edm::InputTag::encode(), edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNamesForType(), edm::ParameterSet::getUntrackedParameter(), hlTriggerResults_, isCustomDatasets_, isCustomStreams_, LogDebug, AlCaRecoCosmics_cfg::name, NEVER, refIndex_, refPath_, refRate_, serviceBy_, streamContents_, and streamNames_.

                                                         :
  hlTriggerResults_ (iConfig.getParameter<edm::InputTag> ("HLTriggerResults")),
  configured_(false),
  nEvents_(0),
  nWasRun_(0),
  nAccept_(0),
  nErrors_(0),
  hlWasRun_(0),
  hltL1s_(0),
  hltPre_(0),
  hlAccept_(0),
  hlAccTot_(0),
  hlErrors_(0),
  posL1s_(0),
  posPre_(0),
  hlNames_(0),
  hlIndex_(0),
  hlAccTotDS_(0),
  datasetNames_(0),
  datasetContents_(0),
  isCustomDatasets_(false),
  dsIndex_(0),
  dsAccTotS_(0),
  streamNames_(0),
  streamContents_(0),
  isCustomStreams_(false),
  refPath_("HLTriggerFinalPath"),
  refIndex_(0),
  refRate_(100.0),
  reportBy_( decode(iConfig.getUntrackedParameter<std::string>("reportBy",  "job")) ),
  resetBy_(  decode(iConfig.getUntrackedParameter<std::string>("resetBy",   "never")) ),
  serviceBy_(decode(iConfig.getUntrackedParameter<std::string>("serviceBy", "never")) ),
  hltConfig_()
{
  const edm::ParameterSet customDatasets(iConfig.getUntrackedParameter<edm::ParameterSet>("CustomDatasets", edm::ParameterSet()));
  isCustomDatasets_ = (customDatasets != edm::ParameterSet());
  if (isCustomDatasets_) {
    datasetNames_ = customDatasets.getParameterNamesForType<std::vector<std::string> >();
    for (std::vector<std::string>::const_iterator name = datasetNames_.begin(); name != datasetNames_.end(); name++) {
      datasetContents_.push_back(customDatasets.getParameter<std::vector<std::string> >(*name));
    }
  }

  const edm::ParameterSet customStreams (iConfig.getUntrackedParameter<edm::ParameterSet>("CustomStreams" , edm::ParameterSet()));
  isCustomStreams_  = (customStreams  != edm::ParameterSet());
  if (isCustomStreams_ ) {
    streamNames_ = customStreams.getParameterNamesForType<std::vector<std::string> >();
    for (std::vector<std::string>::const_iterator name = streamNames_.begin(); name != streamNames_.end(); name++) {
      streamContents_.push_back(customStreams.getParameter<std::vector<std::string> >(*name));
    }
  }

  refPath_ = iConfig.getUntrackedParameter<std::string>("ReferencePath","HLTriggerFinalPath");
  refRate_ = iConfig.getUntrackedParameter<double>("ReferenceRate", 100.0);
  refIndex_= 0;

  LogDebug("HLTrigReport")
    << "HL TiggerResults: " + hlTriggerResults_.encode()
    << " using reference path and rate: " + refPath_ + " " << refRate_;

  if (serviceBy_ != NEVER and edm::Service<HLTrigReportService>()) {
    edm::Service<HLTrigReportService>()->registerModule(this);
  }

}
HLTrigReport::~HLTrigReport ( )

Definition at line 119 of file HLTrigReport.cc.

{ }

Member Function Documentation

void HLTrigReport::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 371 of file HLTrigReport.cc.

References accept(), configured_, datasetCounts(), dsAccTotS_, dsAllTotS_, dsIndex_, edm::InputTag::encode(), EVERY_EVENT, edm::Event::getByLabel(), hlAccept_, hlAccTot_, hlAccTotDS_, hlAllTotDS_, hlErrors_, hlIndex_, hlNames_, hltL1s_, hltPre_, hlTriggerResults_, hlWasRun_, i, edm::EventBase::id(), getHLTprescales::index, LogDebug, edm::EventBase::luminosityBlock(), n, nAccept_, nErrors_, nEvents_, nWasRun_, L1TEmulatorMonitor_cff::p, posL1s_, posPre_, reportBy_, reset(), resetBy_, edm::Event::run(), asciidump::s, serviceBy_, and streamCounts().

{
  // accumulation of statistics event by event

  using namespace std;
  using namespace edm;

  if (resetBy_ == EVERY_EVENT) reset();

  nEvents_++;

  // get hold of TriggerResults
  Handle<TriggerResults> HLTR;
  iEvent.getByLabel(hlTriggerResults_, HLTR);
  if (HLTR.isValid()) {
    if (HLTR->wasrun()) nWasRun_++;
    const bool accept(HLTR->accept());
    LogDebug("HLTrigReport") << "HLT TriggerResults decision: " << accept;
    if (accept) ++nAccept_;
    if (HLTR->error()) nErrors_++;
  } else {
    LogDebug("HLTrigReport") << "HLT TriggerResults with label ["+hlTriggerResults_.encode()+"] not found!";
    nErrors_++;
    return;
  }

  // HLTConfigProvider not configured - cannot produce any detailed statistics
  if (not configured_)
    return;

  // decision for each HL algorithm
  const unsigned int n(hlNames_.size());
  bool acceptedByPrevoiusPaths = false;
  for (unsigned int i=0; i!=n; ++i) {
    if (HLTR->wasrun(i)) hlWasRun_[i]++;
    if (HLTR->accept(i)) {
      acceptedByPrevoiusPaths = true;
      hlAccept_[i]++;
    }
    if (acceptedByPrevoiusPaths) hlAccTot_[i]++;
    if (HLTR->error(i) ) hlErrors_[i]++;
    const int index(static_cast<int>(HLTR->index(i)));
    if (HLTR->accept(i)) {
      if (index >= posL1s_[i]) hltL1s_[i]++;
      if (index >= posPre_[i]) hltPre_[i]++;
    } else {
      if (index >  posL1s_[i]) hltL1s_[i]++;
      if (index >  posPre_[i]) hltPre_[i]++;
    }
  }

  // calculate accumulation of accepted events by a path within a dataset
  std::vector<bool> acceptedByDS(hlIndex_.size(), false);
  for (size_t ds=0; ds<hlIndex_.size(); ++ds) {
    for (size_t p=0; p<hlIndex_[ds].size(); ++p) {
      if (acceptedByDS[ds] or HLTR->accept(hlIndex_[ds][p])) {
        acceptedByDS[ds] = true;
        hlAccTotDS_[ds][p]++;
      }
    }
    if (acceptedByDS[ds]) hlAllTotDS_[ds]++;
  }

  // calculate accumulation of accepted events by a dataset within a stream
  for (size_t s=0; s<dsIndex_.size(); ++s) {
    bool acceptedByS = false;
    for (size_t ds=0; ds<dsIndex_[s].size(); ++ds) {
      if (acceptedByS or acceptedByDS[dsIndex_[s][ds]]) {
        acceptedByS = true;
        dsAccTotS_[s][ds]++;
      }
    }
    if (acceptedByS) dsAllTotS_[s]++;
  }

  if (reportBy_ == EVERY_EVENT) {
    std::stringstream stream;
    stream << "Summary for Run " << iEvent.run() << ", LumiSection " << iEvent.luminosityBlock() << ", Event " << iEvent.id();
  }
  if (serviceBy_ == EVERY_EVENT and edm::Service<HLTrigReportService>()) {
    edm::Service<HLTrigReportService>()->setDatasetCounts(datasetCounts());
    edm::Service<HLTrigReportService>()->setStreamCounts(streamCounts());
  }

}
void HLTrigReport::beginJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 286 of file HLTrigReport.cc.

References EVERY_JOB, reset(), and resetBy_.

                            {
  if (resetBy_ == EVERY_JOB)
    reset();
}
void HLTrigReport::beginLuminosityBlock ( edm::LuminosityBlock const &  lumi,
edm::EventSetup const &  setup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 352 of file HLTrigReport.cc.

References EVERY_RUN, reset(), and resetBy_.

                                                                                                    {
  if (resetBy_ == EVERY_RUN)  reset();
}
void HLTrigReport::beginRun ( edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 302 of file HLTrigReport.cc.

References configured_, dsAccTotS_, dsAllTotS_, dsIndex_, dumpReport(), EVERY_RUN, hlAccept_, hlAccTot_, hlAccTotDS_, hlAllTotDS_, hlErrors_, hlIndex_, hlNames_, hltConfig_, hltL1s_, hltPre_, hlTriggerResults_, hlWasRun_, HLTConfigProvider::init(), nAccept_, nErrors_, nEvents_, nWasRun_, posL1s_, posPre_, edm::InputTag::process(), reset(), and resetBy_.

{
  bool changed = true;
  if (hltConfig_.init(iRun, iSetup, hlTriggerResults_.process(), changed)) {
    configured_ = true;
    if (changed) {
      dumpReport("Summary for this HLT table");
      reset(true);
    }
  } else {
    dumpReport("Summary for this HLT table");
    // cannot initialize the HLT menu - reset and clear all counters and tables
    configured_ = false;
    nEvents_    = 0;
    nWasRun_    = 0;
    nAccept_    = 0;
    nErrors_    = 0;
    hlWasRun_.clear();
    hltL1s_.clear();
    hltPre_.clear();
    hlAccept_.clear();
    hlAccTot_.clear();
    hlErrors_.clear();
    posL1s_.clear();
    posPre_.clear();
    hlNames_.clear();
    hlIndex_.clear();
    hlAccTotDS_.clear();
    hlAllTotDS_.clear();
    dsIndex_.clear();
    dsAccTotS_.clear();
    dsAllTotS_.clear();
  }

  if (resetBy_ == EVERY_RUN) reset();

}
const std::vector< unsigned int > & HLTrigReport::datasetCounts ( ) const

Definition at line 131 of file HLTrigReport.cc.

References hlAllTotDS_.

Referenced by analyze(), endJob(), endLuminosityBlock(), and endRun().

                                                                 {
  return hlAllTotDS_;
}
const std::vector< std::string > & HLTrigReport::datasetNames ( ) const

Definition at line 125 of file HLTrigReport.cc.

References datasetNames_.

                                                             {
  return datasetNames_;
}
HLTrigReport::ReportEvery HLTrigReport::decode ( const std::string &  value) [static]

Definition at line 30 of file HLTrigReport.cc.

References EVERY_EVENT, EVERY_JOB, EVERY_LUMI, EVERY_RUN, Exception, and NEVER.

                                                                    {
  if (value == "never")
    return NEVER;

  if (value == "job")
    return EVERY_JOB;

  if (value == "run")
    return EVERY_RUN;

  if (value == "lumi")
    return EVERY_LUMI;

  if (value == "event")
    return EVERY_EVENT;
  
  throw cms::Exception("Configuration") << "Invalid option value \"" << value << "\". Legal values are \"job\", \"run\", \"lumi\", \"event\" and \"never\".";
}
void HLTrigReport::dumpReport ( std::string const &  header = std::string()) [private]

Definition at line 458 of file HLTrigReport.cc.

References alpha, configured_, datasetNames_, dsAccTotS_, dsAllTotS_, dsIndex_, hlAccept_, hlAccTot_, hlAccTotDS_, hlAllTotDS_, hlErrors_, hlIndex_, hlNames_, hltL1s_, hltPre_, hlWasRun_, i, max(), n, nAccept_, nErrors_, nEvents_, nWasRun_, L1TEmulatorMonitor_cff::p, refIndex_, refRate_, asciidump::s, and streamNames_.

Referenced by beginRun(), endJob(), endLuminosityBlock(), and endRun().

{
  // final printout of accumulated statistics

  using namespace std;
  using namespace edm;
  const unsigned int n(hlNames_.size());

  if ((n==0) and (nEvents_==0)) return;

  LogVerbatim("HLTrigReport") << dec << endl;
  LogVerbatim("HLTrigReport") << "HLT-Report " << "---------- Event  Summary ------------" << endl;
  if (not header.empty())
    LogVerbatim("HLTrigReport") << "HLT-Report " << header << endl;
  LogVerbatim("HLTrigReport") << "HLT-Report"
         << " Events total = " << nEvents_
         << " wasrun = " << nWasRun_
         << " passed = " << nAccept_
         << " errors = " << nErrors_
         << endl;

  // HLTConfigProvider not configured - cannot produce any detailed statistics
  if (not configured_)
    return;

  double scale = hlAccept_[refIndex_]>0 ? refRate_/hlAccept_[refIndex_] : 0.;
  double alpha = 1 - (1.0 - .6854)/2; // for the Clopper-Pearson 68% CI

  LogVerbatim("HLTrigReport") << endl;
  LogVerbatim("HLTrigReport") << "HLT-Report " << "---------- HLTrig Summary ------------" << endl;
  LogVerbatim("HLTrigReport") << "HLT-Report "
         << right << setw(7) << "HLT #" << " "
         << right << setw(7) << "WasRun" << " "
         << right << setw(7) << "L1S" << " "
         << right << setw(7) << "Pre" << " "
         << right << setw(7) << "HLT" << " "
         << right << setw(9) << "%L1sPre" << " "
         << right << setw(7) << "Rate" << " "
         << right << setw(7) << "RateHi" << " "
         << right << setw(7) << "HLTtot" << " "
         << right << setw(7) << "RateTot" << " "
         << right << setw(7) << "Errors" << " "
         << "Name" << endl;

  if (n>0) {
    for (unsigned int i=0; i!=n; ++i) {
      LogVerbatim("HLTrigReport") << "HLT-Report "
           << right << setw(7) << i << " "
           << right << setw(7) << hlWasRun_[i] << " "
           << right << setw(7) << hltL1s_[i] << " "
           << right << setw(7) << hltPre_[i] << " "
           << right << setw(7) << hlAccept_[i] << " "
           << right << setw(9) << fixed << setprecision(5)
           << static_cast<float>(100*hlAccept_[i])/
              static_cast<float>(max(hltPre_[i], 1u)) << " "
           << right << setw(7) << fixed << setprecision(1) << scale*hlAccept_[i] << " "
           << right << setw(7) << fixed << setprecision(1) <<
              ((hlAccept_[refIndex_]-hlAccept_[i] > 0) ? refRate_*ROOT::Math::beta_quantile(alpha, hlAccept_[i]+1, hlAccept_[refIndex_]-hlAccept_[i]) : 0) << " "
           << right << setw(7) << hlAccTot_[i] << " "
           << right << setw(7) << fixed << setprecision(1) << scale*hlAccTot_[i] << " "
           << right << setw(7) << hlErrors_[i] << " "
           << hlNames_[i] << endl;
    }

    // now for each dataset
    for (size_t ds=0; ds<hlIndex_.size(); ++ds) {
      LogVerbatim("HLTrigReport") << endl;
      LogVerbatim("HLTrigReport") << "HLT-Report " << "---------- Dataset Summary: " << datasetNames_[ds] << " ------------" << hlAllTotDS_[ds] << endl;
      LogVerbatim("HLTrigReport") << "HLT-Report "
         << right << setw(7) << "HLT #" << " "
         << right << setw(7) << "WasRun" << " "
         << right << setw(7) << "L1S" << " "
         << right << setw(7) << "Pre" << " "
         << right << setw(7) << "HLT" << " "
         << right << setw(9) << "%L1sPre" << " "
         << right << setw(7) << "Rate" << " "
         << right << setw(7) << "RateHi" << " "
         << right << setw(7) << "HLTtot" << " "
         << right << setw(7) << "RateTot" << " "
         << right << setw(7) << "Errors" << " "
         << "Name" << endl;
      for (size_t p=0; p<hlIndex_[ds].size(); ++p) {
        LogVerbatim("HLTrigReport") << "HLT-Report "
           << right << setw(7) << p << " "
           << right << setw(7) << hlWasRun_[hlIndex_[ds][p]] << " "
           << right << setw(7) << hltL1s_[hlIndex_[ds][p]] << " "
           << right << setw(7) << hltPre_[hlIndex_[ds][p]] << " "
           << right << setw(7) << hlAccept_[hlIndex_[ds][p]] << " "
           << right << setw(9) << fixed << setprecision(5)
           << static_cast<float>(100*hlAccept_[hlIndex_[ds][p]])/
              static_cast<float>(max(hltPre_[hlIndex_[ds][p]], 1u)) << " "
           << right << setw(7) << fixed << setprecision(1) << scale*hlAccept_[hlIndex_[ds][p]] << " "
           << right << setw(7) << fixed << setprecision(1) <<
              ((hlAccept_[refIndex_]-hlAccept_[hlIndex_[ds][p]] > 0) ? refRate_*ROOT::Math::beta_quantile(alpha, hlAccept_[hlIndex_[ds][p]]+1, hlAccept_[refIndex_]-hlAccept_[hlIndex_[ds][p]]) : 0) << " "
           << right << setw(7) << hlAccTotDS_[ds][p] << " "
           << right << setw(7) << fixed << setprecision(1) << scale*hlAccTotDS_[ds][p] << " "
           << right << setw(7) << hlErrors_[hlIndex_[ds][p]] << " "
           << hlNames_[hlIndex_[ds][p]] << endl;
      }
    }

    // now for each stream
    for (size_t s=0; s<dsIndex_.size(); ++s) {
      LogVerbatim("HLTrigReport") << endl;
      LogVerbatim("HLTrigReport") << "HLT-Report " << "---------- Stream Summary: " << streamNames_[s] << " ------------" << dsAllTotS_[s] << endl;
      LogVerbatim("HLTrigReport") << "HLT-Report "
         << right << setw(10) << "Dataset #" << " "
         << right << setw(10) << "Individual" << " "
         << right << setw(10) << "Total" << " "
         << right << setw(10) << "Rate" << " "
         << right << setw(10) << "RateHi" << " "
         << right << setw(10) << "RateTot" << " "
         << "Name" << endl;
      for (size_t ds=0;ds<dsIndex_[s].size(); ++ds) {
        unsigned int acceptedDS = hlAccTotDS_[dsIndex_[s][ds]][hlIndex_[dsIndex_[s][ds]].size()-1];
        LogVerbatim("HLTrigReport") << "HLT-Report "
           << right << setw(10) << ds << " "
           << right << setw(10) << acceptedDS << " "
           << right << setw(10) << dsAccTotS_[s][ds] << " "
           << right << setw(10) << fixed << setprecision(1) << scale*acceptedDS << " "
           << right << setw(10) << fixed << setprecision(1) <<
              ((hlAccept_[refIndex_]-acceptedDS > 0) ? refRate_*ROOT::Math::beta_quantile(alpha, acceptedDS+1, hlAccept_[refIndex_]-acceptedDS) : 0) << " "
           << right << setw(10) << fixed << setprecision(1) << scale*dsAccTotS_[s][ds] << " "
           << datasetNames_[dsIndex_[s][ds]] << endl;
      }
    }

  } else {
    LogVerbatim("HLTrigReport") << "HLT-Report - No HLT paths found!" << endl;
  }

  LogVerbatim("HLTrigReport") << endl;
  LogVerbatim("HLTrigReport") << "HLT-Report end!" << endl;
  LogVerbatim("HLTrigReport") << endl;

  return;
}
void HLTrigReport::endJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 291 of file HLTrigReport.cc.

References datasetCounts(), dumpReport(), EVERY_JOB, reportBy_, serviceBy_, and streamCounts().

void HLTrigReport::endLuminosityBlock ( edm::LuminosityBlock const &  lumi,
edm::EventSetup const &  setup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 356 of file HLTrigReport.cc.

References datasetCounts(), dumpReport(), EVERY_LUMI, edm::LuminosityBlockBase::luminosityBlock(), reportBy_, edm::LuminosityBlockBase::run(), serviceBy_, and streamCounts().

                                                                                                  {
  if (reportBy_ == EVERY_LUMI) {
    std::stringstream stream;
    stream << "Summary for Run " << lumi.run() << ", LumiSection " << lumi.luminosityBlock();
    dumpReport(stream.str());
  }
  if (serviceBy_ == EVERY_LUMI and edm::Service<HLTrigReportService>()) {
    edm::Service<HLTrigReportService>()->setDatasetCounts(datasetCounts());
    edm::Service<HLTrigReportService>()->setStreamCounts(streamCounts());
  }
}
void HLTrigReport::endRun ( edm::Run const &  run,
edm::EventSetup const &  setup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 340 of file HLTrigReport.cc.

References datasetCounts(), dumpReport(), EVERY_RUN, reportBy_, edm::RunBase::run(), serviceBy_, and streamCounts().

                                                                         {
  if (reportBy_ == EVERY_RUN) {
    std::stringstream stream;
    stream << "Summary for Run " << run.run();
    dumpReport(stream.str());
  }
  if (serviceBy_ == EVERY_RUN and edm::Service<HLTrigReportService>()) {
    edm::Service<HLTrigReportService>()->setDatasetCounts(datasetCounts());
    edm::Service<HLTrigReportService>()->setStreamCounts(streamCounts());
  }
}
void HLTrigReport::reset ( bool  changed = false)

Definition at line 138 of file HLTrigReport.cc.

References HLTConfigProvider::datasetContents(), datasetContents_, HLTConfigProvider::datasetNames(), datasetNames_, dsAccTotS_, dsAllTotS_, dsIndex_, hlAccept_, hlAccTot_, hlAccTotDS_, hlAllTotDS_, hlErrors_, hlIndex_, hlNames_, hltConfig_, hltL1s_, hltPre_, hlWasRun_, i, isCustomDatasets_, isCustomStreams_, j, label, HLTConfigProvider::moduleLabels(), HLTConfigProvider::moduleType(), n, nAccept_, nErrors_, nEvents_, NEVER, nWasRun_, L1TEmulatorMonitor_cff::p, posL1s_, posPre_, refIndex_, refPath_, asciidump::s, serviceBy_, findQualityFiles::size, HLTConfigProvider::streamContents(), streamContents_, HLTConfigProvider::streamNames(), streamNames_, HLTConfigProvider::triggerIndex(), and HLTConfigProvider::triggerNames().

Referenced by analyze(), beginJob(), beginLuminosityBlock(), and beginRun().

                                       {

  // reset global counters
  nEvents_ = 0;
  nWasRun_ = 0;
  nAccept_ = 0;
  nErrors_ = 0;

  // update trigger names
  if (changed)
    hlNames_ = hltConfig_.triggerNames();

  const unsigned int n = hlNames_.size();

  if (changed) {
    // resize per-path counters
    hlWasRun_.resize(n);
    hltL1s_.resize(n);
    hltPre_.resize(n);
    hlAccept_.resize(n);
    hlAccTot_.resize(n);
    hlErrors_.resize(n);
    // find the positions of seeding and prescaler modules
    posL1s_.resize(n);
    posPre_.resize(n);
    for (unsigned int i = 0; i < n; ++i) {
      posL1s_[i] = -1;
      posPre_[i] = -1;
      const std::vector<std::string> & moduleLabels(hltConfig_.moduleLabels(i));
      for (unsigned int j = 0; j < moduleLabels.size(); ++j) {
        const std::string & label = hltConfig_.moduleType(moduleLabels[j]);
        if (label == "HLTLevel1GTSeed")
          posL1s_[i] = j;
        else if (label == "HLTPrescaler")
          posPre_[i] = j;
      }
    }
  }

  // reset per-path counters
  for (unsigned int i = 0; i < n; ++i) {
    hlWasRun_[i] = 0;
    hltL1s_[i]   = 0;
    hltPre_[i]   = 0;
    hlAccept_[i] = 0;
    hlAccTot_[i] = 0;
    hlErrors_[i] = 0;
  }

  // if not overridden, reload the datasets and streams
  if (changed and not isCustomDatasets_) {
    datasetNames_    = hltConfig_.datasetNames();
    datasetContents_ = hltConfig_.datasetContents();
  }
  if (changed and not isCustomStreams_) {
    streamNames_     = hltConfig_.streamNames();
    streamContents_  = hltConfig_.streamContents();
  }

  if (changed) {
    // fill the matrices of hlIndex_, hlAccTotDS_
    hlIndex_.clear();
    hlIndex_.resize(datasetNames_.size());
    hlAccTotDS_.clear();
    hlAllTotDS_.clear();
    hlAccTotDS_.resize(datasetNames_.size());
    hlAllTotDS_.resize(datasetNames_.size());
    for (unsigned int ds = 0; ds < datasetNames_.size(); ds++) {
      unsigned int size = datasetContents_[ds].size();
      hlIndex_[ds].reserve(size);
      hlAccTotDS_[ds].reserve(size);
      hlAllTotDS_[ds]=0;
      for (unsigned int p = 0; p < size; ++p) {
        unsigned int i = hltConfig_.triggerIndex(datasetContents_[ds][p]);
        if (i<n) {
          hlIndex_[ds].push_back(i);
          hlAccTotDS_[ds].push_back(0);
        }
      }
    }
  } else {
    // reset the matrix of hlAccTotDS_
    for (unsigned int ds = 0; ds < datasetNames_.size(); ds++) {
      hlAllTotDS_[ds]=0;
      for (unsigned int i = 0; i < hlAccTotDS_[ds].size(); ++i)
          hlAccTotDS_[ds][i] = 0;
    }
  }

  if (changed) {
    // fill the matrices of dsIndex_, dsAccTotS_
    dsIndex_.clear();
    dsIndex_.resize(streamNames_.size());
    dsAccTotS_.clear();
    dsAllTotS_.clear();
    dsAccTotS_.resize(streamNames_.size());
    dsAllTotS_.resize(streamNames_.size());
    for (unsigned int s = 0; s < streamNames_.size(); ++s) {
      unsigned int size = streamContents_[s].size();
      dsIndex_.reserve(size);
      dsAccTotS_.reserve(size);
      dsAllTotS_[s]=0;
      for (unsigned int ds = 0; ds < size; ++ds) {
        unsigned int i = 0;
        for (; i<datasetNames_.size(); i++) if (datasetNames_[i] == streamContents_[s][ds]) 
          break;
        // report only datasets that have at least one path otherwise crash
        if (i < datasetNames_.size() and hlIndex_[i].size() > 0) {
          dsIndex_[s].push_back(i);
          dsAccTotS_[s].push_back(0);
        }
      }
    }
  } else {
    // reset the matrix of dsAccTotS_
    for (unsigned int s = 0; s < streamNames_.size(); ++s) {
      dsAllTotS_[s]=0;
      for (unsigned int i = 0; i < dsAccTotS_[s].size(); ++i)
        dsAccTotS_[s][i] = 0;
    }
  }

  // if needed, update the reference path
  if (changed) {
    refIndex_ = hltConfig_.triggerIndex(refPath_);
    if (refIndex_ >= n) {
      refIndex_ = 0;
      edm::LogWarning("HLTrigReport")
        << "Requested reference path '"+refPath_+"' not in HLT menu. "
        << "Using HLTriggerFinalPath instead.";
      refPath_ = "HLTriggerFinalPath";
      refIndex_ = hltConfig_.triggerIndex(refPath_);
      if (refIndex_ >= n) {
        refIndex_ = 0;
        edm::LogWarning("HLTrigReport")
          << "Requested reference path '"+refPath_+"' not in HLT menu. "
          << "Using first path in table (index=0) instead.";
      }
    }
  }

  if (changed and serviceBy_ != NEVER and edm::Service<HLTrigReportService>()) {
    edm::Service<HLTrigReportService>()->setDatasetNames(datasetNames_);
    edm::Service<HLTrigReportService>()->setStreamNames(streamNames_);
  }

}
const std::vector< unsigned int > & HLTrigReport::streamCounts ( ) const

Definition at line 134 of file HLTrigReport.cc.

References dsAllTotS_.

Referenced by analyze(), endJob(), endLuminosityBlock(), and endRun().

                                                                {
  return dsAllTotS_;
}
const std::vector< std::string > & HLTrigReport::streamNames ( ) const

Definition at line 128 of file HLTrigReport.cc.

References streamNames_.

                                                            {
  return streamNames_;
}

Member Data Documentation

bool HLTrigReport::configured_ [private]

Definition at line 69 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), and dumpReport().

std::vector<std::vector<std::string> > HLTrigReport::datasetContents_ [private]

Definition at line 91 of file HLTrigReport.h.

Referenced by HLTrigReport(), and reset().

std::vector<std::string> HLTrigReport::datasetNames_ [private]

Definition at line 90 of file HLTrigReport.h.

Referenced by datasetNames(), dumpReport(), HLTrigReport(), and reset().

std::vector<std::vector<unsigned int> > HLTrigReport::dsAccTotS_ [private]

Definition at line 94 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::dsAllTotS_ [private]

Definition at line 95 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), reset(), and streamCounts().

std::vector<std::vector<unsigned int> > HLTrigReport::dsIndex_ [private]

Definition at line 93 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::hlAccept_ [private]

Definition at line 79 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::hlAccTot_ [private]

Definition at line 80 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<std::vector<unsigned int> > HLTrigReport::hlAccTotDS_ [private]

Definition at line 88 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::hlAllTotDS_ [private]

Definition at line 89 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), datasetCounts(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::hlErrors_ [private]

Definition at line 81 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<std::vector<unsigned int> > HLTrigReport::hlIndex_ [private]

Definition at line 87 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<std::string> HLTrigReport::hlNames_ [private]

Definition at line 85 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

Definition at line 106 of file HLTrigReport.h.

Referenced by beginRun(), and reset().

std::vector<unsigned int> HLTrigReport::hltL1s_ [private]

Definition at line 77 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<unsigned int> HLTrigReport::hltPre_ [private]

Definition at line 78 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

Definition at line 68 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), and HLTrigReport().

std::vector<unsigned int> HLTrigReport::hlWasRun_ [private]

Definition at line 76 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

Definition at line 92 of file HLTrigReport.h.

Referenced by HLTrigReport(), and reset().

Definition at line 98 of file HLTrigReport.h.

Referenced by HLTrigReport(), and reset().

unsigned int HLTrigReport::nAccept_ [private]

Definition at line 73 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

unsigned int HLTrigReport::nErrors_ [private]

Definition at line 74 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

unsigned int HLTrigReport::nEvents_ [private]

Definition at line 71 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

unsigned int HLTrigReport::nWasRun_ [private]

Definition at line 72 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), dumpReport(), and reset().

std::vector<int> HLTrigReport::posL1s_ [private]

Definition at line 83 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), and reset().

std::vector<int> HLTrigReport::posPre_ [private]

Definition at line 84 of file HLTrigReport.h.

Referenced by analyze(), beginRun(), and reset().

unsigned int HLTrigReport::refIndex_ [private]

Definition at line 100 of file HLTrigReport.h.

Referenced by dumpReport(), HLTrigReport(), and reset().

std::string HLTrigReport::refPath_ [private]

Definition at line 99 of file HLTrigReport.h.

Referenced by HLTrigReport(), and reset().

double HLTrigReport::refRate_ [private]

Definition at line 101 of file HLTrigReport.h.

Referenced by dumpReport(), and HLTrigReport().

Definition at line 103 of file HLTrigReport.h.

Referenced by analyze(), endJob(), endLuminosityBlock(), and endRun().

Definition at line 104 of file HLTrigReport.h.

Referenced by analyze(), beginJob(), beginLuminosityBlock(), and beginRun().

Definition at line 105 of file HLTrigReport.h.

Referenced by analyze(), endJob(), endLuminosityBlock(), endRun(), HLTrigReport(), and reset().

std::vector<std::vector<std::string> > HLTrigReport::streamContents_ [private]

Definition at line 97 of file HLTrigReport.h.

Referenced by HLTrigReport(), and reset().

std::vector<std::string> HLTrigReport::streamNames_ [private]

Definition at line 96 of file HLTrigReport.h.

Referenced by dumpReport(), HLTrigReport(), reset(), and streamNames().