CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/HLTrigger/HLTcore/plugins/HLTPrescaler.cc

Go to the documentation of this file.
00001 
00002 //
00003 // HLTPrescaler
00004 // ------------
00005 //
00006 //           04/25/2008 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
00008 
00009 
00010 #include "HLTrigger/HLTcore/interface/HLTPrescaler.h"
00011 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00012 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00013 
00014 #include "FWCore/Framework/interface/LuminosityBlock.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 
00019 // initialize static member variables
00021 
00022 const unsigned int HLTPrescaler::prescaleSeed_ = 65537;
00023 
00025 // construction/destruction
00027 
00028 //_____________________________________________________________________________
00029 HLTPrescaler::HLTPrescaler(edm::ParameterSet const& iConfig)
00030   : prescaleFactor_(1)
00031   , eventCount_(0)
00032   , acceptCount_(0)
00033   , offsetCount_(0)
00034   , offsetPhase_(iConfig.existsAs<unsigned int>("offset") ? iConfig.getParameter<unsigned int>("offset") : 0)
00035   , prescaleService_(0)
00036   , newLumi_(true)
00037   , gtDigi_ (iConfig.getParameter<edm::InputTag>("L1GtReadoutRecordTag"))
00038 {
00039   if(edm::Service<edm::service::PrescaleService>().isAvailable())
00040     prescaleService_ = edm::Service<edm::service::PrescaleService>().operator->();
00041   else 
00042     LogDebug("NoPrescaleService")<<"PrescaleService unavailable, prescaleFactor=1!";
00043 }
00044 
00045 //_____________________________________________________________________________    
00046 HLTPrescaler::~HLTPrescaler()
00047 {
00048   
00049 }
00050 
00051 
00053 // implementation of member functions
00055 
00056 //______________________________________________________________________________
00057 bool HLTPrescaler::beginLuminosityBlock(edm::LuminosityBlock & lb,
00058                                         edm::EventSetup const& iSetup)
00059 {
00060   newLumi_ = true;
00061 
00062   return true;
00063 }
00064 
00065 
00066 //_____________________________________________________________________________
00067 bool HLTPrescaler::filter(edm::Event& iEvent, const edm::EventSetup&)
00068 {
00069   // during the first event of a LumiSection, read from the GT the prescale index for this
00070   // LumiSection and get the corresponding prescale factor from the PrescaleService
00071   if (newLumi_) {
00072     newLumi_ = false;
00073 
00074     bool needsInit (eventCount_==0);
00075 
00076     if (prescaleService_) {
00077       const unsigned int oldPrescale(prescaleFactor_);
00078 
00079       edm::Handle<L1GlobalTriggerReadoutRecord> handle;
00080       iEvent.getByLabel(gtDigi_ , handle);
00081       if (handle.isValid()) {
00082         unsigned int index = handle->gtFdlWord().gtPrescaleFactorIndexAlgo();
00083         // gtPrescaleFactorIndexTech() is also available
00084         // by construction, they should always return the same index
00085         prescaleFactor_ = prescaleService_->getPrescale(index, *pathName());
00086       } else {
00087         edm::LogWarning("HLT") << "Cannot read prescale column index from GT data: using default as defined by configuration or DAQ";
00088         prescaleFactor_ = prescaleService_->getPrescale(*pathName());
00089       }
00090 
00091       if (prescaleFactor_ != oldPrescale) {
00092         edm::LogInfo("ChangedPrescale")
00093           << "lumiBlockNb="<< iEvent.getLuminosityBlock().id().luminosityBlock() << ", "
00094           << "path="<<*pathName()<<": "
00095           << prescaleFactor_ << " [" <<oldPrescale<<"]";
00096         // reset the prescale counter
00097         needsInit = true;
00098       }
00099     }
00100 
00101     if (needsInit && (prescaleFactor_ != 0)) {
00102       // initialize the prescale counter to the first event number multiplied by a big "seed"
00103       offsetCount_ = ((uint64_t) (iEvent.id().event() + offsetPhase_) * prescaleSeed_) % prescaleFactor_;
00104     }
00105   }
00106 
00107   const bool result ( (prescaleFactor_ == 0) ? 
00108                       false : ((eventCount_ + offsetCount_) % prescaleFactor_ == 0) );
00109 
00110   ++eventCount_;
00111   if (result) ++acceptCount_;
00112   return result;
00113 }
00114 
00115 
00116 //_____________________________________________________________________________
00117 void HLTPrescaler::endJob()
00118 {
00119   edm::LogInfo("PrescaleSummary")
00120     << acceptCount_<< "/" <<eventCount_
00121     << " ("
00122     << 100.*acceptCount_/static_cast<double>(std::max(1u,eventCount_))
00123     << "% of events accepted).";
00124   return;
00125 }