File indexing completed on 2024-04-06 12:14:26
0001
0002
0003
0004
0005
0006 #include <memory>
0007
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "DataFormats/Math/interface/Rounding.h"
0014
0015 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0016
0017 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0018 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0019
0020 #include <iostream>
0021 #include <string>
0022 #include <cmath>
0023 #include <iomanip> // for setw() etc.
0024 #include <vector>
0025
0026 using namespace std;
0027 using namespace cms_rounding;
0028
0029 class DTGeometryAnalyzer : public edm::one::EDAnalyzer<> {
0030 public:
0031 DTGeometryAnalyzer(const edm::ParameterSet& pset);
0032 ~DTGeometryAnalyzer() override;
0033
0034 void beginJob() override {}
0035 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0036 void endJob() override {}
0037
0038 private:
0039 const string& myName() { return myName_; }
0040 const edm::ESGetToken<DTGeometry, MuonGeometryRecord> tokDT_;
0041 const int dashedLineWidth_;
0042 const string dashedLine_;
0043 const string myName_;
0044 double tolerance_;
0045 };
0046
0047 DTGeometryAnalyzer::DTGeometryAnalyzer(const edm::ParameterSet& iConfig)
0048 : tokDT_{esConsumes<DTGeometry, MuonGeometryRecord>(edm::ESInputTag{})},
0049 dashedLineWidth_(104),
0050 dashedLine_(string(dashedLineWidth_, '-')),
0051 myName_("DTGeometryAnalyzer"),
0052 tolerance_(iConfig.getUntrackedParameter<double>("tolerance", 1.e-23)) {}
0053
0054 DTGeometryAnalyzer::~DTGeometryAnalyzer() {}
0055
0056 void DTGeometryAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0057 const auto& pDD = iSetup.getData(tokDT_);
0058
0059 edm::LogVerbatim("DTGeometry") << myName() << ": Analyzer...";
0060 edm::LogVerbatim("DTGeometry") << "start " << dashedLine_;
0061
0062 edm::LogVerbatim("DTGeometry") << " Geometry node for DTGeom is " << &pDD;
0063 edm::LogVerbatim("DTGeometry") << " I have " << pDD.detTypes().size() << " detTypes";
0064 edm::LogVerbatim("DTGeometry") << " I have " << pDD.detUnits().size() << " detUnits";
0065 edm::LogVerbatim("DTGeometry") << " I have " << pDD.dets().size() << " dets";
0066 edm::LogVerbatim("DTGeometry") << " I have " << pDD.layers().size() << " layers";
0067 edm::LogVerbatim("DTGeometry") << " I have " << pDD.superLayers().size() << " superlayers";
0068 edm::LogVerbatim("DTGeometry") << " I have " << pDD.chambers().size() << " chambers";
0069
0070 edm::LogVerbatim("DTGeometry") << myName() << ": Begin iteration over geometry...";
0071 edm::LogVerbatim("DTGeometry") << "iter " << dashedLine_;
0072
0073
0074 for (const auto& det : pDD.detUnits()) {
0075 DetId detId = det->geographicalId();
0076 int id = detId();
0077 const GeomDet* gdet_ = pDD.idToDet(detId);
0078 const GeomDetUnit* gdet = pDD.idToDetUnit(detId);
0079 const DTLayer* lay = dynamic_cast<const DTLayer*>(gdet);
0080 edm::LogVerbatim("DTGeometry") << "GeomDetUnit is of type " << detId.det() << " and raw id = " << id;
0081 assert(det == gdet);
0082 assert(gdet_ == gdet);
0083 assert(gdet_ == lay);
0084 }
0085
0086
0087 edm::LogVerbatim("DTGeometry") << "LAYERS " << dashedLine_;
0088 for (auto det : pDD.layers()) {
0089 const DTTopology& topo = det->specificTopology();
0090 const BoundPlane& surf = det->surface();
0091 edm::LogVerbatim("DTGeometry") << "Layer " << det->id() << " SL " << det->superLayer()->id() << " chamber "
0092 << det->chamber()->id() << " Topology W/H/L: " << topo.cellWidth() << "/"
0093 << topo.cellHeight() << "/" << topo.cellLenght() << " first/last/# wire "
0094 << topo.firstChannel() << "/" << topo.lastChannel() << "/" << topo.channels()
0095 << " Position " << surf.position() << " normVect "
0096 << roundVecIfNear0(surf.normalVector(), tolerance_)
0097 << " bounds W/H/L: " << surf.bounds().width() << "/" << surf.bounds().thickness()
0098 << "/" << surf.bounds().length();
0099 }
0100
0101
0102 edm::LogVerbatim("DTGeometry") << "SUPERLAYERS " << dashedLine_;
0103 for (auto det : pDD.superLayers()) {
0104 const BoundPlane& surf = det->surface();
0105 edm::LogVerbatim("DTGeometry") << "SuperLayer " << det->id() << " chamber " << det->chamber()->id() << " Position "
0106 << surf.position() << " normVect "
0107 << roundVecIfNear0(surf.normalVector(), tolerance_)
0108 << " bounds W/H/L: " << surf.bounds().width() << "/" << surf.bounds().thickness()
0109 << "/" << surf.bounds().length();
0110 }
0111
0112
0113 edm::LogVerbatim("DTGeometry") << "CHAMBERS " << dashedLine_;
0114 for (auto det : pDD.chambers()) {
0115
0116 const BoundPlane& surf = det->surface();
0117
0118 edm::LogVerbatim("DTGeometry") << "Chamber " << det->id() << " Position " << surf.position() << " normVect "
0119 << roundVecIfNear0(surf.normalVector(), tolerance_)
0120 << " bounds W/H/L: " << surf.bounds().width() << "/" << surf.bounds().thickness()
0121 << "/" << surf.bounds().length();
0122 }
0123 edm::LogVerbatim("DTGeometry") << "END " << dashedLine_;
0124
0125
0126 for (int w = -2; w <= 2; ++w) {
0127 for (int st = 1; st <= 4; ++st) {
0128 for (int se = 1; se <= ((st == 4) ? 14 : 12); ++se) {
0129 DTChamberId id(w, st, se);
0130 const DTChamber* ch = pDD.chamber(id);
0131 if (!ch)
0132 edm::LogVerbatim("DTGeometry") << "ERROR ch not found " << id;
0133 else {
0134 if (id != ch->id())
0135 edm::LogVerbatim("DTGeometry")
0136 << "ERROR: got wrong chamber: Cerco camera " << id << " e trovo " << ch->id();
0137
0138 const GeomDet* gdetc = pDD.idToDet(id);
0139 assert(gdetc == ch);
0140
0141 for (int sl = 1; sl <= 3; ++sl) {
0142 if (sl == 2 && st == 4)
0143 continue;
0144 DTSuperLayerId slid(id, sl);
0145 const DTSuperLayer* dtsl = pDD.superLayer(slid);
0146 if (!dtsl)
0147 edm::LogVerbatim("DTGeometry") << "ERROR sl not found " << slid;
0148 else {
0149 if (slid != dtsl->id())
0150 edm::LogVerbatim("DTGeometry") << "ERROR: got wrong sl! Cerco sl " << slid << " e trovo " << dtsl->id();
0151
0152 const GeomDet* gdets = pDD.idToDet(slid);
0153 assert(gdets == dtsl);
0154
0155 for (int l = 1; l <= 4; ++l) {
0156 DTLayerId lid(slid, l);
0157 const DTLayer* lay = pDD.layer(lid);
0158 if (!lay)
0159 edm::LogVerbatim("DTGeometry") << "ERROR lay not found " << lid;
0160 if (lid != lay->id())
0161 edm::LogVerbatim("DTGeometry")
0162 << "ERROR: got wrong layer Cerco lay " << lid << " e trovo " << lay->id();
0163
0164 const GeomDet* gdetl = pDD.idToDet(lid);
0165 assert(gdetl == lay);
0166 }
0167 }
0168 }
0169 }
0170 }
0171 }
0172 }
0173 edm::LogVerbatim("DTGeometry") << "END " << dashedLine_;
0174 }
0175
0176
0177 #include <FWCore/Framework/interface/MakerMacros.h>
0178 DEFINE_FWK_MODULE(DTGeometryAnalyzer);