File indexing completed on 2024-12-08 23:44:09
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <fstream>
0022 #include <iostream>
0023 #include <sstream>
0024 #include <string>
0025 #include <vector>
0026
0027
0028 #include "FWCore/Framework/interface/Frameworkfwd.h"
0029 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/EventSetup.h"
0032 #include "FWCore/Framework/interface/MakerMacros.h"
0033 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0034 #include "FWCore/ParameterSet/interface/FileInPath.h"
0035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0036 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0037 #include "FWCore/Utilities/interface/transform.h"
0038
0039 #include "DataFormats/DetId/interface/DetId.h"
0040 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0041 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0042 #include "Geometry/HGCalCommonData/interface/HGCalGeomUtils.h"
0043 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0044
0045 class HGCalTestCellArea : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0046 public:
0047 explicit HGCalTestCellArea(const edm::ParameterSet &);
0048 ~HGCalTestCellArea() override = default;
0049 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0050
0051 void beginJob() override {}
0052 void beginRun(edm::Run const &, edm::EventSetup const &) override;
0053 void analyze(edm::Event const &iEvent, edm::EventSetup const &) override {}
0054 void endRun(edm::Run const &, edm::EventSetup const &) override {}
0055 void endJob() override {}
0056
0057 private:
0058 const std::vector<std::string> nameDetectors_;
0059 const std::string fileName_;
0060 const std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> tok_hgcal_;
0061 std::vector<const HGCalDDDConstants *> hgcCons_;
0062 std::vector<std::pair<DetId, uint32_t>> detIds_;
0063 };
0064
0065 HGCalTestCellArea::HGCalTestCellArea(const edm::ParameterSet &iC)
0066 : nameDetectors_(iC.getParameter<std::vector<std::string>>("nameDetectors")),
0067 fileName_(iC.getParameter<std::string>("fileName")),
0068 tok_hgcal_{edm::vector_transform(nameDetectors_, [this](const std::string &name) {
0069 return esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name});
0070 })} {
0071 std::ostringstream st1;
0072 for (const auto &name : nameDetectors_)
0073 st1 << " : " << name;
0074 edm::LogVerbatim("HGCGeom") << "Test validity of cells for " << nameDetectors_.size() << " detectors" << st1.str()
0075 << " with inputs from " << fileName_;
0076 if (!fileName_.empty()) {
0077 edm::FileInPath filetmp("Geometry/HGCalCommonData/data/" + fileName_);
0078 std::string fileName = filetmp.fullPath();
0079 std::ifstream fInput(fileName.c_str());
0080 if (!fInput.good()) {
0081 edm::LogVerbatim("HGCGeom") << "Cannot open file " << fileName;
0082 } else {
0083 char buffer[80];
0084 const std::vector<DetId::Detector> dets = {DetId::HGCalEE, DetId::HGCalHSi, DetId::HGCalHSc};
0085 while (fInput.getline(buffer, 80)) {
0086 std::vector<std::string> items = HGCalGeomUtils::splitString(std::string(buffer));
0087 if (items.size() == 8) {
0088 DetId::Detector det = static_cast<DetId::Detector>(std::atoi(items[0].c_str()));
0089 auto itr = std::find(dets.begin(), dets.end(), det);
0090 if (itr != dets.end()) {
0091 uint32_t pos = static_cast<uint32_t>(itr - dets.begin());
0092 DetId id(0);
0093 if ((det == DetId::HGCalEE) || (det == DetId::HGCalHSi)) {
0094 int type = std::atoi(items[1].c_str());
0095 int zside = std::atoi(items[2].c_str());
0096 int layer = std::atoi(items[3].c_str());
0097 int waferU = std::atoi(items[4].c_str());
0098 int waferV = std::atoi(items[5].c_str());
0099 int cellU = std::atoi(items[6].c_str());
0100 int cellV = std::atoi(items[7].c_str());
0101 id = static_cast<DetId>(HGCSiliconDetId(det, zside, type, layer, waferU, waferV, cellU, cellV));
0102 detIds_.emplace_back(id, pos);
0103 }
0104 }
0105 }
0106 }
0107 fInput.close();
0108 }
0109 }
0110 edm::LogVerbatim("HGCGeom") << "Reads " << detIds_.size() << " ID's from " << fileName_;
0111 for (unsigned int k = 0; k < detIds_.size(); ++k) {
0112 edm::LogVerbatim("HGCGeom") << "[" << k << "] " << HGCSiliconDetId(detIds_[k].first) << " from DDConstant "
0113 << (detIds_[k].second);
0114 }
0115 }
0116
0117 void HGCalTestCellArea::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0118 std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive"};
0119 edm::ParameterSetDescription desc;
0120 desc.add<std::vector<std::string>>("nameDetectors", names);
0121 desc.add<std::string>("fileName", "missD88.txt");
0122 descriptions.add("hgcalTestCellArea", desc);
0123 }
0124
0125
0126 void HGCalTestCellArea::beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) {
0127
0128 std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive"};
0129 std::vector<DetId::Detector> dets = {DetId::HGCalEE, DetId::HGCalHSi};
0130 std::map<DetId::Detector, uint32_t> detMap;
0131 for (uint32_t i = 0; i < nameDetectors_.size(); i++) {
0132 edm::LogVerbatim("HGCGeom") << "Tries to initialize HGCalGeometry and HGCalDDDConstants for " << i << ":"
0133 << nameDetectors_[i];
0134 const edm::ESHandle<HGCalDDDConstants> &hgcCons = iSetup.getHandle(tok_hgcal_[i]);
0135 if (hgcCons.isValid()) {
0136 hgcCons_.push_back(hgcCons.product());
0137 } else {
0138 edm::LogWarning("HGCGeom") << "Cannot initiate HGCalDDDConstants for " << nameDetectors_[i] << std::endl;
0139 }
0140 auto ii = std::find(names.begin(), names.end(), nameDetectors_[i]);
0141 if (ii != names.end()) {
0142 uint32_t k = static_cast<uint32_t>(ii - names.begin());
0143 detMap[dets[k]] = i;
0144 }
0145 }
0146 edm::LogVerbatim("HGCGeom") << "Loaded HGCalDDConstants for " << detMap.size() << " detectors";
0147
0148 for (auto itr = detMap.begin(); itr != detMap.end(); ++itr)
0149 edm::LogVerbatim("HGCGeom") << "[" << itr->second << "]: " << nameDetectors_[itr->second] << " for Detector "
0150 << itr->first;
0151
0152 for (unsigned int k = 0; k < detIds_.size(); ++k) {
0153 const HGCalDDDConstants *cons = hgcCons_[detMap[(detIds_[k].first).det()]];
0154 HGCSiliconDetId id(detIds_[k].first);
0155 edm::LogVerbatim("HGCGeom") << "Hit[" << k << "] " << id << " Area " << cons->cellArea(id, false) << " Valid "
0156 << cons->isValidHex8(
0157 id.layer(), id.waferU(), id.waferV(), id.cellU(), id.cellV(), true);
0158 }
0159 }
0160
0161
0162 DEFINE_FWK_MODULE(HGCalTestCellArea);