CaloTowerGeometryAnalyzer

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
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "Geometry/Records/interface/CaloGeometryRecord.h"
#include "Geometry/CaloGeometry/interface/CaloGenericDetId.h"
#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
#include "DataFormats/HcalDetId/interface/HcalDetId.h"

#include <fstream>
#include <iostream>
#include <cmath>

namespace {
  bool AreSame(double a, double b, double epsilon) { return std::fabs(a - b) < epsilon; }

  /*
  inline double
  round( double n, int digits )
  {
    double mult = pow( 10, digits );
    return floor( n * mult ) / mult;
  }
  */

  std::map<int, std::string> hcalmapping = {
      {1, "HB"}, {2, "HE"}, {3, "HO"}, {4, "HF"}, {5, "HT"}, {6, "ZDC"}, {0, "Empty"}};
}  // namespace

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

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

private:
  std::string m_fname;
  double m_epsilon;
  edm::ESGetToken<CaloGeometry, CaloGeometryRecord> tok_geom_;
};

CaloTowerGeometryAnalyzer::CaloTowerGeometryAnalyzer(const edm::ParameterSet& iConfig)
    : m_fname("CaloTower.cells"), m_epsilon(0.004) {
  m_fname = iConfig.getParameter<std::string>("FileName");
  m_epsilon = iConfig.getParameter<double>("Epsilon");
  tok_geom_ = esConsumes<CaloGeometry, CaloGeometryRecord>();
}

CaloTowerGeometryAnalyzer::~CaloTowerGeometryAnalyzer(void) {}

void CaloTowerGeometryAnalyzer::analyze(const edm::Event& /*iEvent*/, const edm::EventSetup& iSetup) {
  std::fstream fAll(std::string(m_fname + ".all").c_str(), std::ios_base::out);
  std::fstream f(std::string(m_fname + ".diff").c_str(), std::ios_base::out);

  const CaloGeometry* caloGeom = &iSetup.getData(tok_geom_);

  const std::vector<DetId>& dha(caloGeom->getSubdetectorGeometry(DetId::Hcal, 1)->getValidDetIds());

  const std::vector<DetId>& ids(caloGeom->getValidDetIds());

  fAll << std::setw(4) << "iEta" << std::setw(4) << "iPhi" << std::setw(4) << "Z" << std::setw(10) << "Sub:Depth"
       << std::setw(6) << "Eta" << std::setw(13) << "Phi" << std::setw(18) << "Corners in Eta\n";

  f << std::setw(4) << "iEta" << std::setw(4) << "iPhi" << std::setw(4) << "Z" << std::setw(10) << "Sub:Depth"
    << std::setw(6) << "Eta" << std::setw(13) << "Phi" << std::setw(18) << "Corners in Eta\n";

  for (auto id : ids) {
    const CaloGenericDetId cgid(id);
    if (cgid.isCaloTower()) {
      const CaloTowerDetId cid(id);
      const int ie(cid.ieta());
      const int ip(cid.iphi());
      const int iz(cid.zside());

      fAll << std::setw(4) << ie << std::setw(4) << ip << std::setw(4) << iz << std::setw(6) << "-";

      auto cell = caloGeom->getGeometry(id);
      assert(cell);
      const GlobalPoint& pos = cell->getPosition();
      double eta = pos.eta();
      double phi = pos.phi();
      fAll << std::setw(10) << eta << std::setw(13) << phi;

      const CaloCellGeometry::CornersVec& corners(cell->getCorners());
      for (unsigned int i(0); i != corners.size(); ++i) {
        const GlobalPoint& cpos = corners[i];
        double ceta = cpos.eta();

        fAll << std::setw(13) << ceta;
      }
      fAll << "\n";

      for (auto ii : dha) {
        const HcalDetId hid(ii);
        const int iie(hid.ieta());
        const int iip(hid.iphi());
        const int iiz(hid.zside());
        const int iid(hid.depth());

        if (ie == iie && ip == iip && iz == iiz) {
          fAll << std::setw(4) << iie << std::setw(4) << iip << std::setw(4) << iiz;

          std::map<int, std::string>::const_iterator iter = hcalmapping.find(hid.subdet());
          if (iter != hcalmapping.end()) {
            fAll << std::setw(4) << iter->second.c_str() << ":" << iid;
          }

          auto hcell = caloGeom->getGeometry(ii);
          assert(hcell);
          const GlobalPoint& hpos = hcell->getPosition();
          double heta = hpos.eta();
          double hphi = hpos.phi();

          fAll << std::setw(10) << heta << std::setw(13) << hphi;

          const CaloCellGeometry::CornersVec& hcorners(hcell->getCorners());
          for (unsigned int i(0); i != hcorners.size(); ++i) {
            const GlobalPoint& hcpos = hcorners[i];
            double hceta = hcpos.eta();

            fAll << std::setw(13) << hceta;
          }

          if (!AreSame(eta, heta, m_epsilon)) {
            fAll << "*DIFFER in Eta*";

            f << std::setw(4) << ie << std::setw(4) << ip << std::setw(4) << iz << std::setw(6) << "-";

            f << std::setw(10) << eta << std::setw(13) << phi;

            for (unsigned int i(0); i != corners.size(); ++i) {
              const GlobalPoint& cpos = corners[i];
              double ceta = cpos.eta();

              f << std::setw(9) << ceta;
            }
            f << "\n";

            f << std::setw(4) << iie << std::setw(4) << iip << std::setw(4) << iiz;

            if (iter != hcalmapping.end()) {
              f << std::setw(4) << iter->second.c_str() << ":" << iid;
            }

            f << std::setw(10) << heta << std::setw(13) << hphi;

            for (unsigned int i(0); i != hcorners.size(); ++i) {
              const GlobalPoint& hcpos = hcorners[i];
              double hceta = hcpos.eta();

              f << std::setw(9) << hceta;
            }
            f << "\n\n";
          }

          if (!AreSame(phi, hphi, m_epsilon))
            fAll << " *DIFFER in Phi*";
          fAll << "\n";
        }
      }
    }
  }

  fAll.close();
  f.close();
}

DEFINE_FWK_MODULE(CaloTowerGeometryAnalyzer);