-
Notifications
You must be signed in to change notification settings - Fork 4.2k
/
HcalLaserEventFilter.cc
372 lines (301 loc) · 14 KB
/
HcalLaserEventFilter.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
// -*- C++ -*-
//
// Package: HcalLaserEventFilter
// Class: HcalLaserEventFilter
//
/**\class HcalLaserEventFilter HcalLaserEventFilter.cc UserCode/HcalLaserEventFilter/src/HcalLaserEventFilter.cc
Description: Filter for removing Hcal laser events
Implementation:
Filter allows users to remove specific (run,event) pairs that have been identified as noise
It also allows users to remove events in which the number of HBHE rechits exceeds a given threshold (5000 by default).
*/
//
// Original Author: Jeff Temple (temple@cern.ch)
// Created: Thu Nov 17 12:44:22 EST 2011
//
//
// system include files
#include <memory>
#include <sstream>
#include <iostream>
#include <string>
// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDFilter.h"
#include "FWCore/Version/interface/GetReleaseVersion.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
// Use for HBHERecHitCollection
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/METReco/interface/HcalNoiseSummary.h"
#include "DataFormats/Provenance/interface/RunLumiEventNumber.h"
//
// class declaration
//
class HcalLaserEventFilter : public edm::global::EDFilter<> {
public:
explicit HcalLaserEventFilter(const edm::ParameterSet&);
~HcalLaserEventFilter() override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
std::vector<int> GetCMSSWVersion(std::string const&) const;
bool IsGreaterThanMinCMSSWVersion(std::vector<int> const&) const;
// ----------member data ---------------------------
// Filter option 1: veto events by run, event number
const bool vetoByRunEventNumber_;
std::vector<std::pair<edm::RunNumber_t,edm::EventNumber_t> > RunEventData_;
// Filter option 2: veto events by HBHE occupancy
const bool vetoByHBHEOccupancy_;
const unsigned int minOccupiedHBHE_;
const bool vetoByLaserMonitor_;
const double minLaserMonitorCharge_;
// Allow for debugging information to be printed
const bool debug_;
// Reverse the filter decision (so instead of selecting only good events, it
// will select only events that fail the filter conditions -- useful for studying
// bad events.)
const bool reverseFilter_;
// InputTag for HBHE rechits
const edm::InputTag hbheInputLabel_;
edm::EDGetTokenT<HBHERecHitCollection> hbheToken_;
const edm::InputTag hcalNoiseSummaryLabel_;
edm::EDGetTokenT<HcalNoiseSummary> hcalNoiseSummaryToken_;
const bool taggingMode_;
// Decide whether to use the HcalNoiseSummary to get the RecHit info, or to use the RecHit Collection itself
bool useHcalNoiseSummary_;
bool forceUseRecHitCollection_;
bool forceUseHcalNoiseSummary_;
std::vector<int> minVersion_;
};
//
// constants, enums and typedefs
//
//
// static data member definitions
//
//
// constructors and destructor
//
HcalLaserEventFilter::HcalLaserEventFilter(const edm::ParameterSet& iConfig)
// Get values from python cfg file
: vetoByRunEventNumber_ (iConfig.getParameter<bool>("vetoByRunEventNumber"))
, vetoByHBHEOccupancy_ (iConfig.getParameter<bool>("vetoByHBHEOccupancy"))
, minOccupiedHBHE_ (iConfig.getParameter<unsigned int>("minOccupiedHBHE"))
, vetoByLaserMonitor_ (iConfig.getParameter<bool>("vetoByLaserMonitor"))
, minLaserMonitorCharge_ (iConfig.getParameter<double>("minLaserMonitorCharge"))
, debug_ (iConfig.getParameter<bool>("debug"))
, reverseFilter_ (iConfig.getParameter<bool>("reverseFilter"))
, hbheInputLabel_ (iConfig.getUntrackedParameter<edm::InputTag>("hbheInputLabel",edm::InputTag("hbhereco")))
, hbheToken_ (mayConsume<HBHERecHitCollection>(hbheInputLabel_))
, hcalNoiseSummaryLabel_ (iConfig.getUntrackedParameter<edm::InputTag>("hcalNoiseSummaryLabel",edm::InputTag("hcalnoise")))
, hcalNoiseSummaryToken_ (mayConsume<HcalNoiseSummary>(hcalNoiseSummaryLabel_))
, taggingMode_ (iConfig.getParameter<bool>("taggingMode"))
, forceUseRecHitCollection_ (iConfig.getParameter<bool>("forceUseRecHitCollection"))
, forceUseHcalNoiseSummary_ (iConfig.getParameter<bool>("forceUseHcalNoiseSummary"))
{
std::vector<unsigned int> temprunevt = iConfig.getParameter<std::vector<unsigned int> >("BadRunEventNumbers");
// Make (run,evt) pairs for storing bad events
// Make this a map for better search performance?
for (unsigned int i=0;i+1<temprunevt.size();i+=2)
{
edm::RunNumber_t run=temprunevt[i];
edm::EventNumber_t evt=temprunevt[i+1];
RunEventData_.push_back(std::make_pair(run,evt));
}
produces<bool>();
// Specify the minimum release that has the rechit counts in the HcalNoiseSummary object.
// If current release >= that release, then HcalNoiseSummary will be used. Otherwise, Rechit collection will be used.
std::string minRelease="\"CMSSW_5_2_0\"";
minVersion_=GetCMSSWVersion(minRelease);
std::vector <int> currentVersion=GetCMSSWVersion(edm::getReleaseVersion());
if (IsGreaterThanMinCMSSWVersion(currentVersion)) // current Version is >= minVersion_
useHcalNoiseSummary_=true;
else
useHcalNoiseSummary_=false;
}
HcalLaserEventFilter::~HcalLaserEventFilter()
{
// do anything here that needs to be done at desctruction time
// (e.g. close files, deallocate resources etc.)
}
//
// member functions
//
// ------------ method called on each new Event ------------
bool
HcalLaserEventFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const
{
using namespace edm;
bool filterDecision=true;
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter> Run = "<<iEvent.id().run()<<" Event = "<<iEvent.id().event();
// Veto events by run/event numbers
if (vetoByRunEventNumber_)
{
for (unsigned int i=0;i<RunEventData_.size();++i)
{
if (iEvent.id().run()==RunEventData_[i].first &&
iEvent.id().event()==RunEventData_[i].second)
{
if (debug_) edm::LogInfo("HcalLaserEventFilter")<<"\t<HcalLaserEventFilter> Filtering bad event; Run "<<iEvent.id().run()<<" Event = "<<iEvent.id().event();
filterDecision=false;
break;
}
}
} // if (vetoByRunEventNumber_)
//Veto events by HBHE rechit collection size
if (vetoByHBHEOccupancy_)
{
// The decision on whether or not to use the noise summary is made based on the CMSSW version.
// As of CMSSW_5_2_0, the HcalNoiseSummary contains the total number of HBHE hits, as well as number of hits > 1.5 GeV and some other info.
// The boolean 'forceUseRecHitCollection_' can be used to override this automatic behavior, and to use the RecHit collection itself, regardless of CMSSW version.
//////////////////////////////////////////////////////////
//
// Apply Filtering based on RecHit information in HBHERecHitcollection
//
////////////////////////////////////////////////////////////
if (useHcalNoiseSummary_==false || forceUseRecHitCollection_==true)
{
edm::Handle<HBHERecHitCollection> hbheRecHits;
if (iEvent.getByToken(hbheToken_,hbheRecHits))
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"Rechit size = "<<hbheRecHits->size()<<" threshold = "<<minOccupiedHBHE_;
if (hbheRecHits->size()>=minOccupiedHBHE_)
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter> Filtering because of large HBHE rechit size; "<<hbheRecHits->size()<<" rechits is greater than or equal to the allowed maximum of "<<minOccupiedHBHE_;
filterDecision=false;
}
}
else
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter::Error> No valid HBHERecHitCollection with label '"<<hbheInputLabel_<<"' found";
}
}
//////////////////////////////////////////////////////////
//
// Apply Filtering based on RecHit information in HcalNoiseSummary object
//
////////////////////////////////////////////////////////////
else if (useHcalNoiseSummary_==true || forceUseHcalNoiseSummary_==true)
{
Handle<HcalNoiseSummary> hSummary;
if (iEvent.getByToken(hcalNoiseSummaryToken_,hSummary)) // get by label, usually with label 'hcalnoise'
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") << " RECHIT SIZE (from HcalNoiseSummary) = "<<hSummary->GetRecHitCount()<<" threshold = "<<minOccupiedHBHE_;
if (hSummary->GetRecHitCount() >= (int)minOccupiedHBHE_)
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter> Filtering because of large HBHE rechit size in HcalNoiseSummary; "<<hSummary->GetRecHitCount()<<" rechits is greater than or equal to the allowed maximum of "<<minOccupiedHBHE_;
filterDecision=false;
}
}
else
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter::Error> No valid HcalNoiseSummary with label '"<<hcalNoiseSummaryLabel_<<"' found";
}
}
}// if (vetoByHBHEOccupancy_)
if( vetoByLaserMonitor_ )
{
//////////////////////////////////////////////////////////
//
// Apply Filtering based on laser monitor information in HcalNoiseSummary object
//
////////////////////////////////////////////////////////////
Handle<HcalNoiseSummary> hSummary;
if (iEvent.getByToken(hcalNoiseSummaryToken_,hSummary)) // get by label, usually with label 'hcalnoise'
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") << " LASERMON CHARGE (from HcalNoiseSummary) = "<<hSummary->GetLaserMonitorCharge()<<" threshold = "<<minLaserMonitorCharge_;
if( hSummary->GetLaserMonitorCharge() > minLaserMonitorCharge_ )
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter> Filtering because of large Laser monitor charge in HcalNoiseSummary; "<<hSummary->GetLaserMonitorCharge()<<" charge is greater than or equal to the allowed maximum of "<<minLaserMonitorCharge_;
filterDecision=false;
}
}
else
{
if (debug_) edm::LogInfo("HcalLaserEventFilter") <<"<HcalLaserEventFilter::Error> No valid HcalNoiseSummary with label '"<<hcalNoiseSummaryLabel_<<"' found";
}
}
// Reverse decision, if specified by user
if (reverseFilter_)
filterDecision=!filterDecision;
iEvent.put(std::make_unique<bool>(filterDecision));
return taggingMode_ || filterDecision;
}
// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void
HcalLaserEventFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
//The following says we do not know what parameters are allowed so do no validation
// Please change this to state exactly what you do use, even if it is no parameters
edm::ParameterSetDescription desc;
desc.add<bool>("vetoByRunEventNumber",false)->
setComment("Enable filtering by run number");
desc.add<bool>("vetoByHBHEOccupancy",true)->
setComment("Enable occupancy filtering");
desc.add<unsigned int>("minOccupiedHBHE",4000)->
setComment("Minimum occupancy to filter events");
desc.add<bool>("vetoByLaserMonitor",true)->
setComment("Enable Laser monitoring filtering");
desc.add<double>("minLaserMonitorCharge",5000.)->
setComment("Set minimum laser monitor charge to filter events");
desc.add<bool>("debug",false)->
setComment("Enable debugging messages");
desc.add<bool>("reverseFilter",false)->
setComment("Invert filter decision");
desc.add<bool>("taggingMode", false)->
setComment("do not filter, just tag the event");
desc.add<bool>("forceUseRecHitCollection",false)->
setComment("force the evaluation using RecHit collection");
desc.add<bool>("forceUseHcalNoiseSummary",false)->
setComment("force the evaluation using Noise Summary");
desc.add<std::vector<unsigned int> >("BadRunEventNumbers",{})->
setComment("vector of bad events to filter");
descriptions.add("hcallaserevent", desc);
}
std::vector<int> HcalLaserEventFilter::GetCMSSWVersion(std::string const& instring) const
{
std::vector <int> temp;
std::string prefix;
std::string v1, v2, v3;
std::istringstream oss(instring);
getline(oss, prefix,'_');
getline(oss, v1,'_');
getline(oss, v2,'_');
getline(oss, v3,'_');
std::stringstream buffer(v1);
int t;
buffer>>t;
temp.push_back(t);
buffer.str();
buffer<<v2;
buffer>>t;
temp.push_back(t);
buffer.str();
buffer<<v3;
buffer>>t;
temp.push_back(t);
//std::cout <<"PREFIX = "<<prefix<<" "<<temp[0]<<" "<<temp[1]<<" "<<temp[2]<<std::endl;
//( ex: PREFIX = "CMSSW 5 5 5 )
return temp;
}
bool HcalLaserEventFilter::IsGreaterThanMinCMSSWVersion(std::vector<int> const& currentVersion) const
{
// Returns false if current version is less than min version
// Otherwise, returns true
// Assumes CMSSW versioning X_Y_Z
// Compare X
if (currentVersion[0]<minVersion_[0]) return false;
if (currentVersion[0]>minVersion_[0]) return true;
// If neither is true, first value of CMSSW versions are the same
// Compare Y
if (currentVersion[1]<minVersion_[1]) return false;
if (currentVersion[1]>minVersion_[1]) return true;
// Compare Z
if (currentVersion[2]<minVersion_[2]) return false;
if (currentVersion[2]>minVersion_[2]) return true;
return true; // versions are identical
}
//define this as a plug-in
DEFINE_FWK_MODULE(HcalLaserEventFilter);