File indexing completed on 2023-10-25 10:03:57
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 "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0025 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0028
0029 #include "FWCore/Framework/interface/Event.h"
0030 #include "FWCore/Framework/interface/EventSetup.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034
0035 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0036 #include "Geometry/Records/interface/HcalSimNumberingRecord.h"
0037 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0038 #include "Geometry/HcalCommonData/interface/HcalDDDSimConstants.h"
0039 #include "Geometry/HcalCommonData/interface/HcalDDDRecConstants.h"
0040 #include "Geometry/HcalCommonData/interface/HcalNumberingFromDDD.h"
0041 #include "Geometry/HcalCommonData/interface/HcalHitRelabeller.h"
0042 #include "SimG4CMS/Calo/interface/HcalNumberingScheme.h"
0043 #include "SimG4CMS/Calo/interface/HcalTestNumberingScheme.h"
0044
0045 class HcalTestNumberingTester : public edm::one::EDAnalyzer<> {
0046 public:
0047 explicit HcalTestNumberingTester(const edm::ParameterSet&);
0048 ~HcalTestNumberingTester() override = default;
0049 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0050
0051 void beginJob() override {}
0052 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0053 void endJob() override {}
0054
0055 private:
0056 const edm::ESGetToken<HcalDDDSimConstants, HcalSimNumberingRecord> tokSim_;
0057 const edm::ESGetToken<HcalDDDRecConstants, HcalRecNumberingRecord> tokReco_;
0058 };
0059
0060 HcalTestNumberingTester::HcalTestNumberingTester(const edm::ParameterSet&)
0061 : tokSim_(esConsumes<HcalDDDSimConstants, HcalSimNumberingRecord>()),
0062 tokReco_(esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord>()) {}
0063
0064 void HcalTestNumberingTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0065 edm::ParameterSetDescription desc;
0066 descriptions.add("hcalTestNumberingTest", desc);
0067 }
0068
0069
0070 void HcalTestNumberingTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0071 const HcalDDDSimConstants* hcs = &iSetup.getData(tokSim_);
0072 const HcalDDDRecConstants* hcr = &iSetup.getData(tokReco_);
0073 HcalNumberingScheme* schme1 = new HcalNumberingScheme();
0074 HcalNumberingScheme* schme2 = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(false));
0075
0076 for (int type = 0; type < 2; ++type) {
0077 HcalSubdetector sub = (type == 0) ? HcalBarrel : HcalEndcap;
0078 for (int zs = 0; zs < 2; ++zs) {
0079 int zside = 2 * zs - 1;
0080 std::pair<int, int> etas = hcr->getEtaRange(type);
0081 for (int eta = etas.first; eta <= etas.second; ++eta) {
0082 std::vector<std::pair<int, double> > phis = hcr->getPhis(sub, eta);
0083 for (unsigned int k = 0; k < phis.size(); ++k) {
0084 int phi = phis[k].first;
0085 int lmin = (type == 1 && eta == 16) ? 8 : 1;
0086 int lmax = (type == 1) ? 19 : ((eta == 16) ? 7 : 17);
0087 for (int lay = lmin; lay <= lmax; ++lay) {
0088 std::pair<int, int> etd = hcs->getEtaDepth(sub, eta, phi, zside, 0, lay);
0089 HcalNumberingFromDDD::HcalID tmp(sub, zs, etd.second, etd.first, phi, phi, lay);
0090 uint32_t id1 = schme1->getUnitID(tmp);
0091 uint32_t id2 = schme2->getUnitID(tmp);
0092 DetId id0 = HcalHitRelabeller::relabel(id2, hcr);
0093 std::string ok = (id1 != id0.rawId()) ? " *** ERROR ***" : " ";
0094 edm::LogVerbatim("HcalSim") << "I/P " << sub << ":" << zside * eta << ":" << phi << ":" << lay << " Normal "
0095 << std::hex << id1 << std::dec << " " << HcalDetId(id1) << " Test " << std::hex
0096 << id2 << std::dec << " " << HcalDetId(id0) << ok;
0097 }
0098 }
0099 }
0100 }
0101 }
0102 }
0103
0104
0105 DEFINE_FWK_MODULE(HcalTestNumberingTester);