Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-29 02:55:25

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/FileInPath.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "FWCore/Utilities/interface/InputTag.h"
0013 #include "FWCore/Utilities/interface/Exception.h"
0014 
0015 #include "DataFormats/ForwardDetId/interface/HGCSiliconDetId.h"
0016 #include "DataFormats/ForwardDetId/interface/HGCScintillatorDetId.h"
0017 
0018 #include "SimDataFormats/CaloHit/interface/PCaloHit.h"
0019 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
0020 #include "SimG4CMS/Calo/interface/CaloSimUtils.h"
0021 #include "SimG4CMS/Calo/interface/HGCGuardRing.h"
0022 
0023 #include "Geometry/CaloTopology/interface/HGCalTopology.h"
0024 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0025 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0026 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0027 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0028 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0029 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0030 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0031 
0032 #include <fstream>
0033 #include <sstream>
0034 #include <string>
0035 #include <map>
0036 
0037 class HGCalTestGuardRing : public edm::one::EDAnalyzer<> {
0038 public:
0039   HGCalTestGuardRing(const edm::ParameterSet& ps);
0040   ~HGCalTestGuardRing() override = default;
0041   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0042 
0043 protected:
0044   void analyze(edm::Event const&, edm::EventSetup const&) override;
0045   void beginJob() override {}
0046   void endJob() override {}
0047 
0048 private:
0049   const std::string nameSense_, waferFile_;
0050   const double guardRingOffset_;
0051   const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> geomToken_;
0052   std::map<HGCSiliconDetId, int> waferID_;
0053 };
0054 
0055 HGCalTestGuardRing::HGCalTestGuardRing(const edm::ParameterSet& ps)
0056     : nameSense_(ps.getParameter<std::string>("nameSense")),
0057       waferFile_(ps.getParameter<std::string>("waferFile")),
0058       guardRingOffset_(ps.getParameter<double>("guardRingOffset")),
0059       geomToken_(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_})) {
0060   DetId::Detector det = (nameSense_ != "HGCalHESiliconSensitive") ? DetId::HGCalEE : DetId::HGCalHSi;
0061   edm::LogVerbatim("HGCalSim") << "Test Guard Ring Offset " << guardRingOffset_ << " for " << nameSense_ << ":" << det
0062                                << " for wafers read from file " << waferFile_;
0063   if (!waferFile_.empty()) {
0064     std::string thick[4] = {"h120", "l200", "l300", "h200"};
0065     int addType[4] = {HGCalTypes::WaferFineThin,
0066                       HGCalTypes::WaferCoarseThin,
0067                       HGCalTypes::WaferCoarseThick,
0068                       HGCalTypes::WaferFineThick};
0069     const int partTypeH[6] = {HGCalTypes::WaferFull,
0070                               HGCalTypes::WaferHalf2,
0071                               HGCalTypes::WaferChopTwoM,
0072                               HGCalTypes::WaferSemi2,
0073                               HGCalTypes::WaferSemi2,
0074                               HGCalTypes::WaferFive2};
0075     const int partTypeL[7] = {HGCalTypes::WaferFull,
0076                               HGCalTypes::WaferHalf,
0077                               HGCalTypes::WaferHalf,
0078                               HGCalTypes::WaferSemi,
0079                               HGCalTypes::WaferSemi,
0080                               HGCalTypes::WaferFive,
0081                               HGCalTypes::WaferThree};
0082     edm::FileInPath filetmp("SimG4CMS/Calo/data/" + waferFile_);
0083     std::string fileName = filetmp.fullPath();
0084     std::ifstream fInput(fileName.c_str());
0085     if (!fInput.good()) {
0086       edm::LogVerbatim("HGCalSim") << "Cannot open file " << fileName;
0087     } else {
0088       char buffer[80];
0089       while (fInput.getline(buffer, 80)) {
0090         std::vector<std::string> items = CaloSimUtils::splitString(std::string(buffer));
0091         if (items.size() > 6) {
0092           int layer = std::atoi(items[0].c_str());
0093           int waferU = std::atoi(items[4].c_str());
0094           int waferV = std::atoi(items[5].c_str());
0095           int thck = static_cast<int>(std::find(thick, thick + 4, items[2]) - thick);
0096           int type = (thck < 4) ? addType[thck] : 0;
0097           HGCSiliconDetId id(det, -1, type, layer, waferU, waferV, 0, 0);
0098           int orient = std::atoi(items[5].c_str());
0099           int part = std::atoi(items[1].c_str());
0100           if (part >= 0) {
0101             if (type == HGCalTypes::WaferFineThin)
0102               part = partTypeH[part];
0103             else
0104               part = partTypeL[part];
0105           }
0106           waferID_[id] = orient * 100 + part;
0107 #ifdef EDM_ML_DEBUG
0108           edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing::Reads " << id << " Orientation:Partial " << orient << ":"
0109                                        << part;
0110 #endif
0111         }
0112       }
0113       edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing::Reads in " << waferID_.size() << " wafer information from "
0114                                    << fileName;
0115       fInput.close();
0116     }
0117   }
0118 }
0119 
0120 void HGCalTestGuardRing::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0121   edm::ParameterSetDescription desc;
0122   desc.add<std::string>("nameSense", "HGCalEESensitive");
0123   desc.add<std::string>("waferFile", "testWafersEE.txt");
0124   desc.add<double>("guardRingOffset", 1.0);
0125   descriptions.add("hgcalTestGuardRingEE", desc);
0126 }
0127 
0128 void HGCalTestGuardRing::analyze(const edm::Event& e, const edm::EventSetup& iS) {
0129   // get HGCalGeometry
0130   const HGCalGeometry* geom = &iS.getData(geomToken_);
0131   const HGCalDDDConstants& hgc = geom->topology().dddConstants();
0132   double waferSize = hgc.waferSize(false);
0133   HGCalCell wafer(waferSize, hgc.getUVMax(0), hgc.getUVMax(1));
0134   // get the hit collection
0135   edm::LogVerbatim("HGCalSim") << "HGCalTestGuardRing: Wafer Szie " << waferSize;
0136 
0137   // Loop over all IDs
0138   int all(0), allSi(0), good(0);
0139   for (std::map<HGCSiliconDetId, int>::const_iterator itr = waferID_.begin(); itr != waferID_.end(); ++itr) {
0140     HGCSiliconDetId id = itr->first;
0141     ++all;
0142     if ((id.det() == DetId::HGCalEE) || (id.det() == DetId::HGCalHSi)) {
0143       ++allSi;
0144       if (((id.det() == DetId::HGCalEE) && (nameSense_ == "HGCalEESensitive")) ||
0145           ((id.det() == DetId::HGCalHSi) && (nameSense_ == "HGCalHESiliconSensitive"))) {
0146         int partial = ((itr->second) % 100);
0147         int orient = (((itr->second) / 100) % 100);
0148         int type = id.type();
0149         int nCells = hgc.getUVMax(type);
0150         for (int u = 0; u < 2 * nCells; ++u) {
0151           for (int v = 0; v < 2 * nCells; ++v) {
0152             if (((v - u) < nCells) && ((u - v) <= nCells)) {
0153               HGCSiliconDetId hid(id.det(), id.zside(), id.type(), id.layer(), id.waferU(), id.waferV(), u, v);
0154               bool valid = (geom->topology()).valid(static_cast<DetId>(hid));
0155               if (valid) {
0156                 ++good;
0157                 int placeIndex = wafer.cellPlacementIndex(1, HGCalTypes::waferFrontBack(0), orient);
0158                 std::pair<double, double> xy = wafer.cellUV2XY1(u, v, placeIndex, type);
0159                 std::vector<std::pair<double, double> > wxy1 =
0160                     HGCalWaferMask::waferXY(partial, orient, -1, waferSize, 0.0, 0.0, 0.0);
0161                 bool check1 = HGCGuardRing::insidePolygon(xy.first, xy.second, wxy1);
0162                 std::ostringstream st1;
0163                 for (unsigned int k1 = 0; k1 < wxy1.size(); ++k1)
0164                   st1 << " (" << wxy1[k1].first << ", " << wxy1[k1].second << ")";
0165                 edm::LogVerbatim("HGCSim")
0166                     << "First " << hid << " Type:Partial:Orient:Place " << type << ":" << partial << ":" << orient
0167                     << ":" << placeIndex << " Boundary with " << wxy1.size() << " points: " << st1.str() << " check "
0168                     << check1 << " for (" << xy.first << ", " << xy.second << ")";
0169 
0170                 std::vector<std::pair<double, double> > wxy2 =
0171                     HGCalWaferMask::waferXY(partial, orient, -1, waferSize, guardRingOffset_, 0.0, 0.0);
0172                 bool check2 = HGCGuardRing::insidePolygon(xy.first, xy.second, wxy2);
0173                 std::ostringstream st2;
0174                 for (unsigned int k1 = 0; k1 < wxy2.size(); ++k1)
0175                   st2 << " (" << wxy2[k1].first << ", " << wxy2[k1].second << ")";
0176                 edm::LogVerbatim("HGCSim") << "Second Offset " << guardRingOffset_ << " Boundary with " << wxy2.size()
0177                                            << " points: " << st2.str() << " check " << check2 << " for (" << xy.first
0178                                            << ", " << xy.second << ")";
0179               }
0180             }
0181           }
0182         }
0183       }
0184     }
0185   }
0186   edm::LogVerbatim("HGCalSim") << "Total hits = " << all << " Good Silicon DetIds = " << allSi << ":" << good;
0187 }
0188 
0189 //define this as a plug-in
0190 DEFINE_FWK_MODULE(HGCalTestGuardRing);