Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:14:25

0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 #include "FWCore/Framework/interface/ESHandle.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0006 #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
0007 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
0008 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0009 
0010 #include <string>
0011 #include <iomanip>  // for setw() etc.
0012 #include <vector>
0013 
0014 class CSCGeometryAsChambers : public edm::one::EDAnalyzer<> {
0015 public:
0016   explicit CSCGeometryAsChambers(const edm::ParameterSet&);
0017   ~CSCGeometryAsChambers() override = default;
0018 
0019   void beginJob() override {}
0020   void analyze(edm::Event const&, edm::EventSetup const&) override;
0021   void endJob() override {}
0022 
0023   const std::string& myName() { return myName_; }
0024 
0025 private:
0026   const int dashedLineWidth_;
0027   const std::string dashedLine_;
0028   const std::string myName_;
0029   const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> tokGeom_;
0030 };
0031 
0032 CSCGeometryAsChambers::CSCGeometryAsChambers(const edm::ParameterSet& iConfig)
0033     : dashedLineWidth_(132),
0034       dashedLine_(std::string(dashedLineWidth_, '-')),
0035       myName_("CSCGeometryAsChambers"),
0036       tokGeom_(esConsumes()) {}
0037 
0038 void CSCGeometryAsChambers::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0039   std::cout << myName() << ": Analyzer..." << std::endl;
0040   std::cout << "start " << dashedLine_ << std::endl;
0041 
0042   const edm::ESHandle<CSCGeometry>& pDD = iSetup.getHandle(tokGeom_);
0043   std::cout << " Geometry node for CSCGeom is  " << &(*pDD) << std::endl;
0044   std::cout << " I have " << pDD->dets().size() << " detectors" << std::endl;
0045   std::cout << " I have " << pDD->detTypes().size() << " types"
0046             << "\n"
0047             << std::endl;
0048   std::cout << " I have " << pDD->detUnits().size() << " detUnits" << std::endl;
0049   std::cout << " I have " << pDD->dets().size() << " dets" << std::endl;
0050   std::cout << " I have " << pDD->layers().size() << " layers" << std::endl;
0051   std::cout << " I have " << pDD->chambers().size() << " chambers" << std::endl;
0052 
0053   std::cout << myName() << ": Begin iteration over geometry..." << std::endl;
0054 
0055   int icount = 0;
0056 
0057   const CSCGeometry::ChamberContainer& vc = pDD->chambers();
0058   std::cout << "No. of chambers stored = " << vc.size() << std::endl;
0059 
0060   std::cout << "\n  #     id(dec)      id(oct)    labels      length       width      thickness   "
0061                "  g(x=0)   g(y=0)   g(z=0)  g(z=-1)  g(z=+1) "
0062                "  phi(0)"
0063             << std::endl;
0064   std::cout << dashedLine_ << std::endl;
0065 
0066   for (auto chamber : vc) {
0067     if (chamber) {
0068       ++icount;
0069 
0070       // What's its DetId?
0071 
0072       CSCDetId detId = chamber->id();
0073       int id = detId();  // or detId.rawId()
0074 
0075       //    std::cout << "GeomDetUnit is of type " << detId.det() << " and raw id = " << id << std::endl;
0076 
0077       // There's going to be a lot of messing with field width (and precision) so
0078       // save input values...
0079       int iw = std::cout.width();      // save current width
0080       int ip = std::cout.precision();  // save current precision
0081 
0082       //    std::cout << "\n" << "Parameters of chamber# " <<
0083       std::cout << std::setw(4) << icount << std::setw(12) << id << std::oct << std::setw(12) << id << std::dec
0084                 << std::setw(iw) << "   E" << detId.endcap() << " S" << detId.station() << " R" << detId.ring() << " C"
0085                 << std::setw(2) << detId.chamber() << std::setw(iw);
0086 
0087       // What's its surface?
0088       // The surface knows how to transform local <-> global
0089 
0090       const Surface& bSurface = chamber->surface();
0091 
0092       //    std::cout << " length=" << bSurface.bounds().length() <<
0093       //                     ", width=" << bSurface.bounds().width() <<
0094       //                     ", thickness=" << bSurface.bounds().thickness() << std::endl;
0095 
0096       std::cout << std::setw(12) << std::setprecision(8) << bSurface.bounds().length() << std::setw(12)
0097                 << std::setprecision(8) << bSurface.bounds().width() << std::setw(12) << std::setprecision(6)
0098                 << bSurface.bounds().thickness();
0099 
0100       // Check global coordinates of centre of CSCChamber, and how
0101       // local z direction relates to global z direction
0102 
0103       LocalPoint lCentre(0., 0., 0.);
0104       GlobalPoint gCentre = bSurface.toGlobal(lCentre);
0105 
0106       LocalPoint lCentre1(0., 0., -1.);
0107       GlobalPoint gCentre1 = bSurface.toGlobal(lCentre1);
0108 
0109       LocalPoint lCentre2(0., 0., 1.);
0110       GlobalPoint gCentre2 = bSurface.toGlobal(lCentre2);
0111 
0112       //    std::cout << "local(0,0,-1) = global " << gCentre1 << std::endl;
0113       //    std::cout << "local(0,0)    = global " << gCentre  << std::endl;
0114       //    std::cout << "local(0,0,+1) = global " << gCentre2 << std::endl;
0115 
0116       double gx = gCentre.x();
0117       double gy = gCentre.y();
0118       double gz = gCentre.z();
0119       double gz1 = gCentre1.z();
0120       double gz2 = gCentre2.z();
0121       if (fabs(gx) < 1.e-06)
0122         gx = 0.;
0123       if (fabs(gy) < 1.e-06)
0124         gy = 0.;
0125       if (fabs(gz) < 1.e-06)
0126         gz = 0.;
0127       if (fabs(gz1) < 1.e-06)
0128         gz1 = 0.;
0129       if (fabs(gz2) < 1.e-06)
0130         gz2 = 0.;
0131 
0132       int now = 9;
0133       int nop = 5;
0134       std::cout << std::setw(now) << std::setprecision(nop) << gx << std::setw(now) << std::setprecision(nop) << gy
0135                 << std::setw(now) << std::setprecision(nop) << gz << std::setw(now) << std::setprecision(nop) << gz1
0136                 << std::setw(now) << std::setprecision(nop) << gz2;
0137 
0138       // Global Phi of centre of CSCChamber
0139 
0140       //@@ CARE The following attempted conversion to degrees can be easily
0141       // subverted by GeometryVector/Phi.h enforcing its range convention!
0142       // Either a) use a separate local double before scaling...
0143       //        double cphi = gCentre.phi();
0144       //        double cphiDeg = cphi * radToDeg;
0145       // or b) use Phi's degree conversion...
0146       double cphiDeg = gCentre.phi().degrees();
0147 
0148       // I want to display in range 0 to 360
0149 
0150       // Handle some occasional ugly precision problems around zero
0151       if (fabs(cphiDeg) < 1.e-06) {
0152         cphiDeg = 0.;
0153       } else if (cphiDeg < 0.) {
0154         cphiDeg += 360.;
0155       } else if (cphiDeg >= 360.) {
0156         std::cout << "WARNING: resetting phi= " << cphiDeg << " to zero." << std::endl;
0157         cphiDeg = 0.;
0158       }
0159 
0160       //        int iphiDeg = static_cast<int>( cphiDeg );
0161       //    std::cout << "phi(0,0,0) = " << iphiDeg << " degrees" << std::endl;
0162 
0163       now = 9;
0164       nop = 4;
0165       std::cout << std::setw(now) << std::setprecision(nop) << cphiDeg << std::endl;
0166 
0167       // Reset the values we changed
0168       std::cout << std::setprecision(ip) << std::setw(iw);
0169 
0170     } else {
0171       std::cout << "*** DANGER WILL ROBINSON! *** Stored a null CSCChamber*?  " << std::endl;
0172     }
0173   }
0174 
0175   std::cout << dashedLine_ << " end" << std::endl;
0176 }
0177 
0178 //define this as a plug-in
0179 DEFINE_FWK_MODULE(CSCGeometryAsChambers);