HcalGeometryTester

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 260 261 262 263 264 265
#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/HcalCommonData/interface/HcalDDDRecConstants.h"
#include "Geometry/Records/interface/CaloGeometryRecord.h"
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
#include "Geometry/HcalTowerAlgo/interface/HcalGeometry.h"
#include "Geometry/HcalTowerAlgo/interface/HcalFlexiHardcodeGeometryLoader.h"
#include "Geometry/HcalTowerAlgo/interface/HcalHardcodeGeometryLoader.h"
#include "Geometry/HcalTowerAlgo/interface/HcalTrigTowerGeometry.h"
#include "DataFormats/HcalDetId/interface/HcalTrigTowerDetId.h"
#include <iostream>
#include <string>

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

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

private:
  void testValidDetIds(CaloSubdetectorGeometry* geom,
                       const HcalTopology& topo,
                       DetId::Detector det,
                       int subdet,
                       const std::string& label);
  void testClosestCells(CaloSubdetectorGeometry* geom, const HcalTopology& top);
  void testClosestCell(const HcalDetId& detId, CaloSubdetectorGeometry* geom);
  void testTriggerGeometry(const HcalTopology& topology);
  void testFlexiValidDetIds(CaloSubdetectorGeometry* geom,
                            const HcalTopology& topology,
                            DetId::Detector det,
                            int subdet,
                            const std::string& label,
                            std::vector<int>& dins);
  void testFlexiGeomHF(CaloSubdetectorGeometry* geom);

  bool useOld_;
  edm::ESGetToken<HcalDDDRecConstants, HcalRecNumberingRecord> tok_ddrec_;
  edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tok_htopo_;
};

HcalGeometryTester::HcalGeometryTester(const edm::ParameterSet& iConfig) {
  useOld_ = iConfig.getParameter<bool>("UseOldLoader");
  tok_ddrec_ = esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord>();
  tok_htopo_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
}

void HcalGeometryTester::analyze(const edm::Event& /*iEvent*/, const edm::EventSetup& iSetup) {
  const HcalDDDRecConstants hcons = iSetup.getData(tok_ddrec_);
  const HcalTopology topology = iSetup.getData(tok_htopo_);
  CaloSubdetectorGeometry* geom(nullptr);
  if (useOld_) {
    HcalHardcodeGeometryLoader m_loader;
    geom = m_loader.load(topology);
  } else {
    HcalFlexiHardcodeGeometryLoader m_loader;
    geom = m_loader.load(topology, hcons);
  }

  testValidDetIds(geom, topology, DetId::Hcal, HcalBarrel, " BARREL ");
  testValidDetIds(geom, topology, DetId::Hcal, HcalEndcap, " ENDCAP ");
  testValidDetIds(geom, topology, DetId::Hcal, HcalOuter, " OUTER ");
  testValidDetIds(geom, topology, DetId::Hcal, HcalForward, " FORWARD ");

  testTriggerGeometry(topology);

  testClosestCells(geom, topology);

  edm::LogVerbatim("HCalGeom") << "HcalGeometryTester::Test SLHC Hcal Geometry";
  std::vector<int> dins;

  testFlexiValidDetIds(geom, topology, DetId::Hcal, HcalBarrel, " BARREL ", dins);
  testFlexiValidDetIds(geom, topology, DetId::Hcal, HcalEndcap, " ENDCAP ", dins);
  testFlexiValidDetIds(geom, topology, DetId::Hcal, HcalOuter, " OUTER ", dins);
  testFlexiValidDetIds(geom, topology, DetId::Hcal, HcalForward, " FORWARD ", dins);

  testFlexiGeomHF(geom);
}

void HcalGeometryTester::testValidDetIds(CaloSubdetectorGeometry* caloGeom,
                                         const HcalTopology& topology,
                                         DetId::Detector det,
                                         int subdet,
                                         const std::string& label) {
  std::stringstream s;
  s << label << " : \n";
  const std::vector<DetId>& idshb = caloGeom->getValidDetIds(det, subdet);

  int counter = 0;
  for (std::vector<DetId>::const_iterator i = idshb.begin(); i != idshb.end(); i++, ++counter) {
    HcalDetId hid = (*i);
    s << counter << ": din " << topology.detId2denseId(*i) << ":" << hid;
    auto cell = caloGeom->getGeometry(*i);
    s << *cell << std::endl;
  }

  s << "=== Total " << counter << " cells in " << label << std::endl;
  edm::LogVerbatim("HCalGeom") << s.str();
}

void HcalGeometryTester::testClosestCells(CaloSubdetectorGeometry* g, const HcalTopology& topology) {
  // make sure each cell is its own closest cell
  HcalDetId barrelDet1(HcalBarrel, 1, 1, 1);
  HcalDetId barrelDet2(HcalBarrel, 16, 50, 1);
  HcalDetId barrelDet3(HcalBarrel, -6, 20, 3);
  HcalDetId barrelDet4(HcalBarrel, -12, 33, 4);
  HcalDetId endcapDet1(HcalEndcap, -17, 72, 1);
  HcalDetId endcapDet2(HcalEndcap, 29, 35, 1);
  HcalDetId endcapDet3(HcalEndcap, 20, 51, 3);
  HcalDetId endcapDet4(HcalEndcap, -24, 63, 5);
  HcalDetId forwardDet1(HcalForward, 30, 71, 1);
  HcalDetId forwardDet3(HcalForward, -40, 71, 1);

  if (topology.valid(barrelDet1))
    testClosestCell(barrelDet1, g);
  if (topology.valid(barrelDet2))
    testClosestCell(barrelDet2, g);
  if (topology.valid(barrelDet3))
    testClosestCell(barrelDet3, g);
  if (topology.valid(barrelDet4))
    testClosestCell(barrelDet4, g);
  if (topology.valid(endcapDet1))
    testClosestCell(endcapDet1, g);
  if (topology.valid(endcapDet2))
    testClosestCell(endcapDet2, g);
  if (topology.valid(endcapDet3))
    testClosestCell(endcapDet3, g);
  if (topology.valid(endcapDet4))
    testClosestCell(endcapDet4, g);
  if (topology.valid(forwardDet1))
    testClosestCell(forwardDet1, g);
  if (topology.valid(forwardDet3))
    testClosestCell(forwardDet3, g);

  const std::vector<DetId>& idsb = g->getValidDetIds(DetId::Hcal, HcalBarrel);
  for (auto id : idsb) {
    testClosestCell(HcalDetId(id), g);
  }

  const std::vector<DetId>& idse = g->getValidDetIds(DetId::Hcal, HcalEndcap);
  for (auto id : idse) {
    testClosestCell(HcalDetId(id), g);
  }
}

void HcalGeometryTester::testClosestCell(const HcalDetId& detId, CaloSubdetectorGeometry* geom) {
  auto cell = geom->getGeometry(detId);
  HcalDetId closest = geom->getClosestCell(cell->getPosition());
  edm::LogVerbatim("HCalGeom") << "i/p " << detId << " position " << cell->getPosition() << " closest " << closest;

  if (closest != detId) {
    edm::LogVerbatim("HCalGeom") << "HcalGeometryTester::Mismatch.  Original HCAL cell is " << detId
                                 << " while nearest is " << closest << " ***ERROR***";
  }
}

void HcalGeometryTester::testTriggerGeometry(const HcalTopology& topology) {
  HcalTrigTowerGeometry trigTowers(&topology);
  edm::LogVerbatim("HCalGeom") << "HCAL trigger tower eta bounds:";
  for (int ieta = 1; ieta <= 32; ++ieta) {
    double eta1, eta2;
    trigTowers.towerEtaBounds(ieta, 0, eta1, eta2);
    edm::LogVerbatim("HCalGeom") << "[" << ieta << "] " << eta1 << " " << eta2;
  }

  // now test some cell mappings
  HcalDetId barrelDet(HcalBarrel, 1, 1, 1);
  HcalDetId endcapDet(HcalEndcap, 29, 1, 1);
  HcalDetId forwardDet1(HcalForward, 29, 71, 1);
  HcalDetId forwardDet2(HcalForward, 29, 71, 2);
  HcalDetId forwardDet3(HcalForward, 40, 71, 1);

  using TowerDets = std::vector<HcalTrigTowerDetId>;
  if (topology.valid(barrelDet)) {
    TowerDets barrelTowers = trigTowers.towerIds(barrelDet);
    edm::LogVerbatim("HCalGeom") << "Trigger Tower Size: Barrel " << barrelTowers.size();
    assert(barrelTowers.size() == 1);
    edm::LogVerbatim("HCalGeom") << "Tower[0] " << barrelTowers[0];
  }
  if (topology.valid(endcapDet)) {
    TowerDets endcapTowers = trigTowers.towerIds(endcapDet);
    edm::LogVerbatim("HCalGeom") << "Trigger Tower Size: Endcap " << endcapTowers.size();
    assert(!endcapTowers.empty());
    for (unsigned int k = 0; k < endcapTowers.size(); ++k)
      edm::LogVerbatim("HCalGeom") << "Tower[" << k << "] " << endcapTowers[k];
  }
  if (topology.valid(forwardDet1)) {
    TowerDets forwardTowers1 = trigTowers.towerIds(forwardDet1);
    edm::LogVerbatim("HCalGeom") << "Trigger Tower Size: Forward1 " << forwardTowers1.size();
    assert(!forwardTowers1.empty());
    for (unsigned int k = 0; k < forwardTowers1.size(); ++k)
      edm::LogVerbatim("HCalGeom") << "Tower[" << k << "] " << forwardTowers1[k];
  }
  if (topology.valid(forwardDet2)) {
    TowerDets forwardTowers2 = trigTowers.towerIds(forwardDet2);
    edm::LogVerbatim("HCalGeom") << "Trigger Tower Size: Forward2 " << forwardTowers2.size();
    assert(!forwardTowers2.empty());
    for (unsigned int k = 0; k < forwardTowers2.size(); ++k)
      edm::LogVerbatim("HCalGeom") << "Tower[" << k << "] " << forwardTowers2[k];
  }
  if (topology.valid(forwardDet3)) {
    TowerDets forwardTowers3 = trigTowers.towerIds(forwardDet3);
    edm::LogVerbatim("HCalGeom") << "Trigger Tower Size: Forward3 " << forwardTowers3.size();
    assert(!forwardTowers3.empty());
    for (unsigned int k = 0; k < forwardTowers3.size(); ++k)
      edm::LogVerbatim("HCalGeom") << "Tower[" << k << "] " << forwardTowers3[k];
  }
}

void HcalGeometryTester::testFlexiValidDetIds(CaloSubdetectorGeometry* caloGeom,
                                              const HcalTopology& topology,
                                              DetId::Detector det,
                                              int subdet,
                                              const std::string& label,
                                              std::vector<int>& dins) {
  std::stringstream s;
  s << label << " : " << std::endl;
  const std::vector<DetId>& idshb = caloGeom->getValidDetIds(det, subdet);

  int counter = 0;
  for (std::vector<DetId>::const_iterator i = idshb.begin(); i != idshb.end(); i++, ++counter) {
    HcalDetId hid = (*i);
    s << counter << ": din " << topology.detId2denseId(*i) << ":" << hid;
    dins.emplace_back(topology.detId2denseId(*i));

    auto cell = caloGeom->getGeometry(*i);
    s << *cell << std::endl;
  }

  std::sort(dins.begin(), dins.end());
  s << "=== Total " << counter << " cells in " << label << std::endl;

  counter = 0;
  for (std::vector<int>::const_iterator i = dins.begin(); i != dins.end(); ++i, ++counter) {
    HcalDetId hid = (topology.denseId2detId(*i));
    HcalDetId ihid = (topology.denseId2detId(dins[counter]));
    s << counter << ": din " << (*i) << " :" << hid << " == " << ihid << std::endl;
  }
  edm::LogVerbatim("HCalGeom") << s.str();
}

void HcalGeometryTester::testFlexiGeomHF(CaloSubdetectorGeometry* caloGeom) {
  std::stringstream s;
  s << "Test HF Geometry : " << std::endl;
  for (int ieta = 29; ieta <= 41; ++ieta) {
    HcalDetId cell3(HcalForward, ieta, 3, 1);
    auto cellGeometry3 = caloGeom->getGeometry(cell3);
    if (cellGeometry3) {
      s << "cell geometry iphi=3 -> ieta=" << ieta << " eta " << cellGeometry3->getPosition().eta() << "+-"
        << std::abs(cellGeometry3->getCorners()[0].eta() - cellGeometry3->getCorners()[2].eta()) / 2 << " phi "
        << cellGeometry3->getPosition().phi() / 3.1415 * 180 << "+-"
        << std::abs(cellGeometry3->getCorners()[0].phi() - cellGeometry3->getCorners()[2].phi()) / 3.1415 * 180 / 2.
        << std::endl;
    }
  }
  edm::LogVerbatim("HCalGeom") << s.str();
}

DEFINE_FWK_MODULE(HcalGeometryTester);