File indexing completed on 2024-07-16 22:52:23
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <array>
0022 #include <fstream>
0023 #include <iostream>
0024 #include <memory>
0025 #include <string>
0026 #include <vector>
0027
0028
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0031
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/EventSetup.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035
0036 #include "DetectorDescription/Core/interface/DDCompactView.h"
0037 #include "DetectorDescription/Core/interface/DDExpandedView.h"
0038 #include "DetectorDescription/Core/interface/DDSpecifics.h"
0039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0040 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0041 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0042 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0043
0044 class HGCalNumberingTester : public edm::one::EDAnalyzer<> {
0045 public:
0046 explicit HGCalNumberingTester(const edm::ParameterSet&);
0047 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0048
0049 void beginJob() override {}
0050 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0051 void endJob() override {}
0052
0053 private:
0054 edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord> dddToken_;
0055 std::string nameSense_, nameDetector_;
0056 std::vector<double> positionX_, positionY_;
0057 int increment_, detType_;
0058 bool reco_;
0059 };
0060
0061 HGCalNumberingTester::HGCalNumberingTester(const edm::ParameterSet& iC) {
0062 nameSense_ = iC.getParameter<std::string>("NameSense");
0063 nameDetector_ = iC.getParameter<std::string>("NameDevice");
0064 positionX_ = iC.getParameter<std::vector<double> >("LocalPositionX");
0065 positionY_ = iC.getParameter<std::vector<double> >("LocalPositionY");
0066 increment_ = iC.getParameter<int>("Increment");
0067 detType_ = iC.getParameter<int>("DetType");
0068 reco_ = iC.getParameter<bool>("Reco");
0069
0070 dddToken_ = esConsumes<HGCalDDDConstants, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_});
0071
0072 std::string unit("mm");
0073 if (reco_) {
0074 for (unsigned int k = 0; k < positionX_.size(); ++k) {
0075 positionX_[k] /= CLHEP::cm;
0076 positionY_[k] /= CLHEP::cm;
0077 }
0078 unit = "cm";
0079 } else {
0080 for (unsigned int k = 0; k < positionX_.size(); ++k) {
0081 positionX_[k] /= CLHEP::mm;
0082 positionY_[k] /= CLHEP::mm;
0083 }
0084 }
0085 edm::LogVerbatim("HGCalGeom") << "Test numbering for " << nameDetector_ << " using constants of " << nameSense_
0086 << " at " << positionX_.size() << " local positions "
0087 << "for every " << increment_ << " layers for DetType " << detType_ << " and RecoFlag "
0088 << reco_;
0089 for (unsigned int k = 0; k < positionX_.size(); ++k)
0090 edm::LogVerbatim("HGCalGeom") << "Position[" << k << "] " << positionX_[k] << " " << unit << ", " << positionY_[k]
0091 << " " << unit;
0092 }
0093
0094 void HGCalNumberingTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0095 std::vector<double> vecxy;
0096 edm::ParameterSetDescription desc;
0097 desc.add<std::string>("NameSense", "HGCalEESensitive");
0098 desc.add<std::string>("NameDevice", "HGCal EE");
0099 desc.add<std::vector<double> >("LocalPositionX", vecxy);
0100 desc.add<std::vector<double> >("LocalPositionY", vecxy);
0101 desc.add<int>("Increment", 19);
0102 desc.add<int>("DetType", 2);
0103 desc.add<bool>("Reco", false);
0104 descriptions.add("hgcalNumberingTesterEE", desc);
0105 }
0106
0107
0108 void HGCalNumberingTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0109 const HGCalDDDConstants& hgdc = iSetup.getData(dddToken_);
0110 edm::LogVerbatim("HGCalGeom") << nameDetector_ << " Layers = " << hgdc.layers(reco_)
0111 << " Sectors = " << hgdc.sectors() << " Minimum Slope = " << hgdc.minSlope();
0112 if (detType_ != 0) {
0113 edm::LogVerbatim("HGCalGeom") << "Minimum Wafer # " << hgdc.waferMin() << " Mamximum Wafer # " << hgdc.waferMax()
0114 << " Wafer counts " << hgdc.waferCount(0) << ":" << hgdc.waferCount(1);
0115 for (unsigned int i = 0; i < hgdc.layers(true); ++i) {
0116 int lay = i + 1;
0117 double z = hgdc.waferZ(lay, reco_);
0118 edm::LogVerbatim("HGCalGeom") << "Layer " << lay << " Wafers " << hgdc.wafers(lay, 0) << ":"
0119 << hgdc.wafers(lay, 1) << ":" << hgdc.wafers(lay, 2) << " Z " << z << " R "
0120 << hgdc.rangeR(z, reco_).first << ":" << hgdc.rangeR(z, reco_).second;
0121 }
0122 }
0123 edm::LogVerbatim("HGCalGeom") << std::endl;
0124 std::pair<float, float> xy;
0125 std::string flg;
0126 int subsec(0);
0127 int loff = hgdc.firstLayer();
0128 double scl = (reco_ ? 10.0 : 1.0);
0129 for (unsigned int k = 0; k < positionX_.size(); ++k) {
0130 float localx(positionX_[k]), localy(positionY_[k]);
0131 for (unsigned int i = 0; i < hgdc.layers(reco_); ++i) {
0132 double zpos = hgdc.waferZ(i + loff, reco_);
0133 int zside = (zpos > 0) ? 1 : -1;
0134 if (detType_ == 1) {
0135 std::pair<int, int> kxy, lxy;
0136 kxy = hgdc.assignCell(localx, localy, i + loff, subsec, reco_);
0137 xy = hgdc.locateCell(kxy.second, i + loff, kxy.first, reco_);
0138 lxy = hgdc.assignCell(xy.first, xy.second, i + loff, 0, reco_);
0139 flg = (kxy == lxy) ? " " : " ***** Error *****";
0140 edm::LogVerbatim("HGCalGeom") << "Input: (" << localx << "," << localy << "," << i + loff << ", " << subsec
0141 << "), assignCell o/p (" << kxy.first << ", " << kxy.second
0142 << ") locateCell o/p (" << xy.first << ", " << xy.second << "),"
0143 << " final (" << lxy.first << ", " << lxy.second << ")" << flg;
0144 kxy = hgdc.assignCell(-localx, -localy, i + loff, subsec, reco_);
0145 xy = hgdc.locateCell(kxy.second, i + loff, kxy.first, reco_);
0146 lxy = hgdc.assignCell(xy.first, xy.second, i + loff, 0, reco_);
0147 flg = (kxy == lxy) ? " " : " ***** Error *****";
0148 edm::LogVerbatim("HGCalGeom") << "Input: (" << -localx << "," << -localy << "," << i + loff << ", " << subsec
0149 << "), assignCell o/p (" << kxy.first << ", " << kxy.second
0150 << ") locateCell o/p (" << xy.first << ", " << xy.second << "), final ("
0151 << lxy.first << ", " << lxy.second << ")" << flg;
0152 } else if (detType_ == 0) {
0153 std::array<int, 3> kxy, lxy;
0154 kxy = hgdc.assignCellTrap(localx, localy, zpos, i + loff, reco_);
0155 xy = hgdc.locateCellTrap(zside, i + loff, kxy[0], kxy[1], reco_, false);
0156 lxy = hgdc.assignCellTrap(xy.first, xy.second, zpos, i + loff, reco_);
0157 flg = (kxy == lxy) ? " " : " ***** Error *****";
0158 edm::LogVerbatim("HGCalGeom") << "Input: (" << localx << "," << localy << "," << zpos << ", " << i + loff
0159 << "), assignCell o/p (" << kxy[0] << ":" << kxy[1] << ":" << kxy[2]
0160 << ") locateCell o/p (" << xy.first << ", " << xy.second << "), final (" << lxy[0]
0161 << ":" << lxy[1] << ":" << lxy[2] << ") Dist "
0162 << hgdc.distFromEdgeTrap(scl * localx, scl * localy, scl * zpos) << " " << flg;
0163 kxy = hgdc.assignCellTrap(-localx, -localy, zpos, i + loff, reco_);
0164 xy = hgdc.locateCellTrap(zside, i + loff, kxy[0], kxy[1], reco_, false);
0165 lxy = hgdc.assignCellTrap(xy.first, xy.second, zpos, i + loff, reco_);
0166 flg = (kxy == lxy) ? " " : " ***** Error *****";
0167 edm::LogVerbatim("HGCalGeom") << "Input: (" << -localx << "," << -localy << "," << zpos << ", " << i + loff
0168 << "), assignCell o/p (" << kxy[0] << ":" << kxy[1] << ":" << kxy[2]
0169 << ") locateCell o/p (" << xy.first << ", " << xy.second << "), final (" << lxy[0]
0170 << ":" << lxy[1] << ":" << lxy[2] << ") Dist "
0171 << hgdc.distFromEdgeTrap(scl * localx, scl * localy, scl * zpos) << " " << flg;
0172 } else {
0173 std::array<int, 5> kxy, lxy;
0174 kxy = hgdc.assignCellHex(localx, localy, zside, i + loff, reco_, false, false);
0175 xy = hgdc.locateCell(zside, i + loff, kxy[0], kxy[1], kxy[3], kxy[4], reco_, true, false, false, false);
0176 lxy = hgdc.assignCellHex(xy.first, xy.second, zside, i + loff, reco_, false, false);
0177 flg = (kxy == lxy) ? " " : " ***** Error *****";
0178 edm::LogVerbatim("HGCalGeom") << "Input: (" << localx << "," << localy << ", " << i + loff
0179 << "), assignCell o/p (" << kxy[0] << ":" << kxy[1] << ":" << kxy[2] << ":"
0180 << kxy[3] << ":" << kxy[4] << ") locateCell o/p (" << xy.first << ", "
0181 << xy.second << "), final (" << lxy[0] << ":" << lxy[1] << ":" << lxy[2] << ":"
0182 << lxy[3] << ":" << lxy[4] << ") Dist "
0183 << hgdc.distFromEdgeHex(scl * localx, scl * localy, scl * zpos) << " " << flg;
0184 kxy = hgdc.assignCellHex(-localx, -localy, zside, i + loff, reco_, false, false);
0185 xy = hgdc.locateCell(zside, i + loff, kxy[0], kxy[1], kxy[3], kxy[4], reco_, true, false, false, false);
0186 lxy = hgdc.assignCellHex(xy.first, xy.second, zside, i + loff, reco_, false, false);
0187 flg = (kxy == lxy) ? " " : " ***** Error *****";
0188 edm::LogVerbatim("HGCalGeom") << "Input: (" << -localx << "," << -localy << ", " << i + loff
0189 << "), assignCell o/p (" << kxy[0] << ":" << kxy[1] << ":" << kxy[2] << ":"
0190 << kxy[3] << ":" << kxy[4] << ") locateCell o/p (" << xy.first << ", "
0191 << xy.second << "), final (" << lxy[0] << ":" << lxy[1] << ":" << lxy[2] << ":"
0192 << lxy[3] << ":" << lxy[4] << ") Dist "
0193 << hgdc.distFromEdgeHex(scl * localx, scl * localy, scl * zpos) << " " << flg;
0194 }
0195 if (k == 0 && i == 0 && detType_ == 1) {
0196 std::vector<int> ncells = hgdc.numberCells(i + 1, reco_);
0197 edm::LogVerbatim("HGCalGeom") << "Layer " << i + 1 << " with " << ncells.size() << " rows";
0198 int ntot(0);
0199 for (unsigned int k = 0; k < ncells.size(); ++k) {
0200 ntot += ncells[k];
0201 edm::LogVerbatim("HGCalGeom") << "Row " << k << " with " << ncells[k] << " cells";
0202 }
0203 edm::LogVerbatim("HGCalGeom") << "Total Cells " << ntot << ":" << hgdc.maxCells(i + 1, reco_);
0204 }
0205 i += increment_;
0206 }
0207 }
0208
0209
0210 if (hgdc.getParameter()->detectorType_ > 0) {
0211 unsigned int kk(0);
0212 for (auto const& zz : hgdc.getParameter()->zLayerHex_) {
0213 std::pair<double, double> rr = hgdc.rangeR(zz, true);
0214 edm::LogVerbatim("HGCalGeom") << "[" << kk << "]\t z = " << zz << "\t rMin = " << rr.first
0215 << "\t rMax = " << rr.second;
0216 ++kk;
0217 }
0218 }
0219
0220
0221 if (hgdc.tileTrapezoid()) {
0222 unsigned int kk(0);
0223 for (auto const& lay : hgdc.getParameter()->layer_) {
0224 auto rRange = hgdc.getREtaRange(lay);
0225 edm::LogVerbatim("HGCalGeom") << "[" << kk << "] Layer " << lay << " R/Eta " << rRange.first << ":"
0226 << rRange.second << " nPhi " << hgdc.getPhiBins(lay);
0227 ++kk;
0228 }
0229 }
0230 }
0231
0232
0233 DEFINE_FWK_MODULE(HGCalNumberingTester);