File indexing completed on 2024-04-06 12:10:42
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDFilter.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030
0031 #include "FWCore/Framework/interface/ESHandle.h"
0032
0033 #include <cmath>
0034 #include <iostream>
0035 #include <sstream>
0036 #include <fstream>
0037
0038 #include "FWCore/Framework/interface/Run.h"
0039 #include "FWCore/Framework/interface/LuminosityBlock.h"
0040
0041 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0042
0043 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0044 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
0045
0046
0047
0048
0049
0050 class HcalLaserHBHEHFFilter2012 : public edm::one::EDFilter<> {
0051 public:
0052 explicit HcalLaserHBHEHFFilter2012(const edm::ParameterSet&);
0053 ~HcalLaserHBHEHFFilter2012() override;
0054
0055 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0056
0057 private:
0058 bool filter(edm::Event&, const edm::EventSetup&) override;
0059 void endJob() override;
0060
0061
0062
0063 bool filterHBHE_;
0064
0065 int minCalibChannelsHBHELaser_;
0066
0067
0068
0069
0070 double minFracDiffHBHELaser_;
0071
0072
0073 double HBHEcalibThreshold_;
0074
0075 std::vector<int> CalibTS_;
0076
0077
0078 bool filterHF_;
0079
0080 int minCalibChannelsHFLaser_;
0081
0082 edm::InputTag digiLabel_;
0083 edm::EDGetTokenT<HcalCalibDigiCollection> tok_calib_;
0084 edm::EDGetTokenT<HBHEDigiCollection> tok_hbhe_;
0085
0086
0087 bool verbose_;
0088
0089 std::string prefix_;
0090 bool WriteBadToFile_;
0091 bool forceFilterTrue_;
0092 std::ofstream outfile_;
0093 };
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106 HcalLaserHBHEHFFilter2012::HcalLaserHBHEHFFilter2012(const edm::ParameterSet& ps) {
0107
0108 filterHBHE_ = ps.getParameter<bool>("filterHBHE");
0109 minCalibChannelsHBHELaser_ = ps.getParameter<int>("minCalibChannelsHBHELaser");
0110 minFracDiffHBHELaser_ = ps.getParameter<double>("minFracDiffHBHELaser");
0111 HBHEcalibThreshold_ = ps.getParameter<double>("HBHEcalibThreshold");
0112 CalibTS_ = ps.getParameter<std::vector<int> >("CalibTS");
0113
0114 filterHF_ = ps.getParameter<bool>("filterHF");
0115 minCalibChannelsHFLaser_ = ps.getParameter<int>("minCalibChannelsHFLaser");
0116
0117 digiLabel_ = ps.getParameter<edm::InputTag>("digiLabel");
0118 tok_calib_ = consumes<HcalCalibDigiCollection>(digiLabel_);
0119 tok_hbhe_ = consumes<HBHEDigiCollection>(digiLabel_);
0120
0121 verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
0122 prefix_ = ps.getUntrackedParameter<std::string>("prefix", "");
0123 WriteBadToFile_ = ps.getUntrackedParameter<bool>("WriteBadToFile", false);
0124 forceFilterTrue_ = ps.getUntrackedParameter<bool>("forceFilterTrue", false);
0125 if (WriteBadToFile_)
0126 outfile_.open("badHcalLaserList_hcalfilter.txt");
0127 }
0128
0129 HcalLaserHBHEHFFilter2012::~HcalLaserHBHEHFFilter2012() {
0130
0131
0132 }
0133
0134
0135
0136
0137
0138
0139 bool HcalLaserHBHEHFFilter2012::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0140
0141
0142 edm::Handle<HcalCalibDigiCollection> calib_digi;
0143 if (!(iEvent.getByToken(tok_calib_, calib_digi))) {
0144 edm::LogWarning("HcalLaserFilter2012") << digiLabel_ << " calib_digi not available";
0145 return true;
0146 }
0147
0148 if (!(calib_digi.isValid())) {
0149 edm::LogWarning("HcalLaserFilter2012") << digiLabel_ << " calib_digi is not valid";
0150 return true;
0151 }
0152
0153 edm::Handle<HBHEDigiCollection> hbhe_digi;
0154 if (!(iEvent.getByToken(tok_hbhe_, hbhe_digi))) {
0155 edm::LogWarning("HcalLaserFilter2012") << digiLabel_ << " hbhe_digi not available";
0156 return true;
0157 }
0158
0159 if (!(hbhe_digi.isValid())) {
0160 edm::LogWarning("HcalLaserHBHEHFFilter2012") << digiLabel_ << " hbhe_digi is not valid";
0161 return true;
0162 }
0163
0164
0165 int ncalibHBHE = 0;
0166 int ncalibHF = 0;
0167
0168 for (HcalCalibDigiCollection::const_iterator Calibiter = calib_digi->begin(); Calibiter != calib_digi->end();
0169 ++Calibiter) {
0170 const HcalCalibDataFrame digi = (const HcalCalibDataFrame)(*Calibiter);
0171 if (digi.zsMarkAndPass())
0172 continue;
0173 HcalCalibDetId myid = (HcalCalibDetId)digi.id();
0174 if (filterHBHE_ && (myid.hcalSubdet() == HcalBarrel || myid.hcalSubdet() == HcalEndcap)) {
0175 if (myid.calibFlavor() == HcalCalibDetId::HOCrosstalk)
0176 continue;
0177
0178
0179 double thischarge = 0;
0180 for (unsigned int i = 0; i < CalibTS_.size(); ++i) {
0181 thischarge += digi[CalibTS_[i]].nominal_fC();
0182 if (thischarge > HBHEcalibThreshold_) {
0183 ++ncalibHBHE;
0184 break;
0185 }
0186 }
0187
0188 if (ncalibHBHE >= minCalibChannelsHBHELaser_) {
0189 if (verbose_)
0190 std::cout << prefix_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event()
0191 << std::endl;
0192 if (WriteBadToFile_)
0193 outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
0194 if (forceFilterTrue_)
0195 return true;
0196 else
0197 return false;
0198 }
0199
0200 } else if (filterHF_ && (myid.hcalSubdet() == HcalForward)) {
0201 ++ncalibHF;
0202 if (ncalibHF >= minCalibChannelsHFLaser_) {
0203 if (verbose_)
0204 std::cout << prefix_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event()
0205 << std::endl;
0206 if (WriteBadToFile_)
0207 outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
0208 if (forceFilterTrue_)
0209 return true;
0210 else
0211 return false;
0212 }
0213 }
0214 }
0215
0216 if (filterHBHE_) {
0217
0218
0219 double badrbxfrac = 0.;
0220 double goodrbxfrac = 0.;
0221 int Nbad = 72 * 3;
0222 int Ngood = 2592 * 2 - Nbad;
0223 for (HBHEDigiCollection::const_iterator hbhe = hbhe_digi->begin(); hbhe != hbhe_digi->end(); ++hbhe) {
0224 const HBHEDataFrame digi = (const HBHEDataFrame)(*hbhe);
0225 HcalDetId myid = (HcalDetId)digi.id();
0226 bool isbad = false;
0227 if (myid.subdet() == HcalBarrel && myid.ieta() < 0) {
0228 if (myid.iphi() >= 15 && myid.iphi() <= 18)
0229 isbad = true;
0230 else if (myid.iphi() >= 27 && myid.iphi() <= 34)
0231 isbad = true;
0232 }
0233 if (isbad == true)
0234 badrbxfrac += 1.;
0235 else
0236 goodrbxfrac += 1.;
0237 }
0238 goodrbxfrac /= Ngood;
0239 badrbxfrac /= Nbad;
0240 if (goodrbxfrac - badrbxfrac > minFracDiffHBHELaser_) {
0241 if (verbose_)
0242 std::cout << prefix_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event()
0243 << std::endl;
0244 if (WriteBadToFile_)
0245 outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
0246
0247 if (forceFilterTrue_)
0248 return true;
0249 else
0250 return false;
0251 }
0252 }
0253
0254 return true;
0255
0256 }
0257
0258
0259 void HcalLaserHBHEHFFilter2012::endJob() {
0260 if (WriteBadToFile_)
0261 outfile_.close();
0262 }
0263
0264
0265 void HcalLaserHBHEHFFilter2012::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0266
0267
0268 edm::ParameterSetDescription desc;
0269 desc.setUnknown();
0270 descriptions.addDefault(desc);
0271 }
0272
0273 DEFINE_FWK_MODULE(HcalLaserHBHEHFFilter2012);