File indexing completed on 2024-04-06 12:25:50
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <iostream>
0021 #include <memory>
0022
0023
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031
0032 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0033 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
0034 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
0035
0036
0037
0038
0039 class HcalSevLvlAnalyzer : public edm::one::EDAnalyzer<> {
0040 public:
0041 explicit HcalSevLvlAnalyzer(const edm::ParameterSet&);
0042 ~HcalSevLvlAnalyzer();
0043
0044 private:
0045 void beginJob() override;
0046 void analyze(const edm::Event&, const edm::EventSetup&) override;
0047 void endJob() override;
0048
0049
0050 const edm::ESGetToken<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd> tokSev_;
0051 };
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064 HcalSevLvlAnalyzer::HcalSevLvlAnalyzer(const edm::ParameterSet& iConfig) : tokSev_(esConsumes()) {
0065
0066
0067
0068 }
0069
0070 HcalSevLvlAnalyzer::~HcalSevLvlAnalyzer() {
0071
0072
0073 }
0074
0075
0076
0077
0078
0079
0080 void HcalSevLvlAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0081
0082 const HcalSeverityLevelComputer* myProd = &iSetup.getData(tokSev_);
0083
0084
0085 std::cout << std::showbase;
0086
0087 HcalGenericDetId myIdHB(1107313548);
0088 uint32_t sampleHBChSt[5] = {0,
0089 4,
0090 64,
0091 1,
0092 32};
0093 uint32_t sampleHBRHFlag[4] = {0,
0094 32,
0095 2,
0096 1};
0097
0098 for (unsigned i = 0; i < 4; i++)
0099 for (unsigned k = 0; k < 5; k++)
0100 {
0101 int theLevel = myProd->getSeverityLevel(myIdHB, sampleHBRHFlag[i], sampleHBChSt[k]);
0102
0103 std::cout << "Status for " << myIdHB << " with RHFlag " << sampleHBRHFlag[i] << " (" << std::hex
0104 << sampleHBRHFlag[i] << ") and ChStFlag " << std::dec << sampleHBChSt[k] << " (" << std::hex
0105 << sampleHBChSt[k] << std::dec << ") is: " << theLevel;
0106 std::cout << std::endl;
0107
0108 bool dropchannel = myProd->dropChannel(sampleHBChSt[k]);
0109 bool recovered = myProd->recoveredRecHit(myIdHB, sampleHBRHFlag[i]);
0110
0111 std::cout << "DropChannel status for " << myIdHB << " with RHFlag " << sampleHBRHFlag[i] << " (" << std::hex
0112 << sampleHBRHFlag[i] << ") and ChStFlag " << std::dec << sampleHBChSt[k] << " (" << std::hex
0113 << sampleHBChSt[k] << std::dec << ") is: " << dropchannel << ", recovered status is: " << recovered
0114 << std::endl;
0115 std::cout << std::endl;
0116 }
0117 std::cout << std::endl;
0118
0119
0120
0121
0122
0123
0124
0125
0126
0127 HcalGenericDetId myIdHE(1140869515);
0128 uint32_t sampleHEChSt[5] = {0,
0129 4,
0130 64,
0131 1,
0132 32};
0133 uint32_t sampleHERHFlag[4] = {0,
0134 32,
0135 2,
0136 1};
0137
0138 for (unsigned i = 0; i < 4; i++)
0139 for (unsigned k = 0; k < 5; k++)
0140 {
0141 int theLevel = myProd->getSeverityLevel(myIdHE, sampleHERHFlag[i], sampleHEChSt[k]);
0142
0143 std::cout << "Status for " << myIdHE << " with RHFlag " << sampleHERHFlag[i] << " (" << std::hex
0144 << sampleHERHFlag[i] << ") and ChStFlag " << std::dec << sampleHEChSt[k] << " (" << std::hex
0145 << sampleHEChSt[k] << std::dec << ") is: " << theLevel;
0146 std::cout << std::endl;
0147 }
0148 std::cout << std::endl;
0149
0150
0151 HcalGenericDetId myIdHF(1207979653);
0152 uint32_t sampleHFChSt[5] = {0,
0153 4,
0154 64,
0155 1,
0156 32};
0157 uint32_t sampleHFRHFlag[4] = {0,
0158 32,
0159 1,
0160 2};
0161
0162 for (unsigned i = 0; i < 4; i++)
0163 for (unsigned k = 0; k < 5; k++)
0164 {
0165 int theLevel = myProd->getSeverityLevel(myIdHF, sampleHFRHFlag[i], sampleHFChSt[k]);
0166
0167 std::cout << "Status for " << myIdHF << " with RHFlag " << sampleHFRHFlag[i] << " (" << std::hex
0168 << sampleHFRHFlag[i] << ") and ChStFlag " << std::dec << sampleHFChSt[k] << " (" << std::hex
0169 << sampleHFChSt[k] << std::dec << ") is: " << theLevel;
0170 std::cout << std::endl;
0171 }
0172 std::cout << std::endl;
0173
0174
0175 HcalGenericDetId myIdHO(1174471682);
0176 uint32_t sampleHOChSt[5] = {0,
0177 4,
0178 64,
0179 1,
0180 32};
0181 uint32_t sampleHORHFlag[3] = {0,
0182 32,
0183 1};
0184
0185 for (unsigned i = 0; i < 3; i++)
0186 for (unsigned k = 0; k < 5; k++)
0187 {
0188 int theLevel = myProd->getSeverityLevel(myIdHO, sampleHORHFlag[i], sampleHOChSt[k]);
0189
0190 std::cout << "Status for " << myIdHO << " with RHFlag " << sampleHORHFlag[i] << " (" << std::hex
0191 << sampleHORHFlag[i] << ") and ChStFlag " << std::dec << sampleHOChSt[k] << " (" << std::hex
0192 << sampleHOChSt[k] << std::dec << ") is: " << theLevel;
0193 std::cout << std::endl;
0194 }
0195 std::cout << std::endl;
0196
0197
0198
0199
0200
0201 std::cout << std::noshowbase;
0202 }
0203
0204
0205 void HcalSevLvlAnalyzer::beginJob() {}
0206
0207
0208 void HcalSevLvlAnalyzer::endJob() {}
0209
0210
0211 DEFINE_FWK_MODULE(HcalSevLvlAnalyzer);