File indexing completed on 2025-07-17 22:23:04
0001 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0002 #include "FWCore/Framework/interface/EventSetup.h"
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0008 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0009 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0010 #include "Geometry/HGCalCommonData/interface/HGCalWaferType.h"
0011 #include "DataFormats/ForwardDetId/interface/HFNoseDetId.h"
0012 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0013 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0014 #include <iostream>
0015
0016 class HGCalWaferSimWt : public edm::one::EDAnalyzer<> {
0017 public:
0018 explicit HGCalWaferSimWt(const edm::ParameterSet&);
0019 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0020
0021 void beginJob() override {}
0022 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0023 void endJob() override {}
0024
0025 private:
0026 const std::vector<std::string> names_;
0027 std::vector<edm::ESGetToken<HGCalGeometry, IdealGeometryRecord>> geomTokens_;
0028 };
0029
0030 HGCalWaferSimWt::HGCalWaferSimWt(const edm::ParameterSet& iC)
0031 : names_(iC.getParameter<std::vector<std::string>>("detectorNames")) {
0032 for (unsigned int k = 0; k < names_.size(); ++k) {
0033 edm::LogVerbatim("HGCalGeomX") << "Study detector [" << k << "] " << names_[k] << std::endl;
0034 geomTokens_.emplace_back(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag{"", names_[k]}));
0035 }
0036 }
0037
0038 void HGCalWaferSimWt::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039 edm::ParameterSetDescription desc;
0040 std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive"};
0041 desc.add<std::vector<std::string>>("detectorNames", names);
0042 descriptions.add("hgcalWaferSimWt", desc);
0043 }
0044
0045 void HGCalWaferSimWt::analyze(const edm::Event& , const edm::EventSetup& iSetup) {
0046 std::vector<std::string> parType = {"Full", "Five", "ChopTwo", "ChopTwoM", "Half", "Semi", "Semi2",
0047 "Three", "Half2", "Five2", "JK10", "LDTop", "LDBottom", "LDLeft",
0048 "LDRight", "LDFive", "LDThree", "JK17", "JK18", "JK19", "JK20",
0049 "HDTop", "HDBottom", "HDLeft", "HDRight", "HDFive"};
0050 std::vector<std::string> detType = {"HD120", "LD200", "LD300", "HD200"};
0051 for (unsigned int k = 0; k < names_.size(); ++k) {
0052 const auto& geomR = iSetup.getData(geomTokens_[k]);
0053 const HGCalGeometry* geom = &geomR;
0054 const std::vector<DetId>& ids = geom->getValidDetIds();
0055 edm::LogVerbatim("HGCalGeomX") << ids.size() << " valid Ids for detector " << names_[k] << std::endl;
0056 int nall(0), ntypes(0);
0057 std::vector<int> idxs;
0058 for (auto id : ids) {
0059 ++nall;
0060 auto cell = geom->getGeometry(id);
0061 HGCSiliconDetId hid(id);
0062 int type = hid.type();
0063 int part = geom->topology().dddConstants().partialWaferType(hid.layer(), hid.waferU(), hid.waferV());
0064 int idx = part * 10 + type;
0065 if (std::find(idxs.begin(), idxs.end(), idx) == idxs.end()) {
0066 ++ntypes;
0067 idxs.push_back(idx);
0068 double xpos = 10.0 * (cell->getPosition().x());
0069 double ypos = 10.0 * (cell->getPosition().y());
0070 int waferU, waferV, cellU, cellV, cellType;
0071 double wt;
0072 geom->topology().dddConstants().waferFromPosition(
0073 xpos, ypos, hid.zside(), hid.layer(), waferU, waferV, cellU, cellV, cellType, wt, false, true);
0074 std::string stype = (type >= 0 && type <= 3) ? detType[type] : ("JK" + std::to_string(type));
0075 std::string spart = (part >= 0 && part <= 25) ? parType[part] : ("JK" + std::to_string(part));
0076 int index = HGCalWaferIndex::waferIndex(hid.layer(), waferU, waferV);
0077 int celltypeX = HGCalWaferType::getType(index, geom->topology().dddConstants().getParameter()->waferInfoMap_);
0078 edm::LogVerbatim("HGCalGeomX") << "[" << ntypes << "] " << stype << " " << spart << " wt " << wt << " for "
0079 << hid << " at " << xpos << ":" << ypos << " Wafer " << waferU << ":" << waferV
0080 << " cell " << cellU << ":" << cellV << " Type " << cellType << ":" << celltypeX;
0081 }
0082 }
0083 edm::LogVerbatim("HGCalGeomX") << "\n\nFinds " << idxs.size() << " different wafer types among " << nall
0084 << " cells of the detector\n";
0085 }
0086 }
0087
0088 DEFINE_FWK_MODULE(HGCalWaferSimWt);