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::one::EDFilter<> {
47 public:
48  explicit PrescalerFHN(const edm::ParameterSet&);
49 
50 private:
51  bool filter(edm::Event&, const edm::EventSetup&) override;
52  // ----------member data ---------------------------
53 
55 
57 
59 
60  std::map<std::string, unsigned int> prescales;
61  std::map<std::string, unsigned int> prescale_counter;
62 
63  std::map<std::string, unsigned int> trigger_indices;
64 };
65 
66 //
67 // constants, enums and typedefs
68 //
69 
70 //
71 // static data member definitions
72 //
73 
74 //
75 // constructors and destructor
76 //
78  tok_trigger = consumes<TriggerResults>(iConfig.getParameter<edm::InputTag>("TriggerResultsTag"));
79  //now do what ever initialization is needed
80  std::vector<edm::ParameterSet> prescales_in(iConfig.getParameter<std::vector<edm::ParameterSet> >("Prescales"));
81 
82  for (std::vector<edm::ParameterSet>::const_iterator cit = prescales_in.begin(); cit != prescales_in.end(); cit++) {
83  std::string name(cit->getParameter<std::string>("HLTName"));
84  unsigned int factor(cit->getParameter<unsigned int>("PrescaleFactor"));
85 
86  // does some exception get thrown if parameters aren't available? should test...
87 
89  prescale_counter[name] = 0;
90  }
91 }
92 
93 //
94 // member functions
95 //
96 
98  trigger_indices.clear();
99 
100  for (std::map<std::string, unsigned int>::const_iterator cit = prescales.begin(); cit != prescales.end(); cit++) {
101  auto index = triggerNames.triggerIndex(cit->first);
102  if (index < result.size()) {
103  trigger_indices[cit->first] = index;
104  } else {
105  // trigger path not found
106  LogDebug("") << "requested HLT path does not exist: " << cit->first;
107  }
108  }
109 }
110 
111 // ------------ method called on each new Event ------------
113  using namespace edm;
114 
115  /* Goal for this skim:
116  Prescaling MET HLT paths
117  Option to turn off HSCP filter
118  - Doing that by treating it as an HLT with prescale 1
119  */
120 
121  // Trying to mirror HLTrigger/HLTfilters/src/HLTHighLevel.cc where possible
122 
123  auto const& trh = iEvent.getHandle(tok_trigger);
124 
125  if (trh.isValid()) {
126  LogDebug("") << "TriggerResults found, number of HLT paths: " << trh->size();
127  } else {
128  LogDebug("") << "TriggerResults product not found - returning result=false!";
129  return false;
130  }
131 
132  const edm::TriggerNames& triggerNames = iEvent.triggerNames(*trh);
133  if (triggerNamesID_ != triggerNames.parameterSetID()) {
134  triggerNamesID_ = triggerNames.parameterSetID();
135  init(*trh, triggerNames);
136  }
137 
138  // Trigger indices are ready at this point
139  // - Begin checking for HLT bits
140 
141  bool accept_event = false;
142  for (std::map<std::string, unsigned int>::const_iterator cit = trigger_indices.begin(); cit != trigger_indices.end();
143  cit++) {
144  if (trh->accept(cit->second)) {
145  prescale_counter[cit->first]++;
146  if (prescale_counter[cit->first] >= prescales[cit->first]) {
147  accept_event = true;
148  prescale_counter[cit->first] = 0;
149  }
150  }
151  }
152 
153  return accept_event;
154 }
155 
156 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int init
Definition: HydjetWrapper.h:66
edm::ParameterSetID triggerNamesID_
Definition: PrescalerFHN.cc:56
void init(const edm::TriggerResults &, const edm::TriggerNames &triggerNames)
Definition: PrescalerFHN.cc:97
int iEvent
Definition: GenABIO.cc:224
std::map< std::string, unsigned int > prescales
Definition: PrescalerFHN.cc:60
std::map< std::string, unsigned int > prescale_counter
Definition: PrescalerFHN.cc:61
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< std::string, unsigned int > trigger_indices
Definition: PrescalerFHN.cc:63
edm::EDGetTokenT< TriggerResults > tok_trigger
Definition: PrescalerFHN.cc:58
bool filter(edm::Event &, const edm::EventSetup &) override
PrescalerFHN(const edm::ParameterSet &)
Definition: PrescalerFHN.cc:77
HLT enums.
#define LogDebug(id)