CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // $Id: PrescalerFHN.cc,v 1.4 2010/02/16 22:24:16 wdd Exp $
17 //
18 //
19 // modified to PrecalerFHN by Grigory Safronov 27/03/09
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
32 
34 
38 
39 #include <string>
40 
41 //
42 // class declaration
43 //
44 
45 using namespace edm;
46 
47 class PrescalerFHN : public edm::EDFilter {
48  public:
49  explicit PrescalerFHN(const edm::ParameterSet&);
50  ~PrescalerFHN();
51 
52  private:
53  virtual void beginJob() ;
54  virtual bool filter(edm::Event&, const edm::EventSetup&);
55  virtual void endJob() ;
56  // ----------member data ---------------------------
57 
58  void init(const edm::TriggerResults &,
59  const edm::TriggerNames & triggerNames);
60 
62 
64 
65  std::map<std::string, unsigned int> prescales;
66  std::map<std::string, unsigned int> prescale_counter;
67 
68  std::map<std::string, unsigned int> trigger_indices;
69 };
70 
71 //
72 // constants, enums and typedefs
73 //
74 
75 //
76 // static data member definitions
77 //
78 
79 //
80 // constructors and destructor
81 //
83  : triggerTag(iConfig.getParameter<edm::InputTag>("TriggerResultsTag"))
84 {
85  //now do what ever initialization is needed
86  std::vector<edm::ParameterSet> prescales_in(iConfig.getParameter<std::vector<edm::ParameterSet> >("Prescales"));
87 
88  for (std::vector<edm::ParameterSet>::const_iterator cit = prescales_in.begin();
89  cit != prescales_in.end(); cit++) {
90 
91  std::string name(cit->getParameter<std::string>("HLTName"));
92  unsigned int factor(cit->getParameter<unsigned int>("PrescaleFactor"));
93 
94  // does some exception get thrown if parameters aren't available? should test...
95 
96  prescales[name] = factor;
98  }
99 
100 }
101 
102 
104 {
105 
106  // do anything here that needs to be done at desctruction time
107  // (e.g. close files, deallocate resources etc.)
108 
109 }
110 
111 
112 //
113 // member functions
114 //
115 
117  const edm::TriggerNames & triggerNames)
118 {
119  trigger_indices.clear();
120 
121  for (std::map<std::string, unsigned int>::const_iterator cit = prescales.begin();
122  cit != prescales.end(); cit++) {
123 
124  trigger_indices[cit->first] = triggerNames.triggerIndex(cit->first);
125 
126  if (trigger_indices[cit->first] >= result.size()) {
127  // trigger path not found
128  LogDebug("") << "requested HLT path does not exist: " << cit->first;
129  trigger_indices.erase(cit->first);
130  }
131 
132  }
133 }
134 
135 // ------------ method called on each new Event ------------
136 bool
138 {
139  using namespace edm;
140 
141  /* Goal for this skim:
142  Prescaling MET HLT paths
143  Option to turn off HSCP filter
144  - Doing that by treating it as an HLT with prescale 1
145  */
146 
147  // Trying to mirror HLTrigger/HLTfilters/src/HLTHighLevel.cc where possible
148 
150  iEvent.getByLabel(triggerTag, trh);
151 
152  if (trh.isValid()) {
153  LogDebug("") << "TriggerResults found, number of HLT paths: " << trh->size();
154  } else {
155  LogDebug("") << "TriggerResults product not found - returning result=false!";
156  return false;
157  }
158 
159  const edm::TriggerNames & triggerNames = iEvent.triggerNames(*trh);
160  if (triggerNamesID_ != triggerNames.parameterSetID()) {
161  triggerNamesID_ = triggerNames.parameterSetID();
162  init(*trh, triggerNames);
163  }
164 
165  // Trigger indices are ready at this point
166  // - Begin checking for HLT bits
167 
168  bool accept_event = false;
169  for (std::map<std::string, unsigned int>::const_iterator cit = trigger_indices.begin();
170  cit != trigger_indices.end(); cit++) {
171  if (trh->accept(cit->second)) {
172  prescale_counter[cit->first]++;
173  if (prescale_counter[cit->first] >= prescales[cit->first]) {
174  accept_event = true;
175  prescale_counter[cit->first] = 0;
176  }
177  }
178  }
179 
180  return accept_event;
181 }
182 
183 // ------------ method called once each job just before starting event loop ------------
184 void
186 {
187 }
188 
189 // ------------ method called once each job just after ending the event loop ------------
190 void
192 }
193 
194 //define this as a plug-in
#define LogDebug(id)
T getParameter(std::string const &) const
virtual void beginJob()
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:207
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int init
Definition: HydjetWrapper.h:63
edm::ParameterSetID triggerNamesID_
Definition: PrescalerFHN.cc:61
void init(const edm::TriggerResults &, const edm::TriggerNames &triggerNames)
void beginJob()
Definition: Breakpoints.cc:15
int iEvent
Definition: GenABIO.cc:243
std::map< std::string, unsigned int > prescales
Definition: PrescalerFHN.cc:65
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
std::map< std::string, unsigned int > prescale_counter
Definition: PrescalerFHN.cc:66
unsigned int size() const
Get number of paths stored.
edm::InputTag triggerTag
Definition: PrescalerFHN.cc:63
tuple result
Definition: query.py:137
std::map< std::string, unsigned int > trigger_indices
Definition: PrescalerFHN.cc:68
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
PrescalerFHN(const edm::ParameterSet &)
Definition: PrescalerFHN.cc:82
virtual void endJob()
virtual bool filter(edm::Event &, const edm::EventSetup &)