File indexing completed on 2024-04-06 11:59:07
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0010
0011 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
0012 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0013 #include "DataFormats/HcalDigi/interface/HcalUMNioDigi.h"
0014 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0015 #include "DataFormats/HcalDigi/interface/HcalCalibrationEventTypes.h"
0016 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0017 #include "Geometry/HcalCommonData/interface/HcalBadLaserChannels.h"
0018 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
0019
0020 #include <array>
0021 #include <atomic>
0022 #include <iostream>
0023 #include <memory>
0024 #include <string>
0025 #include <TH1F.h>
0026
0027 class HcalLaserTest : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0028 public:
0029 explicit HcalLaserTest(const edm::ParameterSet&);
0030 ~HcalLaserTest() override;
0031 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0032
0033 private:
0034 void analyze(edm::Event const&, edm::EventSetup const&) override;
0035 void beginJob() override;
0036 void beginRun(edm::Run const&, edm::EventSetup const&) override {}
0037 void endRun(edm::Run const&, edm::EventSetup const&) override {}
0038 void endJob(void) override;
0039
0040
0041 const edm::EDGetTokenT<HBHEDigiCollection> inputTokenHBHE_;
0042 const edm::EDGetTokenT<QIE10DigiCollection> inputTokenHF_;
0043 const edm::EDGetTokenT<HcalUMNioDigi> ioTokenUMN_;
0044 const double minFracDiffHBHELaser_, minFracHFLaser_;
0045 const int minADCHBHE_, minADCHF_;
0046
0047 const bool testMode_;
0048 mutable std::array<std::atomic<int>, 16> eventsByType_;
0049 mutable std::array<std::atomic<int>, 16> passedEventsByType_;
0050 TH1F *h_hb1_, *h_hb2_, *h_hb3_, *h_hb4_;
0051 TH1F *h_hb5_, *h_hf1_, *h_hf2_;
0052 };
0053
0054 HcalLaserTest::HcalLaserTest(const edm::ParameterSet& config)
0055 : inputTokenHBHE_(consumes<HBHEDigiCollection>(config.getParameter<edm::InputTag>("InputHBHE"))),
0056 inputTokenHF_(consumes<QIE10DigiCollection>(config.getParameter<edm::InputTag>("InputHF"))),
0057 ioTokenUMN_(consumes<HcalUMNioDigi>(config.getParameter<edm::InputTag>("UMNioDigis"))),
0058 minFracDiffHBHELaser_(config.getParameter<double>("minFracDiffHBHELaser")),
0059 minFracHFLaser_(config.getParameter<double>("minFracHFLaser")),
0060 minADCHBHE_(config.getParameter<int>("minADCHBHE")),
0061 minADCHF_(config.getParameter<int>("minADCHF")),
0062 testMode_(config.getUntrackedParameter<bool>("testMode", false)),
0063 eventsByType_(),
0064 passedEventsByType_() {
0065 usesResource(TFileService::kSharedResource);
0066 }
0067
0068 HcalLaserTest::~HcalLaserTest() {}
0069
0070 void HcalLaserTest::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0071 edm::ParameterSetDescription desc;
0072 desc.add<edm::InputTag>("InputHBHE", edm::InputTag("source"));
0073 desc.add<edm::InputTag>("InputHF", edm::InputTag("source"));
0074 desc.add<edm::InputTag>("UMNioDigis", edm::InputTag("UMNioDigis"));
0075 desc.add<int>("minADCHBHE", 10);
0076 desc.add<int>("minADCHF", 10);
0077 desc.add<double>("minFracDiffHBHELaser", 0.3);
0078 desc.add<double>("minFracHFLaser", 0.3);
0079 desc.addUntracked<bool>("testMode", false);
0080 descriptions.add("hcalLaserTest", desc);
0081 }
0082
0083 void HcalLaserTest::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0084 const edm::Handle<HBHEDigiCollection>& hbhe_digi = iEvent.getHandle(inputTokenHBHE_);
0085
0086 const edm::Handle<QIE10DigiCollection>& hf_digi = iEvent.getHandle(inputTokenHF_);
0087
0088
0089 double badrbxfracHBHE(0), goodrbxfracHBHE(0), rbxfracHF(0);
0090 int NbadHBHE = HcalBadLaserChannels::badChannelsHBHE();
0091 int NgoodHBHE = 2592 * 2 - NbadHBHE;
0092 int NallHF = 864 * 4;
0093 int eventType = 0;
0094 int laserType = 0;
0095
0096
0097
0098
0099
0100 const edm::Handle<HcalUMNioDigi>& cumn = iEvent.getHandle(ioTokenUMN_);
0101
0102 eventType = static_cast<int>(cumn->eventType());
0103 laserType = static_cast<int>(cumn->valueUserWord(0));
0104
0105 if (eventType == 14)
0106 eventsByType_.at(laserType)++;
0107
0108 if (testMode_)
0109 edm::LogVerbatim("HcalLaserTest") << "hbhe digi collection size: " << hbhe_digi->size() << "\n"
0110 << "hf digi collection size: " << hf_digi->size() << "\n"
0111 << "Event type: " << eventType << " Laser type: " << laserType << std::endl;
0112
0113 for (auto& hbhe : *hbhe_digi) {
0114 const HBHEDataFrame digi = (const HBHEDataFrame)(hbhe);
0115
0116 int maxdigiHB(0);
0117 for (int i = 0; i < digi.size(); i++)
0118 if (digi.sample(i).adc() > maxdigiHB)
0119 maxdigiHB = digi.sample(i).adc();
0120 bool passCut = (maxdigiHB > minADCHBHE_);
0121
0122 HcalDetId myid = (HcalDetId)digi.id();
0123 bool isbad = HcalBadLaserChannels::badChannelHBHE(myid);
0124
0125 if (isbad)
0126 h_hb2_->Fill(maxdigiHB);
0127 else
0128 h_hb1_->Fill(maxdigiHB);
0129
0130 if (passCut) {
0131 if (isbad)
0132 badrbxfracHBHE += 1.;
0133 else
0134 goodrbxfracHBHE += 1.;
0135 }
0136 }
0137 goodrbxfracHBHE /= NgoodHBHE;
0138 badrbxfracHBHE /= NbadHBHE;
0139 h_hb3_->Fill(goodrbxfracHBHE);
0140 h_hb4_->Fill(badrbxfracHBHE);
0141 h_hb5_->Fill(goodrbxfracHBHE - badrbxfracHBHE);
0142
0143 for (auto& hf : *(hf_digi)) {
0144 const QIE10DataFrame digi = (const QIE10DataFrame)(hf);
0145 bool passCut(false);
0146 int maxdigiHF(0);
0147 for (int i = 0; i < digi.samples(); i++)
0148 if (digi[i].adc() > maxdigiHF)
0149 maxdigiHF = digi[i].adc();
0150 if (maxdigiHF > minADCHF_)
0151 passCut = true;
0152 h_hf1_->Fill(maxdigiHF);
0153
0154 if (passCut)
0155 rbxfracHF += 1.;
0156 }
0157 rbxfracHF /= NallHF;
0158 h_hf2_->Fill(rbxfracHF);
0159
0160 if (testMode_)
0161 edm::LogWarning("HcalLaserTest") << "******************************************************************\n"
0162 << "goodrbxfracHBHE: " << goodrbxfracHBHE << " badrbxfracHBHE: " << badrbxfracHBHE
0163 << " Size " << hbhe_digi->size() << "\n"
0164 << "rbxfracHF: " << rbxfracHF << " Size " << hf_digi->size()
0165 << "\n******************************************************************";
0166
0167 if (((goodrbxfracHBHE - badrbxfracHBHE) < minFracDiffHBHELaser_) || (rbxfracHF < minFracHFLaser_)) {
0168 } else {
0169 passedEventsByType_.at(laserType)++;
0170 }
0171 }
0172
0173 void HcalLaserTest::beginJob() {
0174 for (auto& i : eventsByType_)
0175 i = 0;
0176 for (auto& i : passedEventsByType_)
0177 i = 0;
0178
0179 edm::Service<TFileService> tfile;
0180 if (!tfile.isAvailable())
0181 throw cms::Exception("HcalLaserTest") << "TFileService unavailable: "
0182 << "please add it to config file";
0183 h_hb1_ = tfile->make<TH1F>("hb1", "Maximum ADC in HB (Good)", 5000, 0, 100);
0184 h_hb2_ = tfile->make<TH1F>("hb2", "Maximum ADC in HB (Bad)", 5000, 0, 100);
0185 h_hb3_ = tfile->make<TH1F>("hb3", "Signal Channel fraction (Good)", 1000, 0, 1);
0186 h_hb4_ = tfile->make<TH1F>("hb4", "Signal Channel fraction (Bad)", 1000, 0, 1);
0187 h_hb5_ = tfile->make<TH1F>("hb5", "Signal Channel fraction (Diff)", 1000, -1, 1);
0188 h_hf1_ = tfile->make<TH1F>("hf1", "Maximum ADC in HF", 5000, 0, 100);
0189 h_hf2_ = tfile->make<TH1F>("hf2", "Signal Channel fraction (HF)", 1000, 0, 1);
0190 }
0191
0192 void HcalLaserTest::endJob() {
0193 edm::LogVerbatim("HcalLaserTest") << "Summary of filter decisions (passed/total): \n"
0194 << passedEventsByType_.at(hc_Null) << "/" << eventsByType_.at(hc_Null)
0195 << "(No Calib), " << passedEventsByType_.at(hc_Pedestal) << "/"
0196 << eventsByType_.at(hc_Pedestal) << "(Pedestal), "
0197 << passedEventsByType_.at(hc_RADDAM) << "/" << eventsByType_.at(hc_RADDAM)
0198 << "(RADDAM), " << passedEventsByType_.at(hc_HBHEHPD) << "/"
0199 << eventsByType_.at(hc_HBHEHPD) << "(HBHE/HPD), "
0200 << passedEventsByType_.at(hc_HOHPD) << "/" << eventsByType_.at(hc_HOHPD)
0201 << "(HO/HPD), " << passedEventsByType_.at(hc_HFPMT) << "/"
0202 << eventsByType_.at(hc_HFPMT) << "(HF/PMT), " << passedEventsByType_.at(6) << "/"
0203 << eventsByType_.at(6) << "(ZDC), " << passedEventsByType_.at(7) << "/"
0204 << eventsByType_.at(7) << "(HEPMega)\n"
0205 << passedEventsByType_.at(8) << "/" << eventsByType_.at(8) << "(HEMMega), "
0206 << passedEventsByType_.at(9) << "/" << eventsByType_.at(9) << "(HBPMega), "
0207 << passedEventsByType_.at(10) << "/" << eventsByType_.at(10) << "(HBMMega), "
0208 << passedEventsByType_.at(11) << "/" << eventsByType_.at(11) << "(Undefined), "
0209 << passedEventsByType_.at(12) << "/" << eventsByType_.at(12) << "(CRF), "
0210 << passedEventsByType_.at(13) << "/" << eventsByType_.at(13) << "(Calib), "
0211 << passedEventsByType_.at(14) << "/" << eventsByType_.at(14) << "(Safe), "
0212 << passedEventsByType_.at(15) << "/" << eventsByType_.at(15) << "(Undefined)";
0213 }
0214
0215 #include "FWCore/Framework/interface/MakerMacros.h"
0216
0217 DEFINE_FWK_MODULE(HcalLaserTest);