Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-06-27 22:37:15

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