CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPrescaler.cc
Go to the documentation of this file.
1 //
3 // HLTPrescaler
4 // ------------
5 //
6 // 04/25/2008 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
8 
9 
13 
17 
19 // initialize static member variables
21 
22 const unsigned int HLTPrescaler::prescaleSeed_ = 65537;
23 
25 // construction/destruction
27 
28 //_____________________________________________________________________________
30  : prescaleFactor_(1)
31  , eventCount_(0)
32  , acceptCount_(0)
33  , offsetCount_(0)
34  , offsetPhase_(iConfig.existsAs<unsigned int>("offset") ? iConfig.getParameter<unsigned int>("offset") : 0)
35  , prescaleService_(0)
36  , newLumi_(true)
37  , gtDigi_ (iConfig.getParameter<edm::InputTag>("L1GtReadoutRecordTag"))
38 {
41  else
42  LogDebug("NoPrescaleService")<<"PrescaleService unavailable, prescaleFactor=1!";
43 }
44 
45 //_____________________________________________________________________________
47 {
48 
49 }
50 
51 
53 // implementation of member functions
55 
56 //______________________________________________________________________________
58  edm::EventSetup const& iSetup)
59 {
60  newLumi_ = true;
61 
62  return true;
63 }
64 
65 
66 //_____________________________________________________________________________
68 {
69  // during the first event of a LumiSection, read from the GT the prescale index for this
70  // LumiSection and get the corresponding prescale factor from the PrescaleService
71  if (newLumi_) {
72  newLumi_ = false;
73 
74  bool needsInit (eventCount_==0);
75 
76  if (prescaleService_) {
77  const unsigned int oldPrescale(prescaleFactor_);
78 
80  iEvent.getByLabel(gtDigi_ , handle);
81  if (handle.isValid()) {
82  unsigned int index = handle->gtFdlWord().gtPrescaleFactorIndexAlgo();
83  // gtPrescaleFactorIndexTech() is also available
84  // by construction, they should always return the same index
86  } else {
87  edm::LogWarning("HLT") << "Cannot read prescale column index from GT data: using default as defined by configuration or DAQ";
89  }
90 
91  if (prescaleFactor_ != oldPrescale) {
92  edm::LogInfo("ChangedPrescale")
93  << "lumiBlockNb="<< iEvent.getLuminosityBlock().id().luminosityBlock() << ", "
94  << "path="<<*pathName()<<": "
95  << prescaleFactor_ << " [" <<oldPrescale<<"]";
96  // reset the prescale counter
97  needsInit = true;
98  }
99  }
100 
101  if (needsInit && (prescaleFactor_ != 0)) {
102  // initialize the prescale counter to the first event number multiplied by a big "seed"
104  }
105  }
106 
107  const bool result ( (prescaleFactor_ == 0) ?
108  false : ((eventCount_ + offsetCount_) % prescaleFactor_ == 0) );
109 
110  ++eventCount_;
111  if (result) ++acceptCount_;
112  return result;
113 }
114 
115 
116 //_____________________________________________________________________________
118 {
119  edm::LogInfo("PrescaleSummary")
120  << acceptCount_<< "/" <<eventCount_
121  << " ("
122  << 100.*acceptCount_/static_cast<double>(std::max(1u,eventCount_))
123  << "% of events accepted).";
124  return;
125 }
#define LogDebug(id)
LuminosityBlockID id() const
EventNumber_t event() const
Definition: EventID.h:44
virtual ~HLTPrescaler()
Definition: HLTPrescaler.cc:46
edm::InputTag gtDigi_
GT payload, to extract the prescale column index.
Definition: HLTPrescaler.h:65
edm::service::PrescaleService * prescaleService_
prescale service
Definition: HLTPrescaler.h:59
virtual void endJob()
static const unsigned int prescaleSeed_
&quot;seed&quot; used to initialize the prescale counter
Definition: HLTPrescaler.h:69
const std::string * pathName() const
Definition: HLTadd.h:31
unsigned int acceptCount_
accept counter
Definition: HLTPrescaler.h:52
int iEvent
Definition: GenABIO.cc:243
virtual bool beginLuminosityBlock(edm::LuminosityBlock &lb, edm::EventSetup const &iSetup)
Definition: HLTPrescaler.cc:57
const T & max(const T &a, const T &b)
tuple result
Definition: query.py:137
unsigned int eventCount_
event counter
Definition: HLTPrescaler.h:49
bool newLumi_
check for (re)initialization of the prescale
Definition: HLTPrescaler.h:62
tuple handle
Definition: patZpeak.py:22
LuminosityBlock const & getLuminosityBlock() const
Definition: Event.h:58
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
unsigned long long uint64_t
Definition: Time.h:15
LuminosityBlockNumber_t luminosityBlock() const
unsigned int offsetPhase_
Definition: HLTPrescaler.h:56
edm::EventID id() const
Definition: EventBase.h:56
HLTPrescaler(edm::ParameterSet const &iConfig)
Definition: HLTPrescaler.cc:29
unsigned int offsetCount_
initial offset
Definition: HLTPrescaler.h:55
virtual bool filter(edm::Event &iEvent, edm::EventSetup const &iSetup)
Definition: HLTPrescaler.cc:67
unsigned int prescaleFactor_
accept one in prescaleFactor_; 0 means never to accept an event
Definition: HLTPrescaler.h:46
unsigned int getPrescale(unsigned int lvl1Index, std::string const &prescaledPath)