Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-11-02 03:11:14

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalWaferInFileTest
0004 // Class:      HGCalWaferInFileTest
0005 //
0006 /**\class HGCalWaferInFileTest HGCalWaferInFileTest.cc
0007  test/HGCalWaferInFileTest.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/06/24
0017 //
0018 //
0019 
0020 // system include files
0021 #include <sstream>
0022 #include <string>
0023 #include <vector>
0024 
0025 // user include files
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/EventSetup.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0034 #include "Geometry/HGCalGeometry/interface/HGCalGeometry.h"
0035 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0036 #include "Geometry/HGCalCommonData/interface/HGCalTypes.h"
0037 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0038 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0039 
0040 class HGCalWaferInFileTest : public edm::one::EDAnalyzer<> {
0041 public:
0042   explicit HGCalWaferInFileTest(const edm::ParameterSet&);
0043 
0044   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0045 
0046   void beginJob() override {}
0047   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0048   void endJob() override {}
0049 
0050 private:
0051   std::vector<std::string> getPoints(
0052       double xpos, double ypos, double delX, double delY, double rin, double rout, int lay, int waferU, int waferV);
0053   std::vector<double> getCorners(double xpos, double ypos, double delX, double delY);
0054   std::pair<bool, std::string> getPoints(double xpos,
0055                                          double ypos,
0056                                          double delX,
0057                                          double delY,
0058                                          double rin,
0059                                          double rout,
0060                                          int part,
0061                                          int rotn,
0062                                          int layer,
0063                                          int waferU,
0064                                          int waferV);
0065 
0066   const std::string nameSense_, nameDetector_;
0067   const int verbosity_;
0068   const edm::ESGetToken<HGCalGeometry, IdealGeometryRecord> geomToken_;
0069   double c22_, c27_, c61_, c77_, c88_;
0070 };
0071 
0072 HGCalWaferInFileTest::HGCalWaferInFileTest(const edm::ParameterSet& iC)
0073     : nameSense_(iC.getParameter<std::string>("NameSense")),
0074       nameDetector_(iC.getParameter<std::string>("NameDevice")),
0075       verbosity_(iC.getParameter<int>("Verbosity")),
0076       geomToken_(esConsumes<HGCalGeometry, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_})) {
0077   c22_ = HGCalTypes::c22;
0078   c27_ = HGCalTypes::c27;
0079   c61_ = HGCalTypes::c61;
0080   c77_ = HGCalTypes::c77;
0081   c88_ = HGCalTypes::c88;
0082   edm::LogVerbatim("HGCalGeom") << "Test wafer characteristics for " << nameDetector_ << " using constants of "
0083                                 << nameSense_ << " with verbosity " << verbosity_;
0084 }
0085 
0086 void HGCalWaferInFileTest::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0087   edm::ParameterSetDescription desc;
0088   desc.add<std::string>("NameSense", "HGCalEESensitive");
0089   desc.add<std::string>("NameDevice", "HGCal EE");
0090   desc.add<int>("Verbosity", 0);
0091   descriptions.add("hgcalEEWaferInFileTest", desc);
0092 }
0093 
0094 // ------------ method called to produce the data  ------------
0095 void HGCalWaferInFileTest::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0096   const auto& geomR = iSetup.getData(geomToken_);
0097   const HGCalGeometry* geom = &geomR;
0098   const auto& hgdc = geom->topology().dddConstants();
0099   double delX = 0.5 * hgdc.getParameter()->waferSize_;
0100   double delY = 2.0 * delX / std::sqrt(3.0);
0101   if (hgdc.v17OrLess()) {
0102     c22_ = HGCalTypes::c22O;
0103     c27_ = HGCalTypes::c27O;
0104     c61_ = HGCalTypes::c61O;
0105     c77_ = HGCalTypes::c77O;
0106     c88_ = HGCalTypes::c88O;
0107   }
0108 
0109   edm::LogVerbatim("HGCalGeom") << nameDetector_ << "\nCheck Wafers in file are all valid for " << nameDetector_
0110                                 << std::endl;
0111   if (hgdc.waferHexagon8()) {
0112     DetId::Detector det = (nameSense_ == "HGCalHESiliconSensitive") ? DetId::HGCalHSi : DetId::HGCalEE;
0113     static std::vector<std::string> types = {"F", "b", "g", "gm", "a", "d", "dm", "c", "am", "bm", "X"};
0114     int layers = hgdc.layers(true);
0115     int layerf = hgdc.firstLayer();
0116     std::vector<int> miss(layers, 0);
0117     // See if all entries in the file are valid
0118     int bad1(0);
0119     for (unsigned int k = 0; k < hgdc.waferFileSize(); ++k) {
0120       int indx = hgdc.waferFileIndex(k);
0121       int layer = HGCalWaferIndex::waferLayer(indx);
0122       int waferU = HGCalWaferIndex::waferU(indx);
0123       int waferV = HGCalWaferIndex::waferV(indx);
0124       int type = std::get<0>(hgdc.waferFileInfo(k));
0125       HGCSiliconDetId id(det, 1, type, layer, waferU, waferV, 0, 0);
0126       if (!geom->topology().validModule(id, 3)) {
0127         int part = std::get<1>(hgdc.waferFileInfoFromIndex(indx));
0128         std::string typex = (part < static_cast<int>(types.size())) ? types[part] : "X";
0129         const auto& xy = hgdc.waferPosition(layer, waferU, waferV, true, true);
0130         const auto& rr = hgdc.rangeRLayer(layer, true);
0131         auto points = getPoints(xy.first, xy.second, delX, delY, rr.first, rr.second, layer, waferU, waferV);
0132         auto rpos = getCorners(xy.first, xy.second, delX, delY);
0133         std::ostringstream st1;
0134         st1 << "ID[" << k << "]: (" << (hgdc.getLayerOffset() + layer) << ", " << waferU << ", " << waferV << ", "
0135             << typex << ") at (" << std::setprecision(4) << xy.first << ", " << xy.second << ", "
0136             << hgdc.waferZ(layer, true) << ") not present with " << points.size() << " points:";
0137         for (auto point : points)
0138           st1 << " " << point;
0139         st1 << " in the region " << rr.first << ":" << rr.second << " Corners";
0140         for (auto point : rpos)
0141           st1 << " " << point;
0142         edm::LogVerbatim("HGCalGeom") << st1.str();
0143         ++bad1;
0144         if ((layer - layerf) < layers)
0145           ++miss[layer - layerf];
0146       }
0147     }
0148     edm::LogVerbatim("HGCalGeom") << "\n\nFinds " << bad1 << " invalid wafers among " << hgdc.waferFileSize()
0149                                   << " wafers in the list";
0150     for (unsigned int k = 0; k < miss.size(); ++k) {
0151       if (miss[k] > 0)
0152         edm::LogVerbatim("HGCalGeom") << "Layer[" << k << ":" << (layerf + k) << "] " << miss[k];
0153     }
0154     edm::LogVerbatim("HGCalGeom") << std::endl;
0155 
0156     // Now cross check the content (first type only)
0157     int allG(0), badT(0), badT1(0), badT2(0);
0158     for (unsigned int k = 0; k < hgdc.waferFileSize(); ++k) {
0159       int indx = hgdc.waferFileIndex(k);
0160       int type1 = std::get<0>(hgdc.waferFileInfo(k));
0161       int layer = HGCalWaferIndex::waferLayer(indx);
0162       int waferU = HGCalWaferIndex::waferU(indx);
0163       int waferV = HGCalWaferIndex::waferV(indx);
0164       int type2 = hgdc.waferType(layer, waferU, waferV, false);
0165       HGCSiliconDetId id(det, 1, type2, layer, waferU, waferV, 0, 0);
0166       if (geom->topology().validModule(id, 3)) {
0167         ++allG;
0168         bool typeOK = (type1 == type2);
0169         if (!typeOK) {
0170           ++badT;
0171           if (type1 == 0)
0172             ++badT1;
0173           else if (type2 == 0)
0174             ++badT2;
0175           const auto& xy = hgdc.waferPosition(layer, waferU, waferV, true, false);
0176           edm::LogVerbatim("HGCalGeom") << "ID[" << k << "]: (" << (hgdc.getLayerOffset() + layer) << ", " << waferU
0177                                         << ", " << waferV << ", " << type1 << ":" << type2 << ") at ("
0178                                         << std::setprecision(4) << xy.first << ", " << xy.second << ", "
0179                                         << hgdc.waferZ(layer, true) << ") failure flag " << typeOK;
0180         }
0181       }
0182     }
0183     edm::LogVerbatim("HGCalGeom") << "\n\nFinds " << badT << "[" << badT1 << ":" << badT2 << "] mismatch in type among "
0184                                   << allG << " wafers with the same indices\n";
0185 
0186     // Now cross check the content (partial and orientation)
0187     int allX(0), badG(0), badP(0), badP2(0), badR(0), badR2(0);
0188     std::vector<int> wrongP(layers, 0), wrongP2(layers, 0), wrongR(layers, 0), wrongR2(layers, 0);
0189     for (unsigned int k = 0; k < hgdc.waferFileSize(); ++k) {
0190       int indx = hgdc.waferFileIndex(k);
0191       int part1 = std::get<1>(hgdc.waferFileInfo(k));
0192       int rotn1 = std::get<2>(hgdc.waferFileInfo(k));
0193       int layer = HGCalWaferIndex::waferLayer(indx);
0194       int waferU = HGCalWaferIndex::waferU(indx);
0195       int waferV = HGCalWaferIndex::waferV(indx);
0196       int type2 = hgdc.waferType(layer, waferU, waferV, false);
0197       HGCSiliconDetId id(det, 1, type2, layer, waferU, waferV, 0, 0);
0198       if (geom->topology().validModule(id, 3)) {
0199         ++allX;
0200         int part2 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), false, false).first;
0201         int rotn2 = hgdc.waferTypeRotation(id.layer(), id.waferU(), id.waferV(), false, false).second;
0202         bool partOK = ((part1 == part2) || ((part1 == HGCalTypes::WaferFull) && (part2 == HGCalTypes::WaferOut)));
0203         bool rotnOK = ((rotn1 == rotn2) || (part1 == HGCalTypes::WaferFull) || (part2 == HGCalTypes::WaferFull));
0204         bool partOK2 = (partOK) || (part2 < part1);
0205         bool rotnOK2(true);
0206         std::string miss;
0207         if (!partOK) {
0208           ++badP;
0209           if ((layer - layerf) < layers)
0210             ++wrongP[layer - layerf];
0211         }
0212         const auto& xy = hgdc.waferPosition(layer, waferU, waferV, true, false);
0213         const auto& rr = hgdc.rangeRLayer(layer, true);
0214 
0215         if (!partOK2) {
0216           ++badP2;
0217           if ((layer - layerf) < layers)
0218             ++wrongP2[layer - layerf];
0219           auto result =
0220               getPoints(xy.first, xy.second, delX, delY, rr.first, rr.second, part1, rotn1, layer, waferU, waferV);
0221           rotnOK2 = result.first;
0222           miss = result.second;
0223           if (!rotnOK2) {
0224             ++badR2;
0225             if ((layer - layerf) < layers)
0226               ++wrongR2[layer - layerf];
0227           }
0228         }
0229         if (!rotnOK) {
0230           ++badR;
0231           if ((layer - layerf) < layers)
0232             ++wrongR[layer - layerf];
0233         }
0234         if ((!partOK) || (!rotnOK)) {
0235           ++badG;
0236           if ((verbosity_ > 0) || ((!partOK2) || (!rotnOK2))) {
0237             std::string partx1 = (part1 < static_cast<int>(types.size())) ? types[part1] : "X";
0238             std::string partx2 = (part2 < static_cast<int>(types.size())) ? types[part2] : "X";
0239             auto points = getPoints(xy.first, xy.second, delX, delY, rr.first, rr.second, layer, waferU, waferV);
0240             auto rpos = getCorners(xy.first, xy.second, delX, delY);
0241             std::ostringstream st1;
0242             st1 << "ID[" << k << "]: (" << (hgdc.getLayerOffset() + layer) << ", " << waferU << ", " << waferV << ","
0243                 << type2 << ", " << partx1 << ":" << partx2 << ":" << part1 << ":" << part2 << ", " << rotn1 << ":"
0244                 << rotn2 << ") at (" << std::setprecision(4) << xy.first << ", " << xy.second << ", "
0245                 << hgdc.waferZ(layer, true) << ") failure flags (part = " << partOK << ":" << partOK2 << " rotn "
0246                 << rotnOK << ":" << rotnOK2 << " at " << miss << " with " << points.size() << " points:";
0247             for (auto point : points)
0248               st1 << " " << point;
0249             st1 << " in the region " << rr.first << ":" << rr.second << " Corners";
0250             for (auto point : rpos)
0251               st1 << " " << point;
0252             edm::LogVerbatim("HGCalGeom") << st1.str();
0253           }
0254         }
0255       }
0256     }
0257     edm::LogVerbatim("HGCalGeom") << "\n\nFinds " << badG << " (" << badP << ":" << badP2 << ":" << badR << ":" << badR2
0258                                   << ") mismatch in partial|orientation among " << allX
0259                                   << " wafers with the same indices";
0260     for (int k = 0; k < layers; ++k) {
0261       if ((wrongP[k] > 0) || (wrongR[k] > 0))
0262         edm::LogVerbatim("HGCalGeom") << "Layer[" << k << ":" << (layerf + k) << "] " << wrongP[k] << ":" << wrongP2[k]
0263                                       << ":" << wrongR[k] << ":" << wrongR2[k];
0264     }
0265     edm::LogVerbatim("HGCalGeom") << std::endl;
0266   }
0267 }
0268 
0269 std::vector<std::string> HGCalWaferInFileTest::getPoints(
0270     double xpos, double ypos, double delX, double delY, double rin, double rout, int layer, int waferU, int waferV) {
0271   std::vector<std::string> points;
0272   static const int corners = 6;
0273   static const int base = 10;
0274   static const std::string c0[corners] = {"A0", "B0", "C0", "D0", "E0", "F0"};
0275   double dx0[corners] = {
0276       0.0, HGCalTypes::c10 * delX, HGCalTypes::c10 * delX, 0.0, -HGCalTypes::c10 * delX, -HGCalTypes::c10 * delX};
0277   double dy0[corners] = {-HGCalTypes::c10 * delY,
0278                          -HGCalTypes::c50 * delY,
0279                          HGCalTypes::c50 * delY,
0280                          HGCalTypes::c10 * delY,
0281                          HGCalTypes::c50 * delY,
0282                          -HGCalTypes::c50 * delY};
0283   double xc[corners], yc[corners];
0284   int ncor(0), iok(0);
0285   for (int k = 0; k < corners; ++k) {
0286     xc[k] = xpos + dx0[k];
0287     yc[k] = ypos + dy0[k];
0288     double rpos = sqrt(xc[k] * xc[k] + yc[k] * yc[k]);
0289     if (rpos <= rout && rpos >= rin) {
0290       ++ncor;
0291       iok = iok * base + 1;
0292       points.emplace_back(c0[k]);
0293     } else {
0294       iok *= base;
0295     }
0296   }
0297   if (verbosity_ > 1)
0298     edm::LogVerbatim("HGCalGeom") << "I/p " << layer << ":" << waferU << ":" << waferV << ":" << xpos << ":" << ypos
0299                                   << ":" << delX << ":" << delY << ":" << rin << ":" << rout << " Corners " << ncor
0300                                   << " iok " << iok;
0301 
0302   static const int parts = 3;
0303   static const std::string c1[parts] = {"A1", "A2", "A3"};
0304   double dx1[parts] = {c22_ * delX, HGCalTypes::c50 * delX, c77_ * delX};
0305   double dy1[parts] = {-c88_ * delY, -HGCalTypes::c75 * delY, -c61_ * delY};
0306   if ((((iok / 10000) % 10) == 1) && (((iok / 100000) % 10) == 0)) {
0307     for (int k = 0; k < parts; ++k) {
0308       double xc1 = xpos + dx1[k];
0309       double yc1 = ypos + dy1[k];
0310       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0311       if (rpos <= rout && rpos >= rin) {
0312         points.emplace_back(c1[k]);
0313         break;
0314       }
0315     }
0316   }
0317   if ((((iok / 10000) % 10) == 0) && (((iok / 100000) % 10) == 1)) {
0318     for (int k1 = 0; k1 < parts; ++k1) {
0319       int k = parts - k1 - 1;
0320       double xc1 = xpos + dx1[k];
0321       double yc1 = ypos + dy1[k];
0322       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0323       if (rpos <= rout && rpos >= rin) {
0324         points.emplace_back(c1[k]);
0325         break;
0326       }
0327     }
0328   }
0329 
0330   static const std::string c2[parts] = {"B1", "B2", "B3"};
0331   double dx2[parts] = {HGCalTypes::c10 * delX, HGCalTypes::c10 * delX, HGCalTypes::c10 * delX};
0332   double dy2[parts] = {-c27_ * delY, 0.0, c27_ * delY};
0333   if ((((iok / 1000) % 10) == 1) && (((iok / 10000) % 10) == 0)) {
0334     for (int k = 0; k < parts; ++k) {
0335       double xc1 = xpos + dx2[k];
0336       double yc1 = ypos + dy2[k];
0337       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0338       if (rpos <= rout && rpos >= rin) {
0339         points.emplace_back(c2[k]);
0340         break;
0341       }
0342     }
0343   }
0344   if ((((iok / 1000) % 10) == 0) && (((iok / 10000) % 10) == 1)) {
0345     for (int k1 = 0; k1 < parts; ++k1) {
0346       int k = parts - k1 - 1;
0347       double xc1 = xpos + dx2[k];
0348       double yc1 = ypos + dy2[k];
0349       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0350       if (rpos <= rout && rpos >= rin) {
0351         points.emplace_back(c2[k]);
0352         break;
0353       }
0354     }
0355   }
0356 
0357   static const std::string c3[parts] = {"C1", "C2", "C3"};
0358   double dx3[parts] = {c77_ * delX, HGCalTypes::c50 * delX, c22_ * delX};
0359   double dy3[parts] = {c61_ * delY, HGCalTypes::c75 * delY, c88_ * delY};
0360   if ((((iok / 100) % 10) == 1) && (((iok / 1000) % 10) == 0)) {
0361     for (int k = 0; k < parts; ++k) {
0362       double xc1 = xpos + dx3[k];
0363       double yc1 = ypos + dy3[k];
0364       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0365       if (rpos <= rout && rpos >= rin) {
0366         points.emplace_back(c3[k]);
0367         break;
0368       }
0369     }
0370   }
0371   if ((((iok / 100) % 10) == 0) && (((iok / 1000) % 10) == 1)) {
0372     for (int k1 = 0; k1 < parts; ++k1) {
0373       int k = parts - k1 - 1;
0374       double xc1 = xpos + dx3[k];
0375       double yc1 = ypos + dy3[k];
0376       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0377       if (rpos <= rout && rpos >= rin) {
0378         points.emplace_back(c3[k]);
0379         break;
0380       }
0381     }
0382   }
0383 
0384   static const std::string c4[parts] = {"D1", "D2", "D3"};
0385   double dx4[parts] = {-c22_ * delX, -HGCalTypes::c50 * delX, -c77_ * delX};
0386   double dy4[parts] = {c88_ * delY, HGCalTypes::c75 * delY, c61_ * delY};
0387   if ((((iok / 10) % 10) == 1) && (((iok / 100) % 10) == 0)) {
0388     for (int k = 0; k < parts; ++k) {
0389       double xc1 = xpos + dx4[k];
0390       double yc1 = ypos + dy4[k];
0391       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0392       if (rpos <= rout && rpos >= rin) {
0393         points.emplace_back(c4[k]);
0394         break;
0395       }
0396     }
0397   }
0398   if ((((iok / 10) % 10) == 0) && (((iok / 100) % 10) == 1)) {
0399     for (int k1 = 0; k1 < parts; ++k1) {
0400       int k = parts - k1 - 1;
0401       double xc1 = xpos + dx4[k];
0402       double yc1 = ypos + dy4[k];
0403       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0404       if (rpos <= rout && rpos >= rin) {
0405         points.emplace_back(c4[k]);
0406         break;
0407       }
0408     }
0409   }
0410 
0411   static const std::string c5[parts] = {"E1", "E2", "E3"};
0412   double dx5[parts] = {-delX, -delX, -delX};
0413   double dy5[parts] = {c27_ * delY, 0.0, -c27_ * delY};
0414   if ((((iok / 1) % 10) == 1) && (((iok / 10) % 10) == 0)) {
0415     for (int k = 0; k < parts; ++k) {
0416       double xc1 = xpos + dx5[k];
0417       double yc1 = ypos + dy5[k];
0418       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0419       if (rpos <= rout && rpos >= rin) {
0420         points.emplace_back(c5[k]);
0421         break;
0422       }
0423     }
0424   }
0425   if ((((iok / 1) % 10) == 0) && (((iok / 10) % 10) == 1)) {
0426     for (int k1 = 0; k1 < parts; ++k1) {
0427       int k = parts - k1 - 1;
0428       double xc1 = xpos + dx5[k];
0429       double yc1 = ypos + dy5[k];
0430       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0431       if (rpos <= rout && rpos >= rin) {
0432         points.emplace_back(c5[k]);
0433         break;
0434       }
0435     }
0436   }
0437 
0438   static const std::string c6[parts] = {"F1", "F2", "F3"};
0439   double dx6[parts] = {-c77_ * delX, -HGCalTypes::c50 * delX, -c22_ * delX};
0440   double dy6[parts] = {-c61_ * delY, -HGCalTypes::c75 * delY, -c88_ * delY};
0441   if ((((iok / 100000) % 10) == 1) && (((iok / 1) % 10) == 0)) {
0442     for (int k = 0; k < parts; ++k) {
0443       double xc1 = xpos + dx6[k];
0444       double yc1 = ypos + dy6[k];
0445       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0446       if (rpos <= rout && rpos >= rin) {
0447         points.emplace_back(c6[k]);
0448         break;
0449       }
0450     }
0451   }
0452   if ((((iok / 100000) % 10) == 0) && (((iok / 1) % 10) == 1)) {
0453     for (int k1 = 0; k1 < parts; ++k1) {
0454       int k = parts - k1 - 1;
0455       double xc1 = xpos + dx6[k];
0456       double yc1 = ypos + dy6[k];
0457       double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0458       if (rpos <= rout && rpos >= rin) {
0459         points.emplace_back(c6[k]);
0460         break;
0461       }
0462     }
0463   }
0464 
0465   return points;
0466 }
0467 
0468 std::vector<double> HGCalWaferInFileTest::getCorners(double xpos, double ypos, double delX, double delY) {
0469   std::vector<double> points;
0470   static const int corners = 6;
0471   double dx0[corners] = {
0472       0.0, HGCalTypes::c10 * delX, HGCalTypes::c10 * delX, 0.0, -HGCalTypes::c10 * delX, -HGCalTypes::c10 * delX};
0473   double dy0[corners] = {-HGCalTypes::c10 * delY,
0474                          -HGCalTypes::c50 * delY,
0475                          HGCalTypes::c50 * delY,
0476                          HGCalTypes::c10 * delY,
0477                          HGCalTypes::c50 * delY,
0478                          -HGCalTypes::c50 * delY};
0479   for (int k = 0; k < corners; ++k) {
0480     double xc = xpos + dx0[k];
0481     double yc = ypos + dy0[k];
0482     double rpos = sqrt(xc * xc + yc * yc);
0483     points.emplace_back(rpos);
0484   }
0485   return points;
0486 }
0487 
0488 std::pair<bool, std::string> HGCalWaferInFileTest::getPoints(double xpos,
0489                                                              double ypos,
0490                                                              double delX,
0491                                                              double delY,
0492                                                              double rin,
0493                                                              double rout,
0494                                                              int part,
0495                                                              int rotn,
0496                                                              int layer,
0497                                                              int waferU,
0498                                                              int waferV) {
0499   std::string point;
0500   static const int corners = 6;
0501   static const int corner2 = 12;
0502   static const int base = 10;
0503   static const int base2 = 100;
0504   static const std::string c0[corners] = {"A0", "B0", "C0", "D0", "E0", "F0"};
0505   static const std::string c1[corners] = {"A2", "B2", "C2", "D2", "E2", "F2"};
0506   static const std::string c2[corner2] = {"A1", "A3", "B1", "B3", "C1", "C3", "D1", "D3", "E1", "E3", "F1", "F3"};
0507   double dx0[corners] = {
0508       0.0, HGCalTypes::c10 * delX, HGCalTypes::c10 * delX, 0.0, -HGCalTypes::c10 * delX, -HGCalTypes::c10 * delX};
0509   double dy0[corners] = {-HGCalTypes::c10 * delY,
0510                          -HGCalTypes::c50 * delY,
0511                          HGCalTypes::c50 * delY,
0512                          HGCalTypes::c10 * delY,
0513                          HGCalTypes::c50 * delY,
0514                          -HGCalTypes::c50 * delY};
0515   double dx1[corners] = {HGCalTypes::c50 * delX,
0516                          HGCalTypes::c10 * delX,
0517                          HGCalTypes::c50 * delX,
0518                          -HGCalTypes::c50 * delX,
0519                          -HGCalTypes::c10 * delX,
0520                          -HGCalTypes::c50 * delX};
0521   double dy1[corners] = {
0522       -HGCalTypes::c75 * delY, 0.0, HGCalTypes::c75 * delY, HGCalTypes::c75 * delY, 0.0, -HGCalTypes::c75 * delY};
0523   double dx2[corner2] = {c22_ * delX,
0524                          c77_ * delX,
0525                          HGCalTypes::c10 * delX,
0526                          HGCalTypes::c10 * delX,
0527                          c77_ * delX,
0528                          c22_ * delX,
0529                          -c22_ * delX,
0530                          -c77_ * delX,
0531                          -HGCalTypes::c10 * delX,
0532                          -HGCalTypes::c10 * delX,
0533                          -c77_ * delX,
0534                          -c22_ * delX};
0535   double dy2[corner2] = {-c88_ * delY,
0536                          -c61_ * delY,
0537                          -c27_ * delY,
0538                          c27_ * delY,
0539                          c61_ * delY,
0540                          c88_ * delY,
0541                          c88_ * delY,
0542                          c61_ * delY,
0543                          c27_ * delY,
0544                          -c27_ * delY,
0545                          -c61_ * delY,
0546                          -c88_ * delY};
0547   bool ok(true);
0548   switch (part) {
0549     case (HGCalTypes::WaferThree): {
0550       static const int nc0[corners] = {450, 150, 201, 312, 423, 534};
0551       int nc = nc0[rotn];
0552       for (int k1 = 0; k1 < 3; ++k1) {
0553         int k = nc % base;
0554         double xc1 = xpos + dx0[k];
0555         double yc1 = ypos + dy0[k];
0556         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0557         if (rpos <= rout && rpos >= rin) {
0558           point = c0[k];
0559           ok = false;
0560           break;
0561         }
0562         nc /= base;
0563       }
0564       break;
0565     }
0566     case (HGCalTypes::WaferSemi2): {
0567       static const int nc10[corners] = {450, 150, 201, 312, 423, 534};
0568       static const int nc11[corners] = {700, 902, 1104, 106, 308, 510};
0569       int nc = nc10[rotn];
0570       for (int k1 = 0; k1 < 3; ++k1) {
0571         int k = nc % base;
0572         double xc1 = xpos + dx0[k];
0573         double yc1 = ypos + dy0[k];
0574         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0575         if (rpos <= rout && rpos >= rin) {
0576           point = c0[k];
0577           ok = false;
0578           break;
0579         }
0580         nc /= base;
0581       }
0582       nc = nc11[rotn];
0583       for (int k1 = 0; k1 < 2; ++k1) {
0584         int k = nc % base2;
0585         double xc1 = xpos + dx2[k];
0586         double yc1 = ypos + dy2[k];
0587         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0588         if (rpos <= rout && rpos >= rin) {
0589           point = c2[k];
0590           ok = false;
0591           break;
0592         }
0593         nc /= base2;
0594       }
0595       break;
0596     }
0597     case (HGCalTypes::WaferSemi): {
0598       static const int nc20[corners] = {450, 150, 201, 312, 423, 534};
0599       static const int nc21[corners] = {30, 14, 25, 30, 41, 52};
0600       int nc = nc20[rotn];
0601       for (int k1 = 0; k1 < 3; ++k1) {
0602         int k = nc % base;
0603         double xc1 = xpos + dx0[k];
0604         double yc1 = ypos + dy0[k];
0605         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0606         if (rpos <= rout && rpos >= rin) {
0607           point = c0[k];
0608           ok = false;
0609           break;
0610         }
0611         nc /= base;
0612       }
0613       nc = nc21[rotn];
0614       for (int k1 = 0; k1 < 2; ++k1) {
0615         int k = nc % base;
0616         double xc1 = xpos + dx1[k];
0617         double yc1 = ypos + dy1[k];
0618         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0619         if (rpos <= rout && rpos >= rin) {
0620           point = c1[k];
0621           ok = false;
0622           break;
0623         }
0624         nc /= base;
0625       }
0626       break;
0627     }
0628     case (HGCalTypes::WaferHalf): {
0629       static const int nc3[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
0630       int nc = nc3[rotn];
0631       for (int k1 = 0; k1 < 4; ++k1) {
0632         int k = nc % base;
0633         double xc1 = xpos + dx0[k];
0634         double yc1 = ypos + dy0[k];
0635         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0636         if (rpos <= rout && rpos >= rin) {
0637           point = c0[k];
0638           ok = false;
0639           break;
0640         }
0641         nc /= base;
0642       }
0643       break;
0644     }
0645     case (HGCalTypes::WaferChopTwoM): {
0646       static const int nc40[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
0647       static const int nc41[corners] = {500, 702, 904, 1106, 108, 310};
0648       int nc = nc40[rotn];
0649       for (int k1 = 0; k1 < 4; ++k1) {
0650         int k = nc % base;
0651         double xc1 = xpos + dx0[k];
0652         double yc1 = ypos + dy0[k];
0653         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0654         if (rpos <= rout && rpos >= rin) {
0655           point = c0[k];
0656           ok = false;
0657           break;
0658         }
0659         nc /= base;
0660       }
0661       nc = nc41[rotn];
0662       for (int k1 = 0; k1 < 2; ++k1) {
0663         int k = nc % base2;
0664         double xc1 = xpos + dx2[k];
0665         double yc1 = ypos + dy2[k];
0666         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0667         if (rpos <= rout && rpos >= rin) {
0668           point = c2[k];
0669           ok = false;
0670           break;
0671         }
0672         nc /= base2;
0673       }
0674       break;
0675     }
0676     case (HGCalTypes::WaferChopTwo): {
0677       static const int nc50[corners] = {3450, 1450, 2501, 3012, 4123, 5234};
0678       static const int nc51[corners] = {20, 13, 24, 35, 40, 51};
0679       int nc = nc50[rotn];
0680       for (int k1 = 0; k1 < 4; ++k1) {
0681         int k = nc % base;
0682         double xc1 = xpos + dx0[k];
0683         double yc1 = ypos + dy0[k];
0684         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0685         if (rpos <= rout && rpos >= rin) {
0686           point = c0[k];
0687           ok = false;
0688           break;
0689         }
0690         nc /= base;
0691       }
0692       nc = nc51[rotn];
0693       for (int k1 = 0; k1 < 2; ++k1) {
0694         int k = nc % base;
0695         double xc1 = xpos + dx1[k];
0696         double yc1 = ypos + dy1[k];
0697         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0698         if (rpos <= rout && rpos >= rin) {
0699           point = c1[k];
0700           ok = false;
0701           break;
0702         }
0703         nc /= base;
0704       }
0705       break;
0706     }
0707     case (HGCalTypes::WaferFive): {
0708       static const int nc6[corners] = {23450, 13450, 24501, 35012, 40123, 51234};
0709       int nc = nc6[rotn];
0710       for (int k1 = 0; k1 < 5; ++k1) {
0711         int k = nc % base;
0712         double xc1 = xpos + dx0[k];
0713         double yc1 = ypos + dy0[k];
0714         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0715         if (rpos <= rout && rpos >= rin) {
0716           point = c0[k];
0717           ok = false;
0718           break;
0719         }
0720       }
0721       break;
0722     }
0723     default: {
0724       for (int k = 0; k < corners; ++k) {
0725         double xc1 = xpos + dx0[k];
0726         double yc1 = ypos + dy0[k];
0727         double rpos = sqrt(xc1 * xc1 + yc1 * yc1);
0728         if (rpos <= rout && rpos >= rin) {
0729           point = c0[k];
0730           ok = false;
0731           break;
0732         }
0733       }
0734       break;
0735     }
0736   }
0737   if (verbosity_ > 1)
0738     edm::LogVerbatim("HGCalGeom") << "I/p " << layer << ":" << waferU << ":" << waferV << ":" << xpos << ":" << ypos
0739                                   << ":" << delX << ":" << delY << ":" << rin << ":" << rout << " Results " << ok
0740                                   << " point " << point;
0741   return std::make_pair(ok, point);
0742 }
0743 
0744 // define this as a plug-in
0745 DEFINE_FWK_MODULE(HGCalWaferInFileTest);