CMS 3D CMS Logo

Public Member Functions | Private Attributes | Static Private Attributes

HLTPrescaler Class Reference

#include <HLTPrescaler.h>

Inheritance diagram for HLTPrescaler:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual bool beginLuminosityBlock (edm::LuminosityBlock &lb, edm::EventSetup const &iSetup)
virtual void endJob ()
virtual bool filter (edm::Event &iEvent, edm::EventSetup const &iSetup)
 HLTPrescaler (edm::ParameterSet const &iConfig)
virtual ~HLTPrescaler ()

Private Attributes

unsigned int acceptCount_
 accept counter
unsigned int eventCount_
 event counter
edm::InputTag gtDigi_
 GT payload, to extract the prescale column index.
bool newLumi_
 check for (re)initialization of the prescale
unsigned int offsetCount_
 initial offset
unsigned int offsetPhase_
unsigned int prescaleFactor_
 accept one in prescaleFactor_; 0 means never to accept an event
edm::service::PrescaleServiceprescaleService_
 prescale service
unsigned int prescaleSet_
 l1 prescale set index

Static Private Attributes

static const unsigned int prescaleSeed_ = 65537
 "seed" used to initialize the prescale counter

Detailed Description

This class is an EDFilter implementing an HLT Prescaler module with associated book keeping.

Date:
2012/01/23 10:15:21
Revision:
1.23
Author:
Martin Grunewald
Philipp Schieferdecker

Definition at line 24 of file HLTPrescaler.h.


Constructor & Destructor Documentation

HLTPrescaler::HLTPrescaler ( edm::ParameterSet const &  iConfig) [explicit]

Definition at line 29 of file HLTPrescaler.cc.

References LogDebug, cppFunctionSkipper::operator, and prescaleService_.

                                                         :
  prescaleSet_(0)
  , prescaleFactor_(1)
  , eventCount_(0)
  , acceptCount_(0)
  , offsetCount_(0)
  , offsetPhase_(iConfig.existsAs<unsigned int>("offset") ? iConfig.getParameter<unsigned int>("offset") : 0)
  , prescaleService_(0)
  , newLumi_(true)
  , gtDigi_ (iConfig.getParameter<edm::InputTag>("L1GtReadoutRecordTag"))
{
  if(edm::Service<edm::service::PrescaleService>().isAvailable())
    prescaleService_ = edm::Service<edm::service::PrescaleService>().operator->();
  else 
    LogDebug("NoPrescaleService")<<"PrescaleService unavailable, prescaleFactor=1!";
}
HLTPrescaler::~HLTPrescaler ( ) [virtual]

Definition at line 47 of file HLTPrescaler.cc.

{
  
}

Member Function Documentation

bool HLTPrescaler::beginLuminosityBlock ( edm::LuminosityBlock lb,
edm::EventSetup const &  iSetup 
) [virtual]

Reimplemented from edm::EDFilter.

Definition at line 58 of file HLTPrescaler.cc.

References newLumi_.

{
  newLumi_ = true;

  return true;
}
void HLTPrescaler::endJob ( void  ) [virtual]

Reimplemented from edm::EDFilter.

Definition at line 121 of file HLTPrescaler.cc.

References acceptCount_, eventCount_, and max().

{
  edm::LogInfo("PrescaleSummary")
    << acceptCount_<< "/" <<eventCount_
    << " ("
    << 100.*acceptCount_/static_cast<double>(std::max(1u,eventCount_))
    << "% of events accepted).";
  return;
}
bool HLTPrescaler::filter ( edm::Event iEvent,
edm::EventSetup const &  iSetup 
) [virtual]

Implements edm::EDFilter.

Definition at line 68 of file HLTPrescaler.cc.

References acceptCount_, edm::EDFilter::currentContext(), edm::EventID::event(), eventCount_, edm::Event::getByLabel(), edm::Event::getLuminosityBlock(), edm::service::PrescaleService::getPrescale(), gtDigi_, patZpeak::handle, edm::EventBase::id(), edm::HandleBase::isValid(), newLumi_, offsetCount_, offsetPhase_, EgammaValidation_cff::pathName, prescaleFactor_, prescaleSeed_, prescaleService_, prescaleSet_, query::result, and AlCaHLTBitMon_QueryRunRegistry::string.

{
  // during the first event of a LumiSection, read from the GT the prescale index for this
  // LumiSection and get the corresponding prescale factor from the PrescaleService
  if (newLumi_) {
    newLumi_ = false;

    bool needsInit (eventCount_==0);

    if (prescaleService_) {
      std::string const & pathName = * currentContext()->pathName();
      const unsigned int oldSet(prescaleSet_);
      const unsigned int oldPrescale(prescaleFactor_);

      edm::Handle<L1GlobalTriggerReadoutRecord> handle;
      iEvent.getByLabel(gtDigi_ , handle);
      if (handle.isValid()) {
        prescaleSet_ = handle->gtFdlWord().gtPrescaleFactorIndexAlgo();
        // gtPrescaleFactorIndexTech() is also available
        // by construction, they should always return the same index
        prescaleFactor_ = prescaleService_->getPrescale(prescaleSet_, pathName);
      } else {
        edm::LogWarning("HLT") << "Cannot read prescale column index from GT data: using default as defined by configuration or DAQ";
        prescaleFactor_ = prescaleService_->getPrescale(pathName);
      }

      if (prescaleSet_ != oldSet) {
        edm::LogInfo("ChangedPrescale")
          << "lumiBlockNb = " << iEvent.getLuminosityBlock().id().luminosityBlock()
          << ", set = " << prescaleSet_ << " [" << oldSet <<"]" 
          << ", path = "<< pathName
          << ": " << prescaleFactor_ << " [" <<oldPrescale<<"]";
        // reset the prescale counter
        needsInit = true;
      }
    }

    if (needsInit && (prescaleFactor_ != 0)) {
      // initialize the prescale counter to the first event number multiplied by a big "seed"
      offsetCount_ = ((uint64_t) (iEvent.id().event() + offsetPhase_) * prescaleSeed_) % prescaleFactor_;
    }
  }

  const bool result ( (prescaleFactor_ == 0) ? 
                      false : ((eventCount_ + offsetCount_) % prescaleFactor_ == 0) );

  ++eventCount_;
  if (result) ++acceptCount_;
  return result;
}

Member Data Documentation

unsigned int HLTPrescaler::acceptCount_ [private]

accept counter

Definition at line 58 of file HLTPrescaler.h.

Referenced by endJob(), and filter().

unsigned int HLTPrescaler::eventCount_ [private]

event counter

Definition at line 55 of file HLTPrescaler.h.

Referenced by endJob(), and filter().

GT payload, to extract the prescale column index.

Definition at line 71 of file HLTPrescaler.h.

Referenced by filter().

bool HLTPrescaler::newLumi_ [private]

check for (re)initialization of the prescale

Definition at line 68 of file HLTPrescaler.h.

Referenced by beginLuminosityBlock(), and filter().

unsigned int HLTPrescaler::offsetCount_ [private]

initial offset

Definition at line 61 of file HLTPrescaler.h.

Referenced by filter().

unsigned int HLTPrescaler::offsetPhase_ [private]

Definition at line 62 of file HLTPrescaler.h.

Referenced by filter().

unsigned int HLTPrescaler::prescaleFactor_ [private]

accept one in prescaleFactor_; 0 means never to accept an event

Definition at line 52 of file HLTPrescaler.h.

Referenced by filter().

const unsigned int HLTPrescaler::prescaleSeed_ = 65537 [static, private]

"seed" used to initialize the prescale counter

Definition at line 75 of file HLTPrescaler.h.

Referenced by filter().

prescale service

Definition at line 65 of file HLTPrescaler.h.

Referenced by filter(), and HLTPrescaler().

unsigned int HLTPrescaler::prescaleSet_ [private]

l1 prescale set index

Definition at line 49 of file HLTPrescaler.h.

Referenced by filter().