File indexing completed on 2023-03-17 10:59:52
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 HcalLaserHFFilter2012 : public edm::one::EDFilter<> {
0051 public:
0052 explicit HcalLaserHFFilter2012(const edm::ParameterSet&);
0053 ~HcalLaserHFFilter2012() 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 bool verbose_;
0063 std::string
0064 prefix_;
0065 int minCalibChannelsHFLaser_;
0066 edm::InputTag digiLabel_;
0067 edm::EDGetTokenT<HcalCalibDigiCollection> tok_calib_;
0068
0069 bool WriteBadToFile_;
0070 bool forceFilterTrue_;
0071 std::ofstream outfile_;
0072 };
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085 HcalLaserHFFilter2012::HcalLaserHFFilter2012(const edm::ParameterSet& ps) {
0086
0087 verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
0088 prefix_ = ps.getUntrackedParameter<std::string>("prefix", "");
0089 minCalibChannelsHFLaser_ = ps.getUntrackedParameter<int>("minCalibChannelsHFLaser", 10);
0090 edm::InputTag digi_default("hcalDigis");
0091 digiLabel_ = ps.getUntrackedParameter<edm::InputTag>("digiLabel", digi_default);
0092 tok_calib_ = consumes<HcalCalibDigiCollection>(digiLabel_);
0093 WriteBadToFile_ = ps.getUntrackedParameter<bool>("WriteBadToFile", false);
0094 if (WriteBadToFile_)
0095 outfile_.open("badHcalLaserList_hffilter.txt");
0096 forceFilterTrue_ = ps.getUntrackedParameter<bool>("forceFilterTrue", false);
0097
0098 }
0099
0100 HcalLaserHFFilter2012::~HcalLaserHFFilter2012() {
0101
0102
0103 }
0104
0105
0106
0107
0108
0109
0110 bool HcalLaserHFFilter2012::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0111
0112
0113 edm::Handle<HcalCalibDigiCollection> calib_digi;
0114 if (!(iEvent.getByToken(tok_calib_, calib_digi))) {
0115 edm::LogWarning("HcalLaserHFFilter2012") << digiLabel_ << " calib_digi not available";
0116 return true;
0117 }
0118
0119 if (!(calib_digi.isValid())) {
0120 edm::LogWarning("HcalLaserHFFilter2012") << digiLabel_ << " calib_digi is not valid";
0121 return true;
0122 }
0123
0124
0125 int ncalibHF = 0;
0126
0127 for (HcalCalibDigiCollection::const_iterator Calibiter = calib_digi->begin(); Calibiter != calib_digi->end();
0128 ++Calibiter) {
0129 const HcalCalibDataFrame digi = (const HcalCalibDataFrame)(*Calibiter);
0130 if (digi.zsMarkAndPass())
0131 continue;
0132 HcalCalibDetId myid = (HcalCalibDetId)digi.id();
0133 if (myid.hcalSubdet() != HcalForward)
0134 continue;
0135 ++ncalibHF;
0136 if (ncalibHF >= minCalibChannelsHFLaser_) {
0137 if (verbose_)
0138 std::cout << prefix_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event()
0139 << std::endl;
0140 if (WriteBadToFile_)
0141 outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
0142 if (forceFilterTrue_)
0143 return true;
0144 else
0145 return false;
0146 }
0147 }
0148
0149 return true;
0150 }
0151
0152
0153 void HcalLaserHFFilter2012::endJob() {
0154 if (WriteBadToFile_)
0155 outfile_.close();
0156 }
0157
0158
0159 void HcalLaserHFFilter2012::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0160
0161
0162 edm::ParameterSetDescription desc;
0163 desc.setUnknown();
0164 descriptions.addDefault(desc);
0165 }
0166
0167 DEFINE_FWK_MODULE(HcalLaserHFFilter2012);