File indexing completed on 2024-04-06 12:14:49
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include <memory>
0020 #include <iostream>
0021 #include <fstream>
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/EventSetup.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0033 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0034 #include "Geometry/HcalCommonData/interface/HcalDDDRecConstants.h"
0035
0036 class HcalRecNumberingTester : public edm::one::EDAnalyzer<> {
0037 public:
0038 explicit HcalRecNumberingTester(const edm::ParameterSet&);
0039 ~HcalRecNumberingTester() override;
0040
0041 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0042
0043 private:
0044 edm::ESGetToken<HcalDDDRecConstants, HcalRecNumberingRecord> token_;
0045 };
0046
0047 HcalRecNumberingTester::HcalRecNumberingTester(const edm::ParameterSet&)
0048 : token_{esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord>(edm::ESInputTag{})} {}
0049
0050 HcalRecNumberingTester::~HcalRecNumberingTester() {}
0051
0052
0053 void HcalRecNumberingTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0054 const HcalDDDRecConstants& hdc = iSetup.getData(token_);
0055 for (int i = 0; i < 4; ++i)
0056 edm::LogVerbatim("HCalGeom") << "MaxDepth[" << i << "] = " << hdc.getMaxDepth(i);
0057 edm::LogVerbatim("HCalGeom") << "about to getPhiOff and getPhiBin for 0..2";
0058 int neta = hdc.getNEta();
0059 edm::LogVerbatim("HCalGeom") << neta << " eta bins with phi off set for "
0060 << "barrel = " << hdc.getPhiOff(0) << ", endcap = " << hdc.getPhiOff(1);
0061 for (int i = 0; i < neta; ++i) {
0062 std::pair<double, double> etas = hdc.getEtaLimit(i);
0063 double fbin = hdc.getPhiBin(i);
0064 std::vector<int> depths = hdc.getDepth(i, false);
0065 edm::LogVerbatim("HCalGeom") << "EtaBin[" << i << "]: EtaLimit = (" << etas.first << ":" << etas.second
0066 << ") phiBin = " << fbin << " and " << depths.size() << " depths";
0067 for (unsigned int k = 0; k < depths.size(); ++k) {
0068 edm::LogVerbatim("HCalGeom") << "[" << k << "] " << depths[k];
0069 }
0070 }
0071 for (int type = 0; type < 2; ++type) {
0072 std::pair<int, int> etar = hdc.getEtaRange(type);
0073 edm::LogVerbatim("HCalGeom") << "Detector type: " << type << " with eta ranges " << etar.first << ":"
0074 << etar.second;
0075 for (int eta = etar.first; eta <= etar.second; ++eta) {
0076 std::vector<std::pair<int, double>> phis = hdc.getPhis(type + 1, eta);
0077 for (auto& phi : phis) {
0078 edm::LogVerbatim("HCalGeom") << "Type:Eta:phi " << type << ":" << eta << ":" << phi.first
0079 << " Depth range (+z) " << hdc.getMinDepth(type, eta, phi.first, 1) << ":"
0080 << hdc.getMaxDepth(type, eta, phi.first, 1) << " (-z) "
0081 << hdc.getMinDepth(type, eta, phi.first, -1) << ":"
0082 << hdc.getMaxDepth(type, eta, phi.first, -1);
0083 }
0084 }
0085 }
0086 std::vector<HcalDDDRecConstants::HcalEtaBin> hbar = hdc.getEtaBins(0);
0087 std::vector<HcalDDDRecConstants::HcalEtaBin> hcap = hdc.getEtaBins(1);
0088 edm::LogVerbatim("HCalGeom") << "Topology Mode " << hdc.getTopoMode() << " HB with " << hbar.size()
0089 << " eta sectors and HE with " << hcap.size() << " eta sectors";
0090 std::vector<HcalCellType> hbcell = hdc.HcalCellTypes(HcalBarrel);
0091 edm::LogVerbatim("HCalGeom") << "HB with " << hbcell.size() << " cells";
0092 unsigned int i1(0), i2(0), i3(0), i4(0);
0093 for (const auto& cell : hbcell) {
0094 edm::LogVerbatim("HCalGeom") << "HB[" << i1 << "] det " << cell.detType() << " zside " << cell.zside() << ":"
0095 << cell.halfSize() << " RO " << cell.actualReadoutDirection() << " eta "
0096 << cell.etaBin() << ":" << cell.etaMin() << ":" << cell.etaMax() << " phi "
0097 << cell.nPhiBins() << ":" << cell.nPhiModule() << ":" << cell.phiOffset() << ":"
0098 << cell.phiBinWidth() << ":" << cell.unitPhi() << " depth " << cell.depthSegment()
0099 << ":" << cell.depth() << ":" << cell.depthMin() << ":" << cell.depthMax() << ":"
0100 << cell.depthType();
0101 ++i1;
0102 std::vector<std::pair<int, double>> phis = cell.phis();
0103 edm::LogVerbatim("HCalGeom") << "Phis (" << phis.size() << ") :";
0104 for (const auto& phi : phis)
0105 edm::LogVerbatim("HCalGeom") << " [" << phi.first << ", " << phi.second << "]";
0106 }
0107 std::vector<HcalCellType> hecell = hdc.HcalCellTypes(HcalEndcap);
0108 edm::LogVerbatim("HCalGeom") << "HE with " << hecell.size() << " cells";
0109 for (const auto& cell : hecell) {
0110 edm::LogVerbatim("HCalGeom") << "HE[" << i2 << "] det " << cell.detType() << " zside " << cell.zside() << ":"
0111 << cell.halfSize() << " RO " << cell.actualReadoutDirection() << " eta "
0112 << cell.etaBin() << ":" << cell.etaMin() << ":" << cell.etaMax() << " phi "
0113 << cell.nPhiBins() << ":" << cell.nPhiModule() << ":" << cell.phiOffset() << ":"
0114 << cell.phiBinWidth() << ":" << cell.unitPhi() << " depth " << cell.depthSegment()
0115 << ":" << cell.depth() << ":" << cell.depthMin() << ":" << cell.depthMax() << ":"
0116 << cell.depthType();
0117 ++i2;
0118 std::vector<std::pair<int, double>> phis = cell.phis();
0119 edm::LogVerbatim("HCalGeom") << "Phis (" << phis.size() << ") :";
0120 for (const auto& phi : phis)
0121 edm::LogVerbatim("HCalGeom") << " [" << phi.first << ", " << phi.second << "]";
0122 }
0123 std::vector<HcalCellType> hfcell = hdc.HcalCellTypes(HcalForward);
0124 edm::LogVerbatim("HCalGeom") << "HF with " << hfcell.size() << " cells";
0125 for (const auto& cell : hfcell) {
0126 edm::LogVerbatim("HCalGeom") << "HF[" << i3 << "] det " << cell.detType() << " zside " << cell.zside() << ":"
0127 << cell.halfSize() << " RO " << cell.actualReadoutDirection() << " eta "
0128 << cell.etaBin() << ":" << cell.etaMin() << ":" << cell.etaMax() << " phi "
0129 << cell.nPhiBins() << ":" << cell.nPhiModule() << ":" << cell.phiOffset() << ":"
0130 << cell.phiBinWidth() << ":" << cell.unitPhi() << " depth " << cell.depthSegment()
0131 << ":" << cell.depth() << ":" << cell.depthMin() << ":" << cell.depthMax() << ":"
0132 << cell.depthType();
0133 ++i3;
0134 }
0135 std::vector<HcalCellType> hocell = hdc.HcalCellTypes(HcalOuter);
0136 edm::LogVerbatim("HCalGeom") << "HO with " << hocell.size() << " cells";
0137 for (const auto& cell : hocell) {
0138 edm::LogVerbatim("HCalGeom") << "HO[" << i4 << "] det " << cell.detType() << " zside " << cell.zside() << ":"
0139 << cell.halfSize() << " RO " << cell.actualReadoutDirection() << " eta "
0140 << cell.etaBin() << ":" << cell.etaMin() << ":" << cell.etaMax() << " phi "
0141 << cell.nPhiBins() << ":" << cell.nPhiModule() << ":" << cell.phiOffset() << ":"
0142 << cell.phiBinWidth() << ":" << cell.unitPhi() << " depth " << cell.depthSegment()
0143 << ":" << cell.depth() << ":" << cell.depthMin() << ":" << cell.depthMax() << ":"
0144 << cell.depthType();
0145 ++i4;
0146 }
0147 for (int type = 0; type <= 1; ++type) {
0148 std::vector<HcalDDDRecConstants::HcalActiveLength> act = hdc.getThickActive(type);
0149 edm::LogVerbatim("HCalGeom") << "Hcal type " << type << " has " << act.size() << " eta/depth segments";
0150 for (const auto& active : act) {
0151 edm::LogVerbatim("HCalGeom") << "zside " << active.zside << " ieta " << active.ieta << " depth " << active.depth
0152 << " type " << active.stype << " eta " << active.eta << " active thickness "
0153 << active.thick;
0154 }
0155 }
0156
0157
0158 std::vector<int> phiSp;
0159 HcalSubdetector subdet = HcalSubdetector(hdc.dddConstants()->ldMap()->validDet(phiSp));
0160 if (subdet == HcalBarrel || subdet == HcalEndcap) {
0161 int type = (int)(subdet - 1);
0162 std::pair<int, int> etas = hdc.getEtaRange(type);
0163 for (int eta = etas.first; eta <= etas.second; ++eta) {
0164 for (int k : phiSp) {
0165 int zside = (k > 0) ? 1 : -1;
0166 int iphi = (k > 0) ? k : -k;
0167 #ifdef EDM_ML_DEBUG
0168 edm::LogVerbatim("HCalGeom") << "Look for Subdet " << subdet << " Zside " << zside << " Eta " << eta << " Phi "
0169 << iphi << " depths " << hdc.getMinDepth(type, eta, iphi, zside) << ":"
0170 << hdc.getMaxDepth(type, eta, iphi, zside);
0171 #endif
0172 std::vector<HcalDetId> ids;
0173 for (int depth = hdc.getMinDepth(type, eta, iphi, zside); depth <= hdc.getMaxDepth(type, eta, iphi, zside);
0174 ++depth) {
0175 HcalDetId id(subdet, zside * eta, iphi, depth);
0176 HcalDetId hid = hdc.mergedDepthDetId(id);
0177 hdc.unmergeDepthDetId(hid, ids);
0178 edm::LogVerbatim("HCalGeom") << "Input ID " << id << " Merged ID " << hid << " containing " << ids.size()
0179 << " IDS:";
0180 for (auto id : ids)
0181 edm::LogVerbatim("HCalGeom") << " " << id;
0182 }
0183 }
0184 }
0185 }
0186
0187 for (const auto& cell : hbcell) {
0188 int ieta = cell.etaBin() * cell.zside();
0189 double rz = hdc.getRZ(HcalBarrel, ieta, cell.phis()[0].first, cell.depthSegment());
0190 edm::LogVerbatim("HCalGeom") << "HB (eta=" << ieta << ", phi=" << cell.phis()[0].first
0191 << ", depth=" << cell.depthSegment() << ") r/z = " << rz;
0192 }
0193 for (const auto& cell : hecell) {
0194 int ieta = cell.etaBin() * cell.zside();
0195 double rz = hdc.getRZ(HcalEndcap, ieta, cell.phis()[0].first, cell.depthSegment());
0196 edm::LogVerbatim("HCalGeom") << "HE (eta=" << ieta << ", phi=" << cell.phis()[0].first
0197 << ", depth=" << cell.depthSegment() << ") r/z = " << rz;
0198 }
0199 }
0200
0201
0202 DEFINE_FWK_MODULE(HcalRecNumberingTester);