File indexing completed on 2023-03-17 13:03:45
0001
0002 #include <array>
0003 #include <fstream>
0004 #include <iostream>
0005 #include <memory>
0006 #include <string>
0007 #include <vector>
0008
0009
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0012
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016
0017 #include "DetectorDescription/Core/interface/DDCompactView.h"
0018 #include "DetectorDescription/Core/interface/DDExpandedView.h"
0019 #include "DetectorDescription/Core/interface/DDSpecifics.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0022 #include "Geometry/HGCalTBCommonData/interface/HGCalTBDDDConstants.h"
0023 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0024
0025 class HGCalTBNumberingTester : public edm::one::EDAnalyzer<> {
0026 public:
0027 explicit HGCalTBNumberingTester(const edm::ParameterSet&);
0028 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029
0030 void beginJob() override {}
0031 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0032 void endJob() override {}
0033
0034 private:
0035 edm::ESGetToken<HGCalTBDDDConstants, IdealGeometryRecord> dddToken_;
0036 std::string nameSense_, nameDetector_;
0037 std::vector<double> positionX_, positionY_;
0038 int increment_;
0039 };
0040
0041 HGCalTBNumberingTester::HGCalTBNumberingTester(const edm::ParameterSet& iC) {
0042 nameSense_ = iC.getParameter<std::string>("nameSense");
0043 nameDetector_ = iC.getParameter<std::string>("nameDevice");
0044 positionX_ = iC.getParameter<std::vector<double> >("localPositionX");
0045 positionY_ = iC.getParameter<std::vector<double> >("localPositionY");
0046 increment_ = iC.getParameter<int>("increment");
0047
0048 dddToken_ = esConsumes<HGCalTBDDDConstants, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_});
0049
0050 std::string unit("mm");
0051 for (unsigned int k = 0; k < positionX_.size(); ++k) {
0052 positionX_[k] /= CLHEP::mm;
0053 positionY_[k] /= CLHEP::mm;
0054 }
0055 edm::LogVerbatim("HGCalGeom") << "Test numbering for " << nameDetector_ << " using constants of " << nameSense_
0056 << " at " << positionX_.size() << " local positions "
0057 << "for every " << increment_ << " layers";
0058 for (unsigned int k = 0; k < positionX_.size(); ++k)
0059 edm::LogVerbatim("HGCalGeom") << "Position[" << k << "] " << positionX_[k] << " " << unit << ", " << positionY_[k]
0060 << " " << unit;
0061 }
0062
0063 void HGCalTBNumberingTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0064 std::vector<double> vecxy;
0065 edm::ParameterSetDescription desc;
0066 desc.add<std::string>("nameSense", "HGCalEESensitive");
0067 desc.add<std::string>("nameDevice", "HGCal EE");
0068 desc.add<std::vector<double> >("localPositionX", vecxy);
0069 desc.add<std::vector<double> >("localPositionY", vecxy);
0070 desc.add<int>("increment", 2);
0071 descriptions.add("hgcalTBNumberingTesterEE", desc);
0072 }
0073
0074
0075 void HGCalTBNumberingTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0076 const HGCalTBDDDConstants& hgdc = iSetup.getData(dddToken_);
0077 edm::LogVerbatim("HGCalGeom") << nameDetector_ << " Layers = " << hgdc.layers(false)
0078 << " Sectors = " << hgdc.sectors() << " Minimum Slope = " << hgdc.minSlope();
0079
0080 edm::LogVerbatim("HGCalGeom") << "Minimum Wafer # " << hgdc.waferMin() << " Mamximum Wafer # " << hgdc.waferMax()
0081 << " Wafer counts " << hgdc.waferCount(0) << ":" << hgdc.waferCount(1);
0082 for (unsigned int i = 0; i < hgdc.layers(true); ++i) {
0083 int lay = i + 1;
0084 double z = hgdc.waferZ(lay, false);
0085 edm::LogVerbatim("HGCalGeom") << "Layer " << lay << " Wafers " << hgdc.wafers(lay, 0) << ":" << hgdc.wafers(lay, 1)
0086 << ":" << hgdc.wafers(lay, 2) << " Z " << z << " R " << hgdc.rangeR(z, false).first
0087 << ":" << hgdc.rangeR(z, false).second;
0088 }
0089
0090 edm::LogVerbatim("HGCalGeom") << std::endl;
0091 std::pair<float, float> xy;
0092 std::string flg;
0093 int subsec(0);
0094 int loff = hgdc.firstLayer();
0095 for (unsigned int k = 0; k < positionX_.size(); ++k) {
0096 float localx(positionX_[k]), localy(positionY_[k]);
0097 for (unsigned int i = 0; i < hgdc.layers(false); ++i) {
0098 std::pair<int, int> kxy, lxy;
0099 kxy = hgdc.assignCell(localx, localy, i + loff, subsec, false);
0100 xy = hgdc.locateCell(kxy.second, i + loff, kxy.first, false);
0101 lxy = hgdc.assignCell(xy.first, xy.second, i + loff, 0, false);
0102 flg = (kxy == lxy) ? " " : " ***** Error *****";
0103 edm::LogVerbatim("HGCalGeom") << "Input: (" << localx << "," << localy << "," << i + loff << ", " << subsec
0104 << "), assignCell o/p (" << kxy.first << ", " << kxy.second << ") locateCell o/p ("
0105 << xy.first << ", " << xy.second << "),"
0106 << " final (" << lxy.first << ", " << lxy.second << ")" << flg;
0107 kxy = hgdc.assignCell(-localx, -localy, i + loff, subsec, false);
0108 xy = hgdc.locateCell(kxy.second, i + loff, kxy.first, false);
0109 lxy = hgdc.assignCell(xy.first, xy.second, i + loff, 0, false);
0110 flg = (kxy == lxy) ? " " : " ***** Error *****";
0111 edm::LogVerbatim("HGCalGeom") << "Input: (" << -localx << "," << -localy << "," << i + loff << ", " << subsec
0112 << "), assignCell o/p (" << kxy.first << ", " << kxy.second << ") locateCell o/p ("
0113 << xy.first << ", " << xy.second << "), final (" << lxy.first << ", " << lxy.second
0114 << ")" << flg;
0115
0116 if (k == 0 && i == 0) {
0117 std::vector<int> ncells = hgdc.numberCells(i + 1, false);
0118 edm::LogVerbatim("HGCalGeom") << "Layer " << i + 1 << " with " << ncells.size() << " rows";
0119 int ntot(0);
0120 for (unsigned int k = 0; k < ncells.size(); ++k) {
0121 ntot += ncells[k];
0122 edm::LogVerbatim("HGCalGeom") << "Row " << k << " with " << ncells[k] << " cells";
0123 }
0124 edm::LogVerbatim("HGCalGeom") << "Total Cells " << ntot << ":" << hgdc.maxCells(i + 1, false);
0125 }
0126 i += increment_;
0127 }
0128 }
0129
0130
0131 unsigned int kk(0);
0132 for (auto const& zz : hgdc.getParameter()->zLayerHex_) {
0133 std::pair<double, double> rr = hgdc.rangeR(zz, true);
0134 edm::LogVerbatim("HGCalGeom") << "[" << kk << "]\t z = " << zz << "\t rMin = " << rr.first
0135 << "\t rMax = " << rr.second;
0136 ++kk;
0137 }
0138 }
0139
0140
0141 DEFINE_FWK_MODULE(HGCalTBNumberingTester);