File indexing completed on 2024-04-06 12:14:20
0001 #include <fstream>
0002 #include <iostream>
0003 #include <string>
0004 #include <vector>
0005
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Framework/interface/ESTransientHandle.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015
0016 #include "DetectorDescription/Core/interface/DDCompactView.h"
0017 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0018 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0019 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0020 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0021 #include "DataFormats/HcalDetId/interface/HcalCalibDetId.h"
0022 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0023
0024 class HcalDetId2DenseTester : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0025 public:
0026 explicit HcalDetId2DenseTester(const edm::ParameterSet&);
0027
0028 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029
0030 private:
0031 void analyze(edm::Event const&, edm::EventSetup const&) override;
0032 void beginJob() override {}
0033 void beginRun(edm::Run const&, edm::EventSetup const&) override {}
0034 void endRun(edm::Run const&, edm::EventSetup const&) override {}
0035 void doTestFile(const HcalTopology& topology);
0036 void doTestHcalDetId(const HcalTopology& topology);
0037 void doTestHcalCalibDetId(const HcalTopology& topology);
0038 void doTestOnlyHcalCalibDetId(const HcalTopology& topology);
0039 void doTestOneCell(const HcalTopology&, const std::string&, const HcalCalibDetId&, int&, int&, int&, int&);
0040 std::vector<std::string> splitString(const std::string& fLine);
0041
0042 const std::string fileName_;
0043 const bool testCalib_;
0044 const edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tokTopo_;
0045 };
0046
0047 HcalDetId2DenseTester::HcalDetId2DenseTester(const edm::ParameterSet& iC)
0048 : fileName_(iC.getUntrackedParameter<std::string>("fileName", "")),
0049 testCalib_(iC.getUntrackedParameter<bool>("testCalib", false)),
0050 tokTopo_{esConsumes<HcalTopology, HcalRecNumberingRecord>(edm::ESInputTag{})} {}
0051
0052 void HcalDetId2DenseTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0053 edm::ParameterSetDescription desc;
0054 desc.addUntracked<std::string>("fileName", "");
0055 desc.addUntracked<bool>("testCalib", false);
0056 descriptions.add("hcalDetId2DenseTester", desc);
0057 }
0058
0059 void HcalDetId2DenseTester::analyze(edm::Event const&, edm::EventSetup const& iSetup) {
0060 const auto& topo = iSetup.getData(tokTopo_);
0061 if (!testCalib_) {
0062 doTestFile(topo);
0063 doTestHcalDetId(topo);
0064 doTestHcalCalibDetId(topo);
0065 } else {
0066 doTestOnlyHcalCalibDetId(topo);
0067 }
0068 }
0069
0070 void HcalDetId2DenseTester::doTestFile(const HcalTopology& topology) {
0071 if (!fileName_.empty()) {
0072 std::ifstream fInput(fileName_);
0073 if (!fInput.good()) {
0074 edm::LogVerbatim("HCalGeom") << "Cannot open file " << fileName_;
0075 } else {
0076 char buffer[1024];
0077 unsigned int all(0), total(0), good(0), ok(0), bad(0);
0078 while (fInput.getline(buffer, 1024)) {
0079 ++all;
0080 if (buffer[0] == '#')
0081 continue;
0082 std::vector<std::string> items = splitString(std::string(buffer));
0083 if (items.size() != 4) {
0084 edm::LogVerbatim("HCalGeom") << "Ignore line: " << buffer;
0085 } else {
0086 ++total;
0087 int ieta = std::atoi(items[1].c_str());
0088 ;
0089 int iphi = std::atoi(items[2].c_str());
0090 ;
0091 std::string error;
0092 if ((items[0] == "HB") || (items[0] == "HE") || (items[0] == "HF") || (items[0] == "HO")) {
0093 HcalSubdetector sd(HcalBarrel);
0094 if (items[0] == "HE") {
0095 sd = HcalEndcap;
0096 } else if (items[0] == "HF") {
0097 sd = HcalForward;
0098 } else if (items[0] == "HO") {
0099 sd = HcalOuter;
0100 }
0101 int depth = std::atoi(items[3].c_str());
0102 HcalDetId cell(sd, ieta, iphi, depth);
0103 if (topology.valid(cell)) {
0104 ++ok;
0105 unsigned int dense = topology.detId2denseId(DetId(cell));
0106 DetId id = topology.denseId2detId(dense);
0107 if (cell == HcalDetId(id)) {
0108 error = "";
0109 ++good;
0110 } else {
0111 error = "***ERROR***";
0112 ++bad;
0113 }
0114 edm::LogVerbatim("HCalGeom") << total << " " << cell << " Dense Index " << dense << " gives back "
0115 << HcalDetId(id) << " " << error;
0116 }
0117
0118 } else if (items[0].find("CALIB_") == 0) {
0119 HcalSubdetector sd = HcalOther;
0120 if (items[0].find("HB") != std::string::npos)
0121 sd = HcalBarrel;
0122 else if (items[0].find("HE") != std::string::npos)
0123 sd = HcalEndcap;
0124 else if (items[0].find("HO") != std::string::npos)
0125 sd = HcalOuter;
0126 else if (items[0].find("HF") != std::string::npos)
0127 sd = HcalForward;
0128 int channel = std::atoi(items[3].c_str());
0129 HcalCalibDetId cell(sd, ieta, iphi, channel);
0130 if (topology.validCalib(cell)) {
0131 ++ok;
0132 unsigned int dense = topology.detId2denseIdCALIB(DetId(cell));
0133 DetId id = topology.denseId2detIdCALIB(dense);
0134 if (cell == HcalCalibDetId(id)) {
0135 error = "";
0136 ++good;
0137 } else {
0138 error = "***ERROR***";
0139 ++bad;
0140 }
0141 edm::LogVerbatim("HCalGeom") << total << " " << cell << " Dense Index " << dense << " gives back "
0142 << HcalCalibDetId(id) << " " << error;
0143 }
0144 } else if ((items[0] == "HOX") || (items[0] == "HBX") || (items[0] == "HEX")) {
0145 HcalCalibDetId cell =
0146 ((items[0] == "HOX") ? (HcalCalibDetId(HcalCalibDetId::HOCrosstalk, ieta, iphi))
0147 : ((items[0] == "HBX") ? (HcalCalibDetId(HcalCalibDetId::HBX, ieta, iphi))
0148 : (HcalCalibDetId(HcalCalibDetId::HEX, ieta, iphi))));
0149 if (topology.validCalib(cell)) {
0150 ++ok;
0151 unsigned int dense = topology.detId2denseIdCALIB(DetId(cell));
0152 DetId id = topology.denseId2detIdCALIB(dense);
0153 if (cell == HcalCalibDetId(id)) {
0154 error = "";
0155 ++good;
0156 } else {
0157 error = "***ERROR***";
0158 ++bad;
0159 }
0160 edm::LogVerbatim("HCalGeom") << good << " " << cell << " Dense Index " << dense << " gives back "
0161 << HcalCalibDetId(id) << " " << error;
0162 }
0163 }
0164 }
0165 }
0166 fInput.close();
0167 edm::LogVerbatim("HCalGeom") << "Reads total of " << all << ":" << total << " records with " << good << ":" << ok
0168 << " good and " << bad << " bad DetId's";
0169 }
0170 }
0171 }
0172
0173 void HcalDetId2DenseTester::doTestHcalDetId(const HcalTopology& topology) {
0174 const int n = 48;
0175 HcalSubdetector sds[n] = {HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel,
0176 HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalEndcap, HcalEndcap,
0177 HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap,
0178 HcalEndcap, HcalEndcap, HcalEndcap, HcalOuter, HcalOuter, HcalOuter, HcalOuter,
0179 HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalOuter,
0180 HcalOuter, HcalForward, HcalForward, HcalForward, HcalForward, HcalForward, HcalForward,
0181 HcalForward, HcalForward, HcalForward, HcalForward, HcalForward, HcalForward};
0182 int ietas[n] = {1, 5, 9, 12, 14, 16, -2, -4, -8, -11, -13, -16, 18, 21, 23, 25, 27, 29, -19, -20, -22, -24, -26, -28,
0183 1, 3, 5, 9, 11, 13, -2, -4, -6, -10, -12, -14, 30, 32, 34, 36, 38, 40, -31, -33, -35, -37, -39, -41};
0184 int iphis[n] = {11, 21, 31, 41, 51, 61, 5, 15, 25, 35, 45, 65, 11, 23, 35, 47, 59, 71, 3, 15, 27, 39, 51, 63,
0185 11, 21, 31, 41, 51, 61, 5, 15, 25, 35, 45, 65, 11, 23, 35, 47, 59, 71, 3, 15, 27, 39, 51, 63};
0186 int depth[n] = {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2,
0187 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2};
0188
0189
0190 edm::LogVerbatim("HCalGeom") << "\nCheck on Dense Index for DetId's"
0191 << "\n==================================";
0192 int total(0), good(0), bad(0);
0193 std::string error;
0194 for (int i = 0; i < n; ++i) {
0195 HcalDetId cell(sds[i], ietas[i], iphis[i], depth[i]);
0196 if (topology.valid(cell)) {
0197 ++total;
0198 unsigned int dense = topology.detId2denseId(DetId(cell));
0199 DetId id = topology.denseId2detId(dense);
0200 if (cell == HcalDetId(id)) {
0201 ++good;
0202 error = "";
0203 } else {
0204 ++bad;
0205 error = "**** ERROR *****";
0206 }
0207 edm::LogVerbatim("HCalGeom") << "[" << total << "] " << cell << " Dense " << dense << " o/p " << HcalDetId(id)
0208 << " " << error;
0209 }
0210 }
0211
0212 edm::LogVerbatim("HCalGeom") << "Analyzes total of " << n << ":" << total << " HcalDetIds with " << good
0213 << " good and " << bad << " bad DetId's";
0214 }
0215
0216 void HcalDetId2DenseTester::doTestHcalCalibDetId(const HcalTopology& topology) {
0217 const int n = 48;
0218 HcalCalibDetId::CalibDetType dets[n] = {HcalCalibDetId::CalibrationBox,
0219 HcalCalibDetId::CalibrationBox,
0220 HcalCalibDetId::CalibrationBox,
0221 HcalCalibDetId::CalibrationBox,
0222 HcalCalibDetId::CalibrationBox,
0223 HcalCalibDetId::CalibrationBox,
0224 HcalCalibDetId::CalibrationBox,
0225 HcalCalibDetId::CalibrationBox,
0226 HcalCalibDetId::CalibrationBox,
0227 HcalCalibDetId::CalibrationBox,
0228 HcalCalibDetId::CalibrationBox,
0229 HcalCalibDetId::CalibrationBox,
0230 HcalCalibDetId::CalibrationBox,
0231 HcalCalibDetId::CalibrationBox,
0232 HcalCalibDetId::CalibrationBox,
0233 HcalCalibDetId::CalibrationBox,
0234 HcalCalibDetId::CalibrationBox,
0235 HcalCalibDetId::CalibrationBox,
0236 HcalCalibDetId::CalibrationBox,
0237 HcalCalibDetId::CalibrationBox,
0238 HcalCalibDetId::CalibrationBox,
0239 HcalCalibDetId::CalibrationBox,
0240 HcalCalibDetId::CalibrationBox,
0241 HcalCalibDetId::CalibrationBox,
0242 HcalCalibDetId::HOCrosstalk,
0243 HcalCalibDetId::HOCrosstalk,
0244 HcalCalibDetId::HOCrosstalk,
0245 HcalCalibDetId::HOCrosstalk,
0246 HcalCalibDetId::HOCrosstalk,
0247 HcalCalibDetId::HOCrosstalk,
0248 HcalCalibDetId::HOCrosstalk,
0249 HcalCalibDetId::HOCrosstalk,
0250 HcalCalibDetId::HBX,
0251 HcalCalibDetId::HBX,
0252 HcalCalibDetId::HBX,
0253 HcalCalibDetId::HBX,
0254 HcalCalibDetId::HBX,
0255 HcalCalibDetId::HBX,
0256 HcalCalibDetId::HBX,
0257 HcalCalibDetId::HBX,
0258 HcalCalibDetId::HEX,
0259 HcalCalibDetId::HEX,
0260 HcalCalibDetId::HEX,
0261 HcalCalibDetId::HEX,
0262 HcalCalibDetId::HEX,
0263 HcalCalibDetId::HEX,
0264 HcalCalibDetId::HEX,
0265 HcalCalibDetId::HEX};
0266 HcalSubdetector sds[n] = {HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel,
0267 HcalBarrel, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap, HcalOuter, HcalOuter,
0268 HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalForward, HcalForward, HcalForward,
0269 HcalForward, HcalForward, HcalForward, HcalBarrel, HcalBarrel, HcalBarrel, HcalBarrel,
0270 HcalBarrel, HcalBarrel, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap, HcalEndcap,
0271 HcalEndcap, HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalOuter, HcalOuter,
0272 HcalForward, HcalForward, HcalForward, HcalForward, HcalForward, HcalForward};
0273 int ietas[n] = {1, 1, 1, -1, -1, -1, 1, 1, 1, -1, -1, -1, 1, -1, 0, 2,
0274 -1, -2, 1, 1, -1, -1, 1, -1, 4, 4, -4, -4, 15, 15, -15, -15,
0275 -16, -16, -16, -16, 16, 16, 16, 16, 25, 25, 27, 27, -25, -25, -27, -27};
0276 int iphis[n] = {11, 23, 35, 47, 59, 71, 3, 15, 27, 39, 51, 65, 22, 11, 35, 71, 47, 59, 19, 37, 55, 1, 19, 37,
0277 47, 36, 36, 41, 51, 61, 5, 15, 25, 35, 45, 65, 11, 1, 35, 47, 59, 71, 3, 15, 27, 39, 51, 63};
0278 int depth[n] = {0, 1, 2, 0, 1, 2, 0, 1, 3, 4, 5, 6, 0, 1, 0, 1, 7, 7, 0, 1, 0, 1, 0, 1,
0279 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2};
0280
0281
0282 edm::LogVerbatim("HCalGeom") << "\nCheck on Dense Index for CalibDetId's"
0283 << "\n=======================================";
0284 int total(0), good(0), bad(0);
0285 std::string error;
0286 for (int i = 0; i < n; ++i) {
0287 HcalCalibDetId cell;
0288 if (dets[i] == HcalCalibDetId::CalibrationBox) {
0289 cell = HcalCalibDetId(sds[i], ietas[i], iphis[i], depth[i]);
0290 } else {
0291 cell = HcalCalibDetId(dets[i], ietas[i], iphis[i]);
0292 }
0293 if (topology.validCalib(cell)) {
0294 ++total;
0295 unsigned int dense = topology.detId2denseIdCALIB(DetId(cell));
0296 DetId id = topology.denseId2detIdCALIB(dense);
0297 if (cell == HcalCalibDetId(id)) {
0298 ++good;
0299 error = "";
0300 } else {
0301 ++bad;
0302 error = "**** ERROR *****";
0303 }
0304 edm::LogVerbatim("HCalGeom") << "[" << total << "] " << cell << " Dense " << dense << " o/p "
0305 << HcalCalibDetId(id) << " " << error;
0306 }
0307 }
0308 edm::LogVerbatim("HCalGeom") << "Analyzes total of " << n << ":" << total << " CalibIds with " << good << " good and "
0309 << bad << " bad DetId's";
0310 }
0311
0312 void HcalDetId2DenseTester::doTestOnlyHcalCalibDetId(const HcalTopology& topology) {
0313
0314 edm::LogVerbatim("HCalGeom") << "\nCheck on Dense Index for CalibDetId's"
0315 << "\n=======================================";
0316 int allT(0), totalT(0), goodT(0), badT(0);
0317 int all(0), total(0), good(0), bad(0);
0318
0319 for (int i1 = 0; i1 < 2; ++i1) {
0320 int ieta = 2 * i1 - 1;
0321 for (int i2 = 0; i2 < 3; ++i2) {
0322 int channel = i2;
0323 for (int i3 = 0; i3 < 18; ++i3) {
0324 int iphi = 4 * i3 + 3;
0325 HcalCalibDetId cell(HcalBarrel, ieta, iphi, channel);
0326 doTestOneCell(topology, "HcalBarrel", cell, all, total, good, bad);
0327 }
0328 }
0329 }
0330 edm::LogVerbatim("HCalGeom") << "CalibHB:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0331 << " bad DetId's\n";
0332 allT += all;
0333 totalT += total;
0334 goodT += good;
0335 badT += bad;
0336 all = 0;
0337 total = 0;
0338 good = 0;
0339 bad = 0;
0340
0341
0342 for (int i1 = 0; i1 < 2; ++i1) {
0343 int ieta = 2 * i1 - 1;
0344 for (int i2 = 0; i2 < 7; ++i2) {
0345 int channel = i2;
0346 for (int i3 = 0; i3 < 18; ++i3) {
0347 int iphi = 4 * i3 + 3;
0348 HcalCalibDetId cell(HcalEndcap, ieta, iphi, channel);
0349 doTestOneCell(topology, "HcalEndcap", cell, all, total, good, bad);
0350 }
0351 }
0352 }
0353 edm::LogVerbatim("HCalGeom") << "CalibHE:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0354 << " bad DetId's\n";
0355 allT += all;
0356 totalT += total;
0357 goodT += good;
0358 badT += bad;
0359 all = 0;
0360 total = 0;
0361 good = 0;
0362 bad = 0;
0363
0364
0365 int chanHF[4] = {0, 1, 8, 9};
0366 for (int i1 = 0; i1 < 2; ++i1) {
0367 int ieta = 2 * i1 - 1;
0368 for (int i2 = 0; i2 < 4; ++i2) {
0369 int channel = chanHF[i2];
0370 int phimax = (i2 == 3) ? 1 : 4;
0371 for (int i3 = 0; i3 < phimax; ++i3) {
0372 int iphi = 18 * i3 + 3;
0373 HcalCalibDetId cell(HcalForward, ieta, iphi, channel);
0374 doTestOneCell(topology, "HcalForward", cell, all, total, good, bad);
0375 }
0376 }
0377 }
0378 edm::LogVerbatim("HCalGeom") << "CalibHF:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0379 << " bad DetId's\n";
0380 allT += all;
0381 totalT += total;
0382 goodT += good;
0383 badT += bad;
0384 all = 0;
0385 total = 0;
0386 good = 0;
0387 bad = 0;
0388
0389
0390 int chanHO[3] = {0, 1, 7};
0391 int phiHO[5] = {59, 47, 53, 47, 47};
0392 for (int i1 = 0; i1 < 5; ++i1) {
0393 int ieta = i1 - 2;
0394 for (int i2 = 0; i2 < 3; ++i2) {
0395 int channel = chanHO[i2];
0396 int phimax = (i2 == 2) ? 1 : ((ieta == 0) ? 12 : 6);
0397 for (int i3 = 0; i3 < phimax; ++i3) {
0398 int iphi = (i2 == 2) ? phiHO[i1] : ((ieta == 0) ? 6 * i3 + 5 : 12 * i3 + 11);
0399 HcalCalibDetId cell(HcalOuter, ieta, iphi, channel);
0400 doTestOneCell(topology, "HcalOuter", cell, all, total, good, bad);
0401 }
0402 }
0403 }
0404 edm::LogVerbatim("HCalGeom") << "CalibHO:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0405 << " bad DetId's\n";
0406 allT += all;
0407 totalT += total;
0408 goodT += good;
0409 badT += bad;
0410 all = 0;
0411 total = 0;
0412 good = 0;
0413 bad = 0;
0414
0415
0416 int etaHOX[4] = {4, -4, 15, -15};
0417 for (int i1 = 0; i1 < 4; ++i1) {
0418 int ieta = etaHOX[i1];
0419 int phimax = (i1 < 2) ? 36 : 72;
0420 for (int i3 = 0; i3 < phimax; ++i3) {
0421 int iphi = (i1 >= 2) ? i3 + 1 : ((i3 + 4) % 12 < 6) ? 2 * i3 + 1 : 2 * i3 + 2;
0422 HcalCalibDetId cell(HcalCalibDetId::HOCrosstalk, ieta, iphi);
0423 doTestOneCell(topology, "HOX", cell, all, total, good, bad);
0424 }
0425 }
0426 edm::LogVerbatim("HCalGeom") << "HOX:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0427 << " bad DetId's\n";
0428 allT += all;
0429 totalT += total;
0430 goodT += good;
0431 badT += bad;
0432 all = 0;
0433 total = 0;
0434 good = 0;
0435 bad = 0;
0436
0437
0438 int etaHBX[2] = {16, -16};
0439 for (int i1 = 0; i1 < 2; ++i1) {
0440 int ieta = etaHBX[i1];
0441 for (int i3 = 0; i3 < 72; ++i3) {
0442 int iphi = i3 + 1;
0443 HcalCalibDetId cell(HcalCalibDetId::HBX, ieta, iphi);
0444 doTestOneCell(topology, "HBX", cell, all, total, good, bad);
0445 }
0446 }
0447 edm::LogVerbatim("HCalGeom") << "HBX:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0448 << " bad DetId's\n";
0449 allT += all;
0450 totalT += total;
0451 goodT += good;
0452 badT += bad;
0453 all = 0;
0454 total = 0;
0455 good = 0;
0456 bad = 0;
0457
0458
0459 int etaHEX[4] = {25, -25, 27, -27};
0460 for (int i1 = 0; i1 < 4; ++i1) {
0461 int ieta = etaHEX[i1];
0462 for (int i3 = 0; i3 < 36; ++i3) {
0463 int iphi = 2 * i3 + 1;
0464 HcalCalibDetId cell(HcalCalibDetId::HEX, ieta, iphi);
0465 doTestOneCell(topology, "HEX", cell, all, total, good, bad);
0466 }
0467 }
0468 edm::LogVerbatim("HCalGeom") << "HEX:" << all << ":" << total << " CalibIds with " << good << " good and " << bad
0469 << " bad DetId's\n";
0470 allT += all;
0471 totalT += total;
0472 goodT += good;
0473 badT += bad;
0474
0475 edm::LogVerbatim("HCalGeom") << "\nAnalyzes total of " << allT << ":" << totalT << " CalibIds with " << goodT
0476 << " good and " << badT << " bad DetId's";
0477 }
0478
0479 void HcalDetId2DenseTester::doTestOneCell(const HcalTopology& topology,
0480 const std::string& det,
0481 const HcalCalibDetId& cell,
0482 int& all,
0483 int& total,
0484 int& good,
0485 int& bad) {
0486 std::string error;
0487 ++all;
0488 if (topology.validCalib(cell)) {
0489 ++total;
0490 unsigned int dense = topology.detId2denseIdCALIB(DetId(cell));
0491 DetId id = topology.denseId2detIdCALIB(dense);
0492 if (cell == HcalCalibDetId(id)) {
0493 ++good;
0494 error = "";
0495 } else {
0496 ++bad;
0497 error = "***** ERROR *****";
0498 }
0499 edm::LogVerbatim("HCalGeom") << det << "[" << all << ":" << total << "] " << cell << " Dense " << dense << " o/p "
0500 << HcalCalibDetId(id) << " " << error;
0501 } else {
0502 edm::LogVerbatim("HCalGeom") << det << "[" << all << "] " << cell << "***** INVALID *****";
0503 }
0504 }
0505
0506 std::vector<std::string> HcalDetId2DenseTester::splitString(const std::string& fLine) {
0507 std::vector<std::string> result;
0508 int start = 0;
0509 bool empty = true;
0510 for (unsigned i = 0; i <= fLine.size(); i++) {
0511 if (fLine[i] == ' ' || i == fLine.size()) {
0512 if (!empty) {
0513 std::string item(fLine, start, i - start);
0514 result.push_back(item);
0515 empty = true;
0516 }
0517 start = i + 1;
0518 } else {
0519 if (empty)
0520 empty = false;
0521 }
0522 }
0523 return result;
0524 }
0525
0526
0527 DEFINE_FWK_MODULE(HcalDetId2DenseTester);