File indexing completed on 2025-05-09 22:37:45
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021 #include <algorithm>
0022 #include <iostream>
0023 #include <memory>
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/ParameterSet.h"
0035 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0036
0037 #include "DataFormats/DetId/interface/DetId.h"
0038 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0039 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0040 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0041
0042 class HGCalListValidCells : public edm::one::EDAnalyzer<> {
0043 public:
0044 explicit HGCalListValidCells(const edm::ParameterSet&);
0045 ~HGCalListValidCells() override = default;
0046 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0047
0048 void beginJob() override {}
0049 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0050 void endJob() override {}
0051
0052 private:
0053 const std::string name_;
0054 const int partialType_, verbosity_;
0055 const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> geomToken_;
0056 };
0057
0058 HGCalListValidCells::HGCalListValidCells(const edm::ParameterSet& iC)
0059 : name_(iC.getParameter<std::string>("detector")),
0060 partialType_(iC.getParameter<int>("partialType")),
0061 verbosity_(iC.getParameter<int>("verbosity")),
0062 geomToken_(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag("", name_))) {}
0063
0064 void HGCalListValidCells::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0065 edm::ParameterSetDescription desc;
0066 desc.add<std::string>("detector", "HGCalEESensitive");
0067 desc.add<int>("partialType", 16);
0068 desc.add<int>("verbosity", 0);
0069 descriptions.add("hgcalListValidCellsEE", desc);
0070 }
0071
0072
0073 void HGCalListValidCells::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0074 const auto& geom = &iSetup.getData(geomToken_);
0075 DetId::Detector det = (name_ == "HGCalHESiliconSensitive") ? DetId::HGCalHSi : DetId::HGCalEE;
0076 edm::LogVerbatim("HGCalGeom") << "Perform test for " << name_ << " Detector " << det;
0077
0078 std::string parts[26] = {"Full", "Five", "ChopTwo", "ChopTwoM", "Half", "Semi", "Semi2", "Three", "Half2",
0079 "Five2", "????", "LDTop", "LDBottom", "LDLeft", "LDRight", "LDFive", "LDThree", "????",
0080 "????", "????", "????", "HDTop", "HDBottom", "HDLeft", "HDRight", "HDFive"};
0081 const std::vector<DetId>& ids = geom->getValidDetIds();
0082 edm::LogVerbatim("HGCalGeom") << "Find the list of valid DetIds of type " << partialType_ << ":"
0083 << parts[partialType_] << " among a list of " << ids.size() << " valid ids of "
0084 << geom->cellElement();
0085 std::vector<HGCSiliconDetId> detIds;
0086
0087 for (auto const& id : ids) {
0088 if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
0089 HGCSiliconDetId detId(id);
0090 HGCalParameters::waferInfo info =
0091 geom->topology().dddConstants().waferInfo(detId.layer(), detId.waferU(), detId.waferV());
0092 if (info.part == partialType_) {
0093 if (std::find(detIds.begin(), detIds.end(), detId) == detIds.end())
0094 detIds.emplace_back(detId);
0095 }
0096 } else {
0097 edm::LogVerbatim("HGCalGeom") << "Illegal Det " << id.det() << " in " << std::hex << id.rawId() << std::dec
0098 << " ERROR";
0099 }
0100 }
0101 edm::LogVerbatim("HGCalGeom") << "There are " << detIds.size() << " valid Ids with partial type " << partialType_
0102 << ":" << parts[partialType_];
0103 if (verbosity_ > 0) {
0104 for (auto const detId : detIds)
0105 edm::LogVerbatim("HGCalGeom") << " " << detId;
0106 }
0107
0108 if (detIds.size() > 0) {
0109 std::vector<int> cellPatterns;
0110 for (auto const& detId : detIds) {
0111 int iuv = (100 * detId.cellU() + detId.cellV());
0112 if (std::find(cellPatterns.begin(), cellPatterns.end(), iuv) == cellPatterns.end())
0113 cellPatterns.emplace_back(iuv);
0114 }
0115 std::sort(cellPatterns.begin(), cellPatterns.end());
0116 edm::LogVerbatim("HGCalGeom") << "There are " << cellPatterns.size() << " different cell patterns:";
0117 for (const auto& iuv : cellPatterns) {
0118 int u = ((iuv / 100) % 100);
0119 int v = (iuv % 100);
0120 edm::LogVerbatim("HGCalGeom") << "u = " << std::setw(3) << u << " v = " << std::setw(3) << v;
0121 }
0122 }
0123 }
0124
0125
0126 DEFINE_FWK_MODULE(HGCalListValidCells);