CSCGeometryAnalyzer

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "Geometry/Records/interface/MuonGeometryRecord.h"
#include "Geometry/CSCGeometry/interface/CSCGeometry.h"
#include "Geometry/CSCGeometry/interface/CSCLayer.h"
#include "DataFormats/GeometryVector/interface/Pi.h"
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"

#include <string>
#include <iomanip>  // for setw() etc.
#include <vector>

class CSCGeometryAnalyzer : public edm::one::EDAnalyzer<> {
public:
  explicit CSCGeometryAnalyzer(const edm::ParameterSet&);
  ~CSCGeometryAnalyzer() override;

  void beginJob() override {}
  void analyze(edm::Event const&, edm::EventSetup const&) override;
  void endJob() override {}

  const std::string& myName() { return myName_; }

private:
  const int dashedLineWidth_;
  const std::string dashedLine_;
  const std::string myName_;
  const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> ddToken_;
};

CSCGeometryAnalyzer::CSCGeometryAnalyzer(const edm::ParameterSet& iConfig)
    : dashedLineWidth_(194),
      dashedLine_(std::string(dashedLineWidth_, '-')),
      myName_("CSCGeometryAnalyzer"),
      ddToken_(esConsumes<CSCGeometry, MuonGeometryRecord>()) {}

CSCGeometryAnalyzer::~CSCGeometryAnalyzer() {}

void CSCGeometryAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
  const double dPi = Geom::pi();
  const double radToDeg = 180. / dPi;

  std::cout << myName() << ": Analyzer..." << std::endl;
  std::cout << "start " << dashedLine_ << std::endl;
  std::cout << "pi = " << dPi << ", radToDeg = " << radToDeg << std::endl;

  const edm::ESTransientHandle<CSCGeometry> pDD = iSetup.getTransientHandle(ddToken_);
  std::cout << " Geometry node for CSCGeom is  " << &(*pDD) << std::endl;
  std::cout << " I have " << pDD->detTypes().size() << " detTypes" << std::endl;
  std::cout << " I have " << pDD->detUnits().size() << " detUnits" << std::endl;
  std::cout << " I have " << pDD->dets().size() << " dets" << std::endl;
  std::cout << " I have " << pDD->layers().size() << " layers" << std::endl;
  std::cout << " I have " << pDD->chambers().size() << " chambers" << std::endl;

  std::cout << " Ganged strips? " << pDD->gangedStrips() << std::endl;
  std::cout << " Wires only?    " << pDD->wiresOnly() << std::endl;
  std::cout << " Real wire geometry? " << pDD->realWireGeometry() << std::endl;
  std::cout << " Use offsets to coi? " << pDD->centreTIOffsets() << std::endl;

  std::cout << myName() << ": Begin iteration over geometry..." << std::endl;
  std::cout << "iter " << dashedLine_ << std::endl;

  std::cout << "\n  #     id(dec)      id(oct)                   "
               "  g(x=0)   g(y=0)   g(z=0)  g(z=-1)  g(z=+1)  Ns "
               "  phi(0)  phi(s1)  phi(sN)    dphi    dphi'      ds     off"
               "       uR       uL       lR       lL"
            << std::endl;

  int icount = 0;

  // Check the DetUnits
  for (const auto& it : pDD->detUnits()) {
    // Do we really have a CSC layer?

    auto layer = dynamic_cast<CSCLayer const*>(it);

    if (layer) {
      ++icount;

      // What's its DetId?

      DetId detId = layer->geographicalId();
      int id = detId();  // or detId.rawId()

      //	std::cout << "GeomDetUnit is of type " << detId.det() << " and raw id = " << id << std::endl;

      // There's going to be a lot of messing with field width (and precision) so
      // save input values...
      int iw = std::cout.width();      // save current width
      int ip = std::cout.precision();  // save current precision

      //	std::cout << "\n" << "Parameters of layer# " <<
      std::cout << std::setw(4) << icount << std::setw(12) << id << std::oct << std::setw(12) << id << std::dec
                << std::setw(iw) << "   E" << CSCDetId::endcap(id) << " S" << CSCDetId::station(id) << " R"
                << CSCDetId::ring(id) << " C" << std::setw(2) << CSCDetId::chamber(id) << std::setw(iw) << " L"
                << CSCDetId::layer(id);
      // CSCDetId::layer(id)  << " are:" << std::endl;

      // Check global coordinates of centre of CSCLayer, and how
      // local z direction relates to global z direction

      LocalPoint lCentre(0., 0., 0.);
      GlobalPoint gCentre = layer->toGlobal(lCentre);

      LocalPoint lCentre1(0., 0., -1.);
      GlobalPoint gCentre1 = layer->toGlobal(lCentre1);

      LocalPoint lCentre2(0., 0., 1.);
      GlobalPoint gCentre2 = layer->toGlobal(lCentre2);

      //		std::cout << "\nlocal(0,0,-1) = global " << gCentre1 << std::endl;
      //		std::cout <<   "local(0,0)    = global " << gCentre  << std::endl;
      //		std::cout <<   "local(0,0,+1) = global " << gCentre2 << std::endl;

      double gx = gCentre.x();
      double gy = gCentre.y();
      double gz = gCentre.z();
      double gz1 = gCentre1.z();
      double gz2 = gCentre2.z();
      if (fabs(gx) < 1.e-06)
        gx = 0.;
      if (fabs(gy) < 1.e-06)
        gy = 0.;
      if (fabs(gz) < 1.e-06)
        gz = 0.;
      if (fabs(gz1) < 1.e-06)
        gz1 = 0.;
      if (fabs(gz2) < 1.e-06)
        gz2 = 0.;

      int now = 9;
      int nop = 5;
      std::cout << std::setw(now) << std::setprecision(nop) << gx << std::setw(now) << std::setprecision(nop) << gy
                << std::setw(now) << std::setprecision(nop) << gz << std::setw(now) << std::setprecision(nop) << gz1
                << std::setw(now) << std::setprecision(nop) << gz2;

      // Global Phi of centre of CSCLayer

      //@@ CARE The following attempted conversion to degrees can be easily
      // subverted by GeometryVector/Phi.h enforcing its range convention!
      // Either a) use a separate local double before scaling...
      //        double cphi = gCentre.phi();
      //        double cphiDeg = cphi * radToDeg;
      // or b) use Phi's degree conversion...
      double cphiDeg = gCentre.phi().degrees();

      // I want to display in range 0 to 360
      if (fabs(cphiDeg) < 1.e-06) {
        cphiDeg = 0.;
      } else if (cphiDeg < 0.) {
        cphiDeg += 360.;
      } else if (cphiDeg >= 360.) {
        std::cout << "WARNING: resetting phi= " << cphiDeg << " to zero." << std::endl;
        cphiDeg = 0.;
      }

      // Handle some occasional ugly precision problems around zero
      if (fabs(cphiDeg) < 1.e-06)
        cphiDeg = 0.;
      //        int iphiDeg = static_cast<int>( cphiDeg );
      //	std::cout << "phi(0,0,0) = " << iphiDeg << " degrees" << std::endl;

      int nStrips = layer->geometry()->numberOfStrips();
      std::cout << std::setw(4) << nStrips;

      double cstrip1 = layer->centerOfStrip(1).phi();
      double cstripN = layer->centerOfStrip(nStrips).phi();

      double phiwid = layer->geometry()->stripPhiPitch();
      double stripwid = layer->geometry()->stripPitch();
      double stripoff = layer->geometry()->stripOffset();
      double phidif = fabs(cstrip1 - cstripN);

      // May have one strip at 180-epsilon and other at -180+epsilon
      // If so the raw difference is 360-(phi extent of chamber)
      // Want to reset that to (phi extent of chamber):
      if (phidif > dPi)
        phidif = fabs(phidif - 2. * dPi);
      double phiwid_check = phidif / double(nStrips - 1);

      // Clean up some stupid floating decimal aesthetics
      cstrip1 = cstrip1 * radToDeg;
      if (fabs(cstrip1) < 1.e-06)
        cstrip1 = 0.;
      else if (cstrip1 < 0.)
        cstrip1 += 360.;

      cstripN = cstripN * radToDeg;
      if (fabs(cstripN) < 1.e-06)
        cstripN = 0.;
      else if (cstripN < 0.)
        cstripN += 360.;

      if (fabs(stripoff) < 1.e-06)
        stripoff = 0.;

      now = 9;
      nop = 4;
      std::cout << std::setw(now) << std::setprecision(nop) << cphiDeg << std::setw(now) << std::setprecision(nop)
                << cstrip1 << std::setw(now) << std::setprecision(nop) << cstripN << std::setw(now)
                << std::setprecision(nop) << phiwid << std::setw(now) << std::setprecision(nop) << phiwid_check
                << std::setw(now) << std::setprecision(nop) << stripwid << std::setw(now) << std::setprecision(nop)
                << stripoff;
      //          << std::setw(8) << (layer->getOwner())->sector() ; //@@ No sector yet!

      // Layer geometry:  layer corner phi's...

      std::array<const float, 4> const& parameters = layer->geometry()->parameters();
      // these parameters are half-lengths, due to GEANT
      float hBottomEdge = parameters[0];
      float hTopEdge = parameters[1];
      float hThickness = parameters[2];
      float hApothem = parameters[3];

      //@@ CHECK
      //	std::cout << "hB=" << hBottomEdge << " hT=" << hTopEdge << " hTh=" << hThickness <<
      //	  " hAp=" << hApothem << std::endl;

      // first the face nearest the interaction
      // get the other face by using positive hThickness
      LocalPoint upperRightLocal(hTopEdge, hApothem, -hThickness);
      LocalPoint upperLeftLocal(-hTopEdge, hApothem, -hThickness);
      LocalPoint lowerRightLocal(hBottomEdge, -hApothem, -hThickness);
      LocalPoint lowerLeftLocal(-hBottomEdge, -hApothem, -hThickness);

      GlobalPoint upperRightGlobal = layer->toGlobal(upperRightLocal);
      GlobalPoint upperLeftGlobal = layer->toGlobal(upperLeftLocal);
      GlobalPoint lowerRightGlobal = layer->toGlobal(lowerRightLocal);
      GlobalPoint lowerLeftGlobal = layer->toGlobal(lowerLeftLocal);

      float uRGp = upperRightGlobal.phi().degrees();
      float uLGp = upperLeftGlobal.phi().degrees();
      float lRGp = lowerRightGlobal.phi().degrees();
      float lLGp = lowerLeftGlobal.phi().degrees();

      now = 9;
      std::cout << std::setw(now) << uRGp << std::setw(now) << uLGp << std::setw(now) << lRGp << std::setw(now) << lLGp
                << std::endl;

      // Reset the values we changed
      std::cout << std::setprecision(ip) << std::setw(iw);

      // Check idToDetUnit
      const GeomDetUnit* gdu = pDD->idToDetUnit(detId);
      assert(gdu == layer);
      // Check idToDet
      const GeomDet* gd = pDD->idToDet(detId);
      assert(gd == layer);
    } else {
      std::cout << "Could not dynamic_cast to a CSCLayer* " << std::endl;
    }
  }
  std::cout << dashedLine_ << " end" << std::endl;
}

//define this as a plug-in
DEFINE_FWK_MODULE(CSCGeometryAnalyzer);