CMS 3D CMS Logo

PrescaleEventFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CalibTracker/SiStripCommon
4 // Class: PrescaleEventFilter
5 //
13 //
14 // Original Author: Marco Musich
15 // Created: Wed, 29 Nov 2017 15:27:07 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
32 
33 namespace prescale {
34  struct Efficiency {
36  mutable std::atomic<unsigned int> eventCount_;
37  mutable std::atomic<unsigned int> acceptCount_;
38  };
39 }
40 
41 //
42 // class declaration
43 //
44 
45 class PrescaleEventFilter : public edm::stream::EDFilter<edm::GlobalCache<prescale::Efficiency> > {
46  public:
48  ~PrescaleEventFilter() override;
49 
50  static std::unique_ptr<prescale::Efficiency> initializeGlobalCache(edm::ParameterSet const&) {
51  return std::unique_ptr<prescale::Efficiency>(new prescale::Efficiency());
52  };
53 
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56  static void globalEndJob(const prescale::Efficiency* efficiency);
57 
58  private:
59  void beginStream(edm::StreamID) override;
60  bool filter(edm::Event&, const edm::EventSetup&) override;
61  void endStream() override;
62 
63  // ----------member data ---------------------------
64 
66  unsigned int prescaleFactor_;
67 
69  unsigned int eventCount_;
70 
72  unsigned int acceptCount_;
73 
75  unsigned int offsetCount_;
76  unsigned int offsetPhase_;
77 
79  bool newLumi_;
80 
82  static const unsigned int prescaleSeed_ = 65537;
83 
84 };
85 
86 //
87 // constructors and destructor
88 //
90  prescaleFactor_(iConfig.getParameter<unsigned int>("prescale")),
91  eventCount_(0),
92  acceptCount_(0),
93  offsetCount_(0),
94  offsetPhase_(iConfig.getParameter<unsigned int>("offset"))
95 {
96  //now do what ever initialization is needed
97 
98 }
99 
100 
102 {
103 }
104 
105 //
106 // member functions
107 //
108 
109 // ------------ method called on each new Event ------------
110 bool
112 {
113  using namespace edm;
114 
115  bool needsInit (eventCount_==0);
116 
117  if (needsInit && (prescaleFactor_ != 0)) {
118  // initialize the prescale counter to the first event number multiplied by a big "seed"
120  }
121 
122  const bool result ( (prescaleFactor_ == 0) ?
123  false : ((eventCount_ + offsetCount_) % prescaleFactor_ == 0) );
124 
125  ++eventCount_;
126  if (result) ++acceptCount_;
127  return result;
128 
129 }
130 
131 // ------------ method called once each stream before processing any runs, lumis or events ------------
132 void
134 {
135 }
136 
137 //_____________________________________________________________________________
138 void
140 {
141  //since these are std::atomic, it is safe to increment them
142  // even if multiple endStreams are being called.
143  globalCache()->eventCount_ += eventCount_;
144  globalCache()->acceptCount_ += acceptCount_;
145  return;
146 }
147 
148 //_____________________________________________________________________________
149 void
151 {
152  unsigned int accept(efficiency->acceptCount_);
153  unsigned int event (efficiency->eventCount_);
154  edm::LogInfo("PrescaleSummary")
155  << accept << "/" << event
156  << " ("
157  << 100.*accept/static_cast<double>(std::max(1u,event))
158  << "% of events accepted).";
159  return;
160 }
161 
162 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
163 void
165 
167  desc.add<unsigned int>("prescale",1);
168  desc.add<unsigned int>("offset",0);
169  descriptions.add("prescaleEvent", desc);
170 
171 }
172 //define this as a plug-in
EventNumber_t event() const
Definition: EventID.h:41
bool newLumi_
check for (re)initialization of the prescale
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
PrescaleEventFilter(edm::ParameterSet const &iConfig, const prescale::Efficiency *efficiency)
unsigned int acceptCount_
accept counter
static std::unique_ptr< prescale::Efficiency > initializeGlobalCache(edm::ParameterSet const &)
unsigned int prescaleFactor_
accept one in prescaleFactor_; 0 means never to accept an event
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool filter(edm::Event &, const edm::EventSetup &) override
static const unsigned int prescaleSeed_
"seed" used to initialize the prescale counter
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
unsigned int offsetCount_
initial offset
void beginStream(edm::StreamID) override
unsigned int eventCount_
event counter
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned long long uint64_t
Definition: Time.h:15
static void globalEndJob(const prescale::Efficiency *efficiency)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
std::atomic< unsigned int > acceptCount_
std::atomic< unsigned int > eventCount_