Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:15:12

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalWaferCell
0004 // Class:      HGCalWaferCell
0005 //
0006 /**\class HGCalWaferCell HGCalWaferCell.cc
0007  test/HGCalWaferCell.cc
0008 
0009  Description: <one line class summary>
0010 
0011  Implementation:
0012      <Notes on implementation>
0013 */
0014 //
0015 // Original Author:  Sunanda Banerjee
0016 //         Created:  Mon 2020/09/15
0017 //
0018 //
0019 
0020 // system include files
0021 #include <fstream>
0022 #include <iostream>
0023 #include <memory>
0024 #include <string>
0025 #include <vector>
0026 
0027 // user include files
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 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0037 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0038 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0039 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0040 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0041 
0042 class HGCalWaferCell : public edm::one::EDAnalyzer<> {
0043 public:
0044   explicit HGCalWaferCell(const edm::ParameterSet&);
0045 
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 nameSense_, nameDetector_;
0054   const bool debug_;
0055   const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> geomToken_;
0056 };
0057 
0058 HGCalWaferCell::HGCalWaferCell(const edm::ParameterSet& iC)
0059     : nameSense_(iC.getParameter<std::string>("NameSense")),
0060       nameDetector_(iC.getParameter<std::string>("NameDevice")),
0061       debug_(iC.getParameter<bool>("Verbosity")),
0062       geomToken_(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_})) {
0063   edm::LogVerbatim("HGCalGeomX") << "Check # of cells for " << nameDetector_ << " using constants of " << nameSense_;
0064 }
0065 
0066 void HGCalWaferCell::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0067   edm::ParameterSetDescription desc;
0068   desc.add<std::string>("NameSense", "HGCalEESensitive");
0069   desc.add<std::string>("NameDevice", "HGCal EE");
0070   desc.add<bool>("Verbosity", false);
0071   descriptions.add("hgcalEEWaferCell", desc);
0072 }
0073 
0074 // ------------ method called to produce the data  ------------
0075 void HGCalWaferCell::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0076   const auto& geomR = iSetup.getData(geomToken_);
0077   const HGCalGeometry* geom = &geomR;
0078   const auto& hgdc = geom->topology().dddConstants();
0079 
0080   if (hgdc.waferHexagon8()) {
0081     // Find all valid wafers
0082     std::vector<DetId> ids = geom->getValidGeomDetIds();
0083     edm::LogVerbatim("HGCalGeomX") << "\n\nCheck Wafers for " << nameDetector_ << " with " << ids.size() << " wafers\n";
0084     DetId::Detector det = (nameSense_ == "HGCalHESiliconSensitive") ? DetId::HGCalHSi : DetId::HGCalEE;
0085     int bad(0);
0086     std::map<int, int> waferMap;
0087     for (unsigned int k = 0; k < hgdc.waferFileSize(); ++k) {
0088       int indx = hgdc.waferFileIndex(k);
0089       int layer = HGCalWaferIndex::waferLayer(indx);
0090       int waferU = HGCalWaferIndex::waferU(indx);
0091       int waferV = HGCalWaferIndex::waferV(indx);
0092       int type = std::get<0>(hgdc.waferFileInfo(k));
0093       int part = std::get<1>(hgdc.waferFileInfo(k));
0094       int rotn = std::get<2>(hgdc.waferFileInfo(k));
0095       HGCSiliconDetId id(det, 1, type, layer, waferU, waferV, 0, 0);
0096       int kndx = ((rotn * 10 + part) * 10 + type);
0097       if (debug_)
0098         edm::LogVerbatim("HGCalGeomX") << std::hex << id.rawId() << std::dec << " " << id << " Index:" << kndx;
0099       if (waferMap.find(kndx) == waferMap.end()) {
0100         waferMap[kndx] = 1;
0101       } else {
0102         ++waferMap[kndx];
0103       }
0104       if ((type < 0) || (type > 2) || (part < 0) || (part > 7) || (rotn < 0) || (rotn > 5))
0105         ++bad;
0106     }
0107     edm::LogVerbatim("HGCalGeomX") << bad << " wafers of unknown types among " << hgdc.waferFileSize() << " wafers\n";
0108 
0109     // Now print out the summary
0110     static const std::vector<int> itype = {0, 1, 2};
0111     static const std::vector<int> itypc = {0, 1, 2, 3, 4, 5};
0112     static const std::vector<int> itypp = {0, 1, 2, 3, 4, 5, 6, 7};
0113     static const std::vector<std::string> typep = {"F", "b", "g", "gm", "a", "d", "dm", "c"};
0114     for (const auto& type : itype) {
0115       for (unsigned int k = 0; k < itypp.size(); ++k) {
0116         int part = itypp[k];
0117         for (const auto& rotn : itypc) {
0118           int kndx = ((rotn * 10 + part) * 10 + type);
0119           auto itr = waferMap.find(kndx);
0120           if (itr != waferMap.end())
0121             edm::LogVerbatim("HGCalGeomX") << "Type:" << type << " Partial:" << typep[k] << " Orientation:" << rotn
0122                                            << " with " << (itr->second) << " wafers";
0123         }
0124       }
0125     }
0126 
0127     edm::LogVerbatim("HGCalGeomX") << "\n\nSummary of Cells\n================";
0128     for (const auto& type : itype) {
0129       int N = (type == 0) ? hgdc.getParameter()->nCellsFine_ : hgdc.getParameter()->nCellsCoarse_;
0130       for (unsigned int k = 0; k < itypp.size(); ++k) {
0131         int part = itypp[k];
0132         for (const auto& rotn : itypc) {
0133           int num(0);
0134           for (int cellU = 0; cellU < 2 * N; ++cellU) {
0135             for (int cellV = 0; cellV < 2 * N; ++cellV) {
0136               if (((cellV - cellU) < N) && (cellU - cellV) <= N) {
0137                 if (HGCalWaferMask::goodCell(cellU, cellV, N, part, rotn))
0138                   ++num;
0139               }
0140             }
0141           }
0142           edm::LogVerbatim("HGCalGeomX") << "Type:" << type << " Partial:" << typep[k] << " Orientation:" << rotn
0143                                          << " with " << num << " cells";
0144         }
0145       }
0146     }
0147     edm::LogVerbatim("HGCalGeomX") << "\n\n";
0148   }
0149 }
0150 
0151 // define this as a plug-in
0152 DEFINE_FWK_MODULE(HGCalWaferCell);