CMS 3D CMS Logo

PrescalerFHN.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HCALNoiseAlCaReco
4 // Class: HCALNoiseAlCaReco
5 //
13 //
14 // Original Author: Kenneth Case Rossato
15 // Created: Wed Mar 25 13:05:10 CET 2008
16 //
17 //
18 // modified to PrecalerFHN by Grigory Safronov 27/03/09
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
33 
37 
38 #include <string>
39 
40 //
41 // class declaration
42 //
43 
44 using namespace edm;
45 
46 class PrescalerFHN : public edm::EDFilter {
47 public:
48  explicit PrescalerFHN(const edm::ParameterSet&);
49  ~PrescalerFHN() override;
50 
51 private:
52  void beginJob() override;
53  bool filter(edm::Event&, const edm::EventSetup&) override;
54  void endJob() override;
55  // ----------member data ---------------------------
56 
57  void init(const edm::TriggerResults&, const edm::TriggerNames& triggerNames);
58 
60 
62 
63  std::map<std::string, unsigned int> prescales;
64  std::map<std::string, unsigned int> prescale_counter;
65 
66  std::map<std::string, unsigned int> trigger_indices;
67 };
68 
69 //
70 // constants, enums and typedefs
71 //
72 
73 //
74 // static data member definitions
75 //
76 
77 //
78 // constructors and destructor
79 //
81  tok_trigger = consumes<TriggerResults>(iConfig.getParameter<edm::InputTag>("TriggerResultsTag"));
82  //now do what ever initialization is needed
83  std::vector<edm::ParameterSet> prescales_in(iConfig.getParameter<std::vector<edm::ParameterSet> >("Prescales"));
84 
85  for (std::vector<edm::ParameterSet>::const_iterator cit = prescales_in.begin(); cit != prescales_in.end(); cit++) {
86  std::string name(cit->getParameter<std::string>("HLTName"));
87  unsigned int factor(cit->getParameter<unsigned int>("PrescaleFactor"));
88 
89  // does some exception get thrown if parameters aren't available? should test...
90 
91  prescales[name] = factor;
92  prescale_counter[name] = 0;
93  }
94 }
95 
97  // do anything here that needs to be done at desctruction time
98  // (e.g. close files, deallocate resources etc.)
99 }
100 
101 //
102 // member functions
103 //
104 
106  trigger_indices.clear();
107 
108  for (std::map<std::string, unsigned int>::const_iterator cit = prescales.begin(); cit != prescales.end(); cit++) {
109  trigger_indices[cit->first] = triggerNames.triggerIndex(cit->first);
110 
111  if (trigger_indices[cit->first] >= result.size()) {
112  // trigger path not found
113  LogDebug("") << "requested HLT path does not exist: " << cit->first;
114  trigger_indices.erase(cit->first);
115  }
116  }
117 }
118 
119 // ------------ method called on each new Event ------------
121  using namespace edm;
122 
123  /* Goal for this skim:
124  Prescaling MET HLT paths
125  Option to turn off HSCP filter
126  - Doing that by treating it as an HLT with prescale 1
127  */
128 
129  // Trying to mirror HLTrigger/HLTfilters/src/HLTHighLevel.cc where possible
130 
132  iEvent.getByToken(tok_trigger, trh);
133 
134  if (trh.isValid()) {
135  LogDebug("") << "TriggerResults found, number of HLT paths: " << trh->size();
136  } else {
137  LogDebug("") << "TriggerResults product not found - returning result=false!";
138  return false;
139  }
140 
141  const edm::TriggerNames& triggerNames = iEvent.triggerNames(*trh);
142  if (triggerNamesID_ != triggerNames.parameterSetID()) {
143  triggerNamesID_ = triggerNames.parameterSetID();
144  init(*trh, triggerNames);
145  }
146 
147  // Trigger indices are ready at this point
148  // - Begin checking for HLT bits
149 
150  bool accept_event = false;
151  for (std::map<std::string, unsigned int>::const_iterator cit = trigger_indices.begin(); cit != trigger_indices.end();
152  cit++) {
153  if (trh->accept(cit->second)) {
154  prescale_counter[cit->first]++;
155  if (prescale_counter[cit->first] >= prescales[cit->first]) {
156  accept_event = true;
157  prescale_counter[cit->first] = 0;
158  }
159  }
160  }
161 
162  return accept_event;
163 }
164 
165 // ------------ method called once each job just before starting event loop ------------
167 
168 // ------------ method called once each job just after ending the event loop ------------
170 
171 //define this as a plug-in
#define LogDebug(id)
T getParameter(std::string const &) const
void beginJob() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool accept() const
Has at least one path accepted the event?
int init
Definition: HydjetWrapper.h:67
edm::ParameterSetID triggerNamesID_
Definition: PrescalerFHN.cc:59
void init(const edm::TriggerResults &, const edm::TriggerNames &triggerNames)
void beginJob()
Definition: Breakpoints.cc:14
ParameterSetID const & parameterSetID() const
Definition: TriggerNames.cc:33
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< std::string, unsigned int > prescales
Definition: PrescalerFHN.cc:63
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:24
std::map< std::string, unsigned int > prescale_counter
Definition: PrescalerFHN.cc:64
unsigned int size() const
Get number of paths stored.
void endJob() override
std::map< std::string, unsigned int > trigger_indices
Definition: PrescalerFHN.cc:66
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< TriggerResults > tok_trigger
Definition: PrescalerFHN.cc:61
~PrescalerFHN() override
Definition: PrescalerFHN.cc:96
bool filter(edm::Event &, const edm::EventSetup &) override
PrescalerFHN(const edm::ParameterSet &)
Definition: PrescalerFHN.cc:80
HLT enums.
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:256