File indexing completed on 2024-04-06 12:14:39
0001 #include <iostream>
0002 #include <sstream>
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
0015 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0016 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0017 #include "DataFormats/ForwardDetId/interface/ForwardSubdetector.h"
0018 #include "DataFormats/ForwardDetId/interface/HFNoseDetId.h"
0019 #include "CoralBase/Exception.h"
0020
0021 class HFNoseGeometryTester : public edm::one::EDAnalyzer<> {
0022 public:
0023 explicit HFNoseGeometryTester(const edm::ParameterSet&);
0024 ~HFNoseGeometryTester() override {}
0025
0026 void beginJob() override {}
0027 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0028 void endJob() override {}
0029
0030 private:
0031 void doTestWafer(const HGCalGeometry* geom);
0032
0033 const std::string name_;
0034 const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> tokGeom_;
0035 };
0036
0037 HFNoseGeometryTester::HFNoseGeometryTester(const edm::ParameterSet& iC)
0038 : name_(iC.getParameter<std::string>("Detector")),
0039 tokGeom_(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag{"", name_})) {}
0040
0041 void HFNoseGeometryTester::analyze(const edm::Event&, const edm::EventSetup& iSetup) {
0042 const HGCalGeometry* geom = &iSetup.getData(tokGeom_);
0043 if (geom->topology().isHFNose()) {
0044 doTestWafer(geom);
0045 } else {
0046 edm::LogVerbatim("ForwardGeom") << name_ << " is not a valid name for HFNose Detecor";
0047 }
0048 }
0049
0050 void HFNoseGeometryTester::doTestWafer(const HGCalGeometry* geom) {
0051 const std::vector<DetId>& ids = geom->getValidDetIds();
0052 edm::LogVerbatim("ForwardGeom") << "doTestWafer:: " << ids.size() << " valid ids for " << geom->cellElement();
0053 int layers[] = {1, 4, 7};
0054 int zsides[] = {1, -1};
0055 int cells[] = {1, 4, 7};
0056 int wafers[] = {7, 5, 3, -3, -5, -7};
0057 for (int zside : zsides) {
0058 for (int layer : layers) {
0059 for (int waferU : wafers) {
0060 for (int waferV : wafers) {
0061 int type = geom->topology().dddConstants().getTypeHex(layer, waferU, waferV);
0062 edm::LogVerbatim("ForwardGeom")
0063 << "zside " << zside << " layer " << layer << " wafer " << waferU << ":" << waferV << " type " << type;
0064 for (int cellU : cells) {
0065 for (int cellV : cells) {
0066 edm::LogVerbatim("ForwardGeom") << " cell " << cellU << ":" << cellV;
0067 DetId id1 = (DetId)(HFNoseDetId(zside, type, layer, waferU, waferV, cellU, cellV));
0068 edm::LogVerbatim("ForwardGeom") << HFNoseDetId(id1);
0069 if (geom->topology().valid(id1)) {
0070 auto icell1 = geom->getGeometry(id1);
0071 GlobalPoint global1 = geom->getPosition(id1);
0072 DetId idc1 = geom->getClosestCell(global1);
0073 std::ostringstream st1;
0074 st1 << "DetId (" << zside << ":" << type << ":" << layer << ":" << waferU << ":" << waferV << ":"
0075 << cellU << ":" << cellV << ") Geom " << icell1 << " position (" << global1.x() << ", "
0076 << global1.y() << ", " << global1.z() << ") ids " << std::hex << id1.rawId() << ":" << idc1.rawId()
0077 << std::dec << ":" << HFNoseDetId(id1) << ":" << HFNoseDetId(idc1)
0078 << " parameter[3] = " << icell1->param()[2] << ":" << icell1->param()[2];
0079 if (id1.rawId() != idc1.rawId())
0080 st1 << "***** ERROR *****";
0081 edm::LogVerbatim("ForwardGeom") << st1.str();
0082 std::vector<GlobalPoint> corners = geom->getCorners(idc1);
0083 std::ostringstream st2;
0084 st2 << corners.size() << " corners";
0085 for (auto const& cor : corners)
0086 st2 << " [" << cor.x() << "," << cor.y() << "," << cor.z() << "]";
0087 edm::LogVerbatim("ForwardGeom") << st2.str();
0088 }
0089 }
0090 }
0091 }
0092 }
0093 }
0094 }
0095 }
0096
0097
0098 DEFINE_FWK_MODULE(HFNoseGeometryTester);