Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-16 22:52:23

0001 // system include files
0002 #include <memory>
0003 #include <string>
0004 #include <vector>
0005 #include <cmath>
0006 
0007 // user include files
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0018 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0019 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0020 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0021 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0022 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0023 
0024 class HGCalPartialWaferTester : public edm::one::EDAnalyzer<> {
0025 public:
0026   explicit HGCalPartialWaferTester(const edm::ParameterSet&);
0027   ~HGCalPartialWaferTester() override = default;
0028   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029 
0030   void beginJob() override {}
0031   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0032   void endJob() override {}
0033 
0034 private:
0035   const std::string nameSense_;
0036   const std::vector<int> orientations_, partialTypes_;
0037   const int nTrials_;
0038   const edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord> dddToken_;
0039 };
0040 
0041 HGCalPartialWaferTester::HGCalPartialWaferTester(const edm::ParameterSet& iC)
0042     : nameSense_(iC.getParameter<std::string>("nameSense")),
0043       orientations_(iC.getParameter<std::vector<int>>("waferOrientations")),
0044       partialTypes_(iC.getParameter<std::vector<int>>("partialTypes")),
0045       nTrials_(iC.getParameter<int>("numberOfTrials")),
0046       dddToken_(esConsumes<HGCalDDDConstants, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_})) {
0047   edm::LogVerbatim("HGCalGeom") << "Test positions for " << partialTypes_.size() << " partial wafer types "
0048                                 << " and " << orientations_.size() << " Orientations for " << nameSense_;
0049 }
0050 
0051 void HGCalPartialWaferTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0052   edm::ParameterSetDescription desc;
0053   std::vector<int> orients = {HGCalTypes::WaferOrient0,
0054                               HGCalTypes::WaferOrient1,
0055                               HGCalTypes::WaferOrient2,
0056                               HGCalTypes::WaferOrient3,
0057                               HGCalTypes::WaferOrient4,
0058                               HGCalTypes::WaferOrient5};
0059   std::vector<int> types = {HGCalTypes::WaferLDTop,
0060                             HGCalTypes::WaferLDBottom,
0061                             HGCalTypes::WaferLDLeft,
0062                             HGCalTypes::WaferLDRight,
0063                             HGCalTypes::WaferLDFive,
0064                             HGCalTypes::WaferLDThree,
0065                             HGCalTypes::WaferHDTop,
0066                             HGCalTypes::WaferHDBottom,
0067                             HGCalTypes::WaferHDLeft,
0068                             HGCalTypes::WaferHDRight,
0069                             HGCalTypes::WaferHDFive};
0070   desc.add<std::string>("nameSense", "HGCalHESiliconSensitive");
0071   desc.add<std::vector<int>>("waferOrientations", orients);
0072   desc.add<std::vector<int>>("partialTypes", types);
0073   desc.add<int>("numberOfTrials", 1000);
0074   descriptions.add("hgcalPartialWaferTester", desc);
0075 }
0076 
0077 void HGCalPartialWaferTester::analyze(const edm::Event&, const edm::EventSetup& iSetup) {
0078   const HGCalDDDConstants& hgdc = iSetup.getData(dddToken_);
0079   const bool reco(true), all(true), norot(false), debug1(false), debug2(false), extend(true);
0080   for (const auto& partialType : partialTypes_) {
0081     for (const auto& orientation : orientations_) {
0082       int indx(-1), type(-1);
0083       for (auto itr = hgdc.getParameter()->waferInfoMap_.begin(); itr != hgdc.getParameter()->waferInfoMap_.end();
0084            ++itr) {
0085         if (((itr->second).part == partialType) && ((itr->second).orient == orientation)) {
0086           indx = itr->first;
0087           type = (itr->second).type;
0088           break;
0089         }
0090       }
0091 
0092       if (indx > 0) {
0093         int alltry(0), error(0);
0094         int layer = HGCalWaferIndex::waferLayer(indx);
0095         int waferU = HGCalWaferIndex::waferU(indx);
0096         int waferV = HGCalWaferIndex::waferV(indx);
0097         auto xy = hgdc.waferPosition(layer, waferU, waferV, true, false);
0098         edm::LogVerbatim("HGCalGeom") << "\n\nPartial Type " << partialType << " Orientation " << orientation
0099                                       << " Wafer " << waferU << ":" << waferV << " in layer " << layer << " at "
0100                                       << xy.first << ":" << xy.second << "\n\n";
0101         int nCells = (type == 0) ? HGCSiliconDetId::HGCalHighDensityN : HGCSiliconDetId::HGCalLowDensityN;
0102         for (int i = 0; i < nTrials_; i++) {
0103           int ui = std::floor(nCells * 0.0002 * (rand() % 10000));
0104           int vi = std::floor(nCells * 0.0002 * (rand() % 10000));
0105           if ((ui < 2 * nCells) && (vi < 2 * nCells) && ((vi - ui) < nCells) && ((ui - vi) <= nCells) &&
0106               HGCalWaferMask::goodCell(ui, vi, partialType)) {
0107             ++alltry;
0108             double zpos = hgdc.waferZ(layer, reco);
0109             int zside = (zpos > 0) ? 1 : -1;
0110             auto xy = hgdc.locateCell(zside, layer, waferU, waferV, ui, vi, reco, all, norot, false, debug1);
0111             int lay(layer), cU(0), cV(0), wType(-1), wU(0), wV(0);
0112             double wt(0);
0113             hgdc.waferFromPosition(HGCalParameters::k_ScaleToDDD * xy.first,
0114                                    HGCalParameters::k_ScaleToDDD * xy.second,
0115                                    zside,
0116                                    lay,
0117                                    wU,
0118                                    wV,
0119                                    cU,
0120                                    cV,
0121                                    wType,
0122                                    wt,
0123                                    extend,
0124                                    debug2);
0125             bool ok =
0126                 ((wType == type) && (layer == lay) && (waferU == wU) && (waferV == wV) && (ui == cU) && (vi == cV));
0127             std::string comment = (ok) ? "" : " ***** ERROR *****";
0128             edm::LogVerbatim("HGCalGeom")
0129                 << "Layer " << layer << ":" << lay << " waferU " << waferU << ":" << wU << " waferV " << waferV << ":"
0130                 << wV << " Type " << type << ":" << wType << " cellU " << ui << ":" << cU << " cellV " << vi << ":"
0131                 << cV << " position " << xy.first << ":" << xy.second << comment;
0132             if (!ok)
0133               ++error;
0134           }
0135         }
0136         edm::LogVerbatim("HGCalGeom") << "\n\nFound " << error << " errors among " << alltry << ":" << nTrials_
0137                                       << " trials for partial type " << partialType << " orientation " << orientation
0138                                       << " of " << nameSense_;
0139       } else {
0140         edm::LogVerbatim("HGCalGeom") << "\n\nCannot find a wafer of type " << partialType << " and orientation "
0141                                       << orientation << " for " << nameSense_;
0142       }
0143     }
0144   }
0145 }
0146 
0147 // define this as a plug-in
0148 DEFINE_FWK_MODULE(HGCalPartialWaferTester);