Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 13:03:42

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalWaferInFileOrientation
0004 // Class:      HGCalWaferInFileOrientation
0005 //
0006 /**\class HGCalWaferInFileOrientation HGCalWaferInFileOrientation.cc
0007  test/HGCalWaferInFileOrientation.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 2021/10/07
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 
0031 #include "FWCore/Framework/interface/Event.h"
0032 #include "FWCore/Framework/interface/EventSetup.h"
0033 #include "FWCore/Framework/interface/MakerMacros.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0036 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0037 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0038 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0039 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0040 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0041 
0042 class HGCalWaferInFileOrientation : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0043 public:
0044   explicit HGCalWaferInFileOrientation(const edm::ParameterSet&);
0045 
0046   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0047 
0048   void beginRun(edm::Run const&, edm::EventSetup const&) override;
0049   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override {}
0050   void endRun(edm::Run const&, edm::EventSetup const&) override {}
0051 
0052 private:
0053   const std::string nameDetector_;
0054   const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> geomToken_;
0055   const std::vector<int> layers_;
0056   const std::vector<int> waferU_, waferV_, types_;
0057 };
0058 
0059 HGCalWaferInFileOrientation::HGCalWaferInFileOrientation(const edm::ParameterSet& iC)
0060     : nameDetector_(iC.getParameter<std::string>("detectorName")),
0061       geomToken_(esConsumes<HGCalGeometry, IdealGeometryRecord, edm::Transition::BeginRun>(
0062           edm::ESInputTag{"", nameDetector_})),
0063       layers_(iC.getParameter<std::vector<int>>("layers")),
0064       waferU_(iC.getParameter<std::vector<int>>("waferUs")),
0065       waferV_(iC.getParameter<std::vector<int>>("waferVs")),
0066       types_(iC.getParameter<std::vector<int>>("types")) {
0067   edm::LogVerbatim("HGCalGeom") << "Test Orientation for " << layers_.size() << " ID's in " << nameDetector_;
0068 }
0069 
0070 void HGCalWaferInFileOrientation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0071   edm::ParameterSetDescription desc;
0072   std::vector<int> layer = {33, 33, 33, 33, 33, 33};
0073   std::vector<int> waferU = {12, 12, -3, -3, -9, -9};
0074   std::vector<int> waferV = {3, 9, 9, -12, 3, -12};
0075   std::vector<int> type = {2, 2, 2, 2, 2, 2};
0076   desc.add<std::string>("detectorName", "HGCalHESiliconSensitive");
0077   desc.add<std::vector<int>>("layers", layer);
0078   desc.add<std::vector<int>>("waferUs", waferU);
0079   desc.add<std::vector<int>>("waferVs", waferV);
0080   desc.add<std::vector<int>>("types", type);
0081   descriptions.add("hgcalWaferInFileOrientation", desc);
0082 }
0083 
0084 void HGCalWaferInFileOrientation::beginRun(edm::Run const&, edm::EventSetup const& iSetup) {
0085   const auto& geomR = iSetup.getData(geomToken_);
0086   const HGCalGeometry* geom = &geomR;
0087   const auto& hgdc = geom->topology().dddConstants();
0088 
0089   edm::LogVerbatim("HGCalGeom") << "Check Wafer orientations in file for " << nameDetector_ << "\n";
0090   if (hgdc.waferHexagon8()) {
0091     DetId::Detector det = (nameDetector_ == "HGCalHESiliconSensitive") ? DetId::HGCalHSi : DetId::HGCalEE;
0092     static std::vector<std::string> types = {"F", "b", "g", "gm", "a", "d", "dm", "c", "am", "bm", "X"};
0093     int allG(0), badP(0), badR(0), badG(0);
0094     int layerOff = hgdc.getLayerOffset();
0095     for (unsigned int k = 0; k < layers_.size(); ++k) {
0096       int layer = layers_[k] - layerOff;
0097       int indx = HGCalWaferIndex::waferIndex(layer, waferU_[k], waferV_[k]);
0098       int part1 = std::get<1>(hgdc.waferFileInfoFromIndex(indx));
0099       int rotn1 = std::get<2>(hgdc.waferFileInfoFromIndex(indx));
0100       HGCSiliconDetId id(det, 1, types_[k], layer, waferU_[k], waferV_[k], 0, 0);
0101       if (geom->topology().validModule(id, 3)) {
0102         ++allG;
0103         int part2 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), false, false).first;
0104         int rotn2 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), false, false).second;
0105         int part3 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), true, false).first;
0106         int rotn3 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), true, false).second;
0107         bool partOK = (part1 == part2);
0108         bool rotnOK = (rotn1 == rotn2);
0109         if (!partOK)
0110           ++badP;
0111         if (!rotnOK)
0112           ++badR;
0113         if ((!partOK) || (!rotnOK)) {
0114           ++badG;
0115           std::string partx1 = (part1 < static_cast<int>(types.size())) ? types[part1] : "X";
0116           std::string partx2 = (part2 < static_cast<int>(types.size())) ? types[part2] : "X";
0117           std::string partx3 = (part3 < static_cast<int>(types.size())) ? types[part3] : "X";
0118           const auto& xy = hgdc.waferPosition(layer, waferU_[k], waferV_[k], true, false);
0119           edm::LogVerbatim("HGCalGeom") << "ID[" << k << "]: " << id << " (" << partx1 << ":" << partx2 << ":" << partx3
0120                                         << ", " << rotn1 << ":" << rotn2 << ":" << rotn3 << ") at ("
0121                                         << std::setprecision(4) << xy.first << ", " << xy.second << ", "
0122                                         << hgdc.waferZ(layer, true) << ") failure flag " << partOK << ":" << rotnOK;
0123         }
0124       }
0125     }
0126     edm::LogVerbatim("HGCalGeom") << "\n\nFinds " << badG << " (" << badP << ":" << badR << ") mismatch among " << allG
0127                                   << ":" << layers_.size() << " wafers\n";
0128   }
0129 }
0130 
0131 // define this as a plug-in
0132 DEFINE_FWK_MODULE(HGCalWaferInFileOrientation);