File indexing completed on 2023-10-25 10:01:00
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022 #include <memory>
0023 #include <sstream>
0024 #include <iostream>
0025 #include <string>
0026
0027
0028 #include "FWCore/Framework/interface/Frameworkfwd.h"
0029 #include "FWCore/Framework/interface/global/EDFilter.h"
0030 #include "FWCore/Version/interface/GetReleaseVersion.h"
0031
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/MakerMacros.h"
0034
0035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0037
0038
0039 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0040
0041 #include "DataFormats/METReco/interface/HcalNoiseSummary.h"
0042
0043 #include "DataFormats/Provenance/interface/RunLumiEventNumber.h"
0044
0045
0046
0047
0048 class HcalLaserEventFilter : public edm::global::EDFilter<> {
0049 public:
0050 explicit HcalLaserEventFilter(const edm::ParameterSet&);
0051 ~HcalLaserEventFilter() override;
0052
0053 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0054
0055 private:
0056 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0057
0058 std::vector<int> GetCMSSWVersion(std::string const&) const;
0059 bool IsGreaterThanMinCMSSWVersion(std::vector<int> const&) const;
0060
0061
0062
0063
0064 const bool vetoByRunEventNumber_;
0065 std::vector<std::pair<edm::RunNumber_t, edm::EventNumber_t> > RunEventData_;
0066
0067
0068 const bool vetoByHBHEOccupancy_;
0069 const unsigned int minOccupiedHBHE_;
0070
0071 const bool vetoByLaserMonitor_;
0072 const double minLaserMonitorCharge_;
0073
0074
0075 const bool debug_;
0076
0077
0078
0079 const bool reverseFilter_;
0080
0081
0082 const edm::InputTag hbheInputLabel_;
0083 edm::EDGetTokenT<HBHERecHitCollection> hbheToken_;
0084 const edm::InputTag hcalNoiseSummaryLabel_;
0085 edm::EDGetTokenT<HcalNoiseSummary> hcalNoiseSummaryToken_;
0086
0087 const bool taggingMode_;
0088
0089
0090 bool useHcalNoiseSummary_;
0091 bool forceUseRecHitCollection_;
0092 bool forceUseHcalNoiseSummary_;
0093 std::vector<int> minVersion_;
0094 };
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107 HcalLaserEventFilter::HcalLaserEventFilter(const edm::ParameterSet& iConfig)
0108
0109
0110 : vetoByRunEventNumber_(iConfig.getParameter<bool>("vetoByRunEventNumber")),
0111 vetoByHBHEOccupancy_(iConfig.getParameter<bool>("vetoByHBHEOccupancy")),
0112 minOccupiedHBHE_(iConfig.getParameter<unsigned int>("minOccupiedHBHE")),
0113 vetoByLaserMonitor_(iConfig.getParameter<bool>("vetoByLaserMonitor")),
0114 minLaserMonitorCharge_(iConfig.getParameter<double>("minLaserMonitorCharge")),
0115 debug_(iConfig.getParameter<bool>("debug")),
0116 reverseFilter_(iConfig.getParameter<bool>("reverseFilter")),
0117 hbheInputLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hbheInputLabel", edm::InputTag("hbhereco"))),
0118 hbheToken_(mayConsume<HBHERecHitCollection>(hbheInputLabel_)),
0119 hcalNoiseSummaryLabel_(
0120 iConfig.getUntrackedParameter<edm::InputTag>("hcalNoiseSummaryLabel", edm::InputTag("hcalnoise"))),
0121 hcalNoiseSummaryToken_(mayConsume<HcalNoiseSummary>(hcalNoiseSummaryLabel_)),
0122 taggingMode_(iConfig.getParameter<bool>("taggingMode")),
0123 forceUseRecHitCollection_(iConfig.getParameter<bool>("forceUseRecHitCollection")),
0124 forceUseHcalNoiseSummary_(iConfig.getParameter<bool>("forceUseHcalNoiseSummary"))
0125
0126 {
0127 std::vector<unsigned int> temprunevt = iConfig.getParameter<std::vector<unsigned int> >("BadRunEventNumbers");
0128
0129
0130
0131 for (unsigned int i = 0; i + 1 < temprunevt.size(); i += 2) {
0132 edm::RunNumber_t run = temprunevt[i];
0133 edm::EventNumber_t evt = temprunevt[i + 1];
0134 RunEventData_.push_back(std::make_pair(run, evt));
0135 }
0136 produces<bool>();
0137
0138
0139
0140 std::string minRelease = "\"CMSSW_5_2_0\"";
0141
0142 minVersion_ = GetCMSSWVersion(minRelease);
0143 std::vector<int> currentVersion = GetCMSSWVersion(edm::getReleaseVersion());
0144
0145 if (IsGreaterThanMinCMSSWVersion(currentVersion))
0146 useHcalNoiseSummary_ = true;
0147 else
0148 useHcalNoiseSummary_ = false;
0149 }
0150
0151 HcalLaserEventFilter::~HcalLaserEventFilter() {
0152
0153
0154 }
0155
0156
0157
0158
0159
0160
0161 bool HcalLaserEventFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0162 using namespace edm;
0163
0164 bool filterDecision = true;
0165
0166 if (debug_)
0167 edm::LogInfo("HcalLaserEventFilter") << "<HcalLaserEventFilter> Run = " << iEvent.id().run()
0168 << " Event = " << iEvent.id().event();
0169
0170
0171 if (vetoByRunEventNumber_) {
0172 for (unsigned int i = 0; i < RunEventData_.size(); ++i) {
0173 if (iEvent.id().run() == RunEventData_[i].first && iEvent.id().event() == RunEventData_[i].second) {
0174 if (debug_)
0175 edm::LogInfo("HcalLaserEventFilter") << "\t<HcalLaserEventFilter> Filtering bad event; Run "
0176 << iEvent.id().run() << " Event = " << iEvent.id().event();
0177 filterDecision = false;
0178 break;
0179 }
0180 }
0181 }
0182
0183
0184 if (vetoByHBHEOccupancy_) {
0185
0186
0187
0188
0189
0190
0191
0192
0193
0194
0195 if (useHcalNoiseSummary_ == false || forceUseRecHitCollection_ == true) {
0196 edm::Handle<HBHERecHitCollection> hbheRecHits;
0197 if (iEvent.getByToken(hbheToken_, hbheRecHits)) {
0198 if (debug_)
0199 edm::LogInfo("HcalLaserEventFilter")
0200 << "Rechit size = " << hbheRecHits->size() << " threshold = " << minOccupiedHBHE_;
0201 if (hbheRecHits->size() >= minOccupiedHBHE_) {
0202 if (debug_)
0203 edm::LogInfo("HcalLaserEventFilter")
0204 << "<HcalLaserEventFilter> Filtering because of large HBHE rechit size; " << hbheRecHits->size()
0205 << " rechits is greater than or equal to the allowed maximum of " << minOccupiedHBHE_;
0206 filterDecision = false;
0207 }
0208 } else {
0209 if (debug_)
0210 edm::LogInfo("HcalLaserEventFilter")
0211 << "<HcalLaserEventFilter::Error> No valid HBHERecHitCollection with label '" << hbheInputLabel_
0212 << "' found";
0213 }
0214 }
0215
0216
0217
0218
0219
0220
0221 else if (useHcalNoiseSummary_ == true || forceUseHcalNoiseSummary_ == true) {
0222 Handle<HcalNoiseSummary> hSummary;
0223 if (iEvent.getByToken(hcalNoiseSummaryToken_, hSummary))
0224 {
0225 if (debug_)
0226 edm::LogInfo("HcalLaserEventFilter")
0227 << " RECHIT SIZE (from HcalNoiseSummary) = " << hSummary->GetRecHitCount()
0228 << " threshold = " << minOccupiedHBHE_;
0229 if (hSummary->GetRecHitCount() >= (int)minOccupiedHBHE_) {
0230 if (debug_)
0231 edm::LogInfo("HcalLaserEventFilter")
0232 << "<HcalLaserEventFilter> Filtering because of large HBHE rechit size in HcalNoiseSummary; "
0233 << hSummary->GetRecHitCount() << " rechits is greater than or equal to the allowed maximum of "
0234 << minOccupiedHBHE_;
0235 filterDecision = false;
0236 }
0237 } else {
0238 if (debug_)
0239 edm::LogInfo("HcalLaserEventFilter") << "<HcalLaserEventFilter::Error> No valid HcalNoiseSummary with label '"
0240 << hcalNoiseSummaryLabel_ << "' found";
0241 }
0242 }
0243 }
0244 if (vetoByLaserMonitor_) {
0245
0246
0247
0248
0249
0250 Handle<HcalNoiseSummary> hSummary;
0251 if (iEvent.getByToken(hcalNoiseSummaryToken_, hSummary))
0252 {
0253 if (debug_)
0254 edm::LogInfo("HcalLaserEventFilter")
0255 << " LASERMON CHARGE (from HcalNoiseSummary) = " << hSummary->GetLaserMonitorCharge()
0256 << " threshold = " << minLaserMonitorCharge_;
0257 if (hSummary->GetLaserMonitorCharge() > minLaserMonitorCharge_) {
0258 if (debug_)
0259 edm::LogInfo("HcalLaserEventFilter")
0260 << "<HcalLaserEventFilter> Filtering because of large Laser monitor charge in HcalNoiseSummary; "
0261 << hSummary->GetLaserMonitorCharge() << " charge is greater than or equal to the allowed maximum of "
0262 << minLaserMonitorCharge_;
0263 filterDecision = false;
0264 }
0265 } else {
0266 if (debug_)
0267 edm::LogInfo("HcalLaserEventFilter") << "<HcalLaserEventFilter::Error> No valid HcalNoiseSummary with label '"
0268 << hcalNoiseSummaryLabel_ << "' found";
0269 }
0270 }
0271
0272
0273 if (reverseFilter_)
0274 filterDecision = !filterDecision;
0275
0276 iEvent.put(std::make_unique<bool>(filterDecision));
0277
0278 return taggingMode_ || filterDecision;
0279 }
0280
0281
0282 void HcalLaserEventFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0283
0284
0285 edm::ParameterSetDescription desc;
0286
0287 desc.add<bool>("vetoByRunEventNumber", false)->setComment("Enable filtering by run number");
0288 desc.add<bool>("vetoByHBHEOccupancy", true)->setComment("Enable occupancy filtering");
0289 desc.add<unsigned int>("minOccupiedHBHE", 4000)->setComment("Minimum occupancy to filter events");
0290 desc.add<bool>("vetoByLaserMonitor", true)->setComment("Enable Laser monitoring filtering");
0291 desc.add<double>("minLaserMonitorCharge", 5000.)->setComment("Set minimum laser monitor charge to filter events");
0292 desc.add<bool>("debug", false)->setComment("Enable debugging messages");
0293 desc.add<bool>("reverseFilter", false)->setComment("Invert filter decision");
0294 desc.add<bool>("taggingMode", false)->setComment("do not filter, just tag the event");
0295 desc.add<bool>("forceUseRecHitCollection", false)->setComment("force the evaluation using RecHit collection");
0296 desc.add<bool>("forceUseHcalNoiseSummary", false)->setComment("force the evaluation using Noise Summary");
0297 desc.add<std::vector<unsigned int> >("BadRunEventNumbers", {})->setComment("vector of bad events to filter");
0298
0299 descriptions.add("hcallaserevent", desc);
0300 }
0301
0302 std::vector<int> HcalLaserEventFilter::GetCMSSWVersion(std::string const& instring) const {
0303 std::vector<int> temp;
0304
0305 std::string prefix;
0306 std::string v1, v2, v3;
0307
0308 std::istringstream oss(instring);
0309 getline(oss, prefix, '_');
0310 getline(oss, v1, '_');
0311 getline(oss, v2, '_');
0312 getline(oss, v3, '_');
0313
0314 std::stringstream buffer(v1);
0315 int t;
0316 buffer >> t;
0317 temp.push_back(t);
0318 buffer.str();
0319 buffer << v2;
0320 buffer >> t;
0321 temp.push_back(t);
0322 buffer.str();
0323 buffer << v3;
0324 buffer >> t;
0325 temp.push_back(t);
0326
0327
0328
0329 return temp;
0330 }
0331
0332 bool HcalLaserEventFilter::IsGreaterThanMinCMSSWVersion(std::vector<int> const& currentVersion) const {
0333
0334
0335
0336
0337
0338 if (currentVersion[0] < minVersion_[0])
0339 return false;
0340 if (currentVersion[0] > minVersion_[0])
0341 return true;
0342
0343
0344
0345 if (currentVersion[1] < minVersion_[1])
0346 return false;
0347 if (currentVersion[1] > minVersion_[1])
0348 return true;
0349
0350
0351 if (currentVersion[2] < minVersion_[2])
0352 return false;
0353 if (currentVersion[2] > minVersion_[2])
0354 return true;
0355
0356 return true;
0357 }
0358
0359
0360 DEFINE_FWK_MODULE(HcalLaserEventFilter);