Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalCellPartialCellTester
0004 // Class:      HGCalCellPartialCellTester
0005 //
0006 /**\class HGCalCellPartialCellTester HGCalCellPartialCellTester.cc
0007  test/HGCalCellPartialCellTester.cc
0008 
0009  Description: <one line class summary>
0010 
0011  Implementation:
0012      <Notes on implementation>
0013 */
0014 //
0015 // Original Author:  Sunanda Banerjee, Pruthvi Suryadevara
0016 //         Created:  Mon 2022/06/10
0017 //
0018 //
0019 
0020 // system include files
0021 #include <fstream>
0022 #include <iostream>
0023 #include <memory>
0024 #include <string>
0025 #include <vector>
0026 #include <stdlib.h>
0027 #include <cmath>
0028 //#include <chrono>
0029 
0030 // user include files
0031 #include "FWCore/Framework/interface/Frameworkfwd.h"
0032 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0033 
0034 #include "FWCore/Framework/interface/Event.h"
0035 #include "FWCore/Framework/interface/EventSetup.h"
0036 #include "FWCore/Framework/interface/MakerMacros.h"
0037 
0038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0040 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0041 #include "Geometry/HGCalCommonData/interface/HGCalDDDConstants.h"
0042 #include "Geometry/HGCalCommonData/interface/HGCalParameters.h"
0043 #include "Geometry/HGCalCommonData/interface/HGCalCellUV.h"
0044 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0045 #include "Geometry/HGCalCommonData/interface/HGCalWaferMask.h"
0046 
0047 class HGCalPartialCellTester : public edm::one::EDAnalyzer<> {
0048 public:
0049   explicit HGCalPartialCellTester(const edm::ParameterSet&);
0050   ~HGCalPartialCellTester() override = default;
0051   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0052 
0053   void beginJob() override {}
0054   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0055   void endJob() override {}
0056 
0057 private:
0058   const double waferSize_;
0059   const int waferType_;
0060   const int placeIndex_;
0061   const int partialType_;
0062   const int nTrials_;
0063   const bool v17OrLess_;
0064   const int modeUV_;
0065 };
0066 
0067 HGCalPartialCellTester::HGCalPartialCellTester(const edm::ParameterSet& iC)
0068     : waferSize_(iC.getParameter<double>("waferSize")),
0069       waferType_(iC.getParameter<int>("waferType")),
0070       placeIndex_(iC.getParameter<int>("cellPlacementIndex")),
0071       partialType_(iC.getParameter<int>("partialType")),
0072       nTrials_(iC.getParameter<int>("numbberOfTrials")),
0073       v17OrLess_(iC.getParameter<bool>("v17OrLess")),
0074       modeUV_(iC.getParameter<int>("modeUV")) {
0075   edm::LogVerbatim("HGCalGeom") << "Test positions for wafer of size " << waferSize_ << " Type " << waferType_
0076                                 << " Partial Type " << partialType_ << " Placement Index " << placeIndex_ << " mode "
0077                                 << modeUV_ << " with " << nTrials_ << " trials"
0078                                 << " v17OrLess " << v17OrLess_;
0079 }
0080 
0081 void HGCalPartialCellTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0082   edm::ParameterSetDescription desc;
0083   desc.add<double>("waferSize", 166.4408);
0084   desc.add<int>("waferType", 0);
0085   desc.add<int>("cellPlacementIndex", 3);
0086   desc.add<int>("partialType", 25);
0087   desc.add<int>("numbberOfTrials", 1000);
0088   desc.add<bool>("v17OrLess", false);
0089   desc.add<int>("modeUV", 0);
0090   descriptions.add("hgcalPartialCellTester", desc);
0091 }
0092 
0093 // ------------ method called to produce the data  ------------
0094 void HGCalPartialCellTester::analyze(const edm::Event&, const edm::EventSetup&) {
0095   const int nFine(12), nCoarse(8);
0096   double r2 = 0.5 * waferSize_;
0097   double R2 = 2 * r2 / sqrt(3);
0098   int nCells = (waferType_ == 0) ? nFine : nCoarse;
0099   HGCalCellUV wafer(waferSize_, 0.0, nFine, nCoarse);
0100   HGCalCell wafer2(waferSize_, nFine, nCoarse);
0101   edm::LogVerbatim("HGCalGeom") << "\nHGCalPartialCellTester:: nCells " << nCells << " and placement index "
0102                                 << placeIndex_ << "\n\n";
0103   auto start_t = std::chrono::high_resolution_clock::now();
0104 
0105   if (modeUV_ <= 0) {
0106     for (int i = 0; i < nTrials_; i++) {
0107       double xi = (2 * r2 * static_cast<double>(rand()) / RAND_MAX) - r2;
0108       double yi = (2 * R2 * static_cast<double>(rand()) / RAND_MAX) - R2;
0109       bool goodPoint = true;
0110       int ug = 0;
0111       int vg = 0;
0112       if (partialType_ == 11 || partialType_ == 13 || partialType_ == 15 || partialType_ == 21 || partialType_ == 23 ||
0113           partialType_ == 25) {
0114         ug = 0;
0115         vg = 0;
0116       } else if (partialType_ == 12 || partialType_ == 14 || partialType_ == 16 || partialType_ == 22 ||
0117                  partialType_ == 24) {
0118         ug = nCells + 1;
0119         vg = 2 * (nCells - 1);
0120       }
0121       std::pair<double, double> xyg = wafer2.cellUV2XY2(ug, vg, placeIndex_, waferType_);
0122       std::vector<std::pair<double, double> > wxy =
0123           HGCalWaferMask::waferXY(partialType_, placeIndex_, waferSize_, 0.0, 0.0, 0.0, v17OrLess_);
0124       for (unsigned int i = 0; i < (wxy.size() - 1); ++i) {
0125         double xp1 = wxy[i].first;
0126         double yp1 = wxy[i].second;
0127         double xp2 = wxy[i + 1].first;
0128         double yp2 = wxy[i + 1].second;
0129         if ((((xi - xp1) / (xp2 - xp1)) - ((yi - yp1) / (yp2 - yp1))) *
0130                 (((xyg.first - xp1) / (xp2 - xp1)) - ((xyg.second - yp1) / (yp2 - yp1))) <=
0131             0) {
0132           goodPoint = false;
0133         }
0134       }
0135       if (goodPoint) {  //Only allowing (x, y) inside a partial wafer 11, placement index 2
0136         std::pair<int32_t, int32_t> uv1 = wafer.cellUVFromXY1(xi, yi, placeIndex_, waferType_, true, false);
0137         std::pair<int32_t, int32_t> uv5 =
0138             wafer.cellUVFromXY1(xi, yi, placeIndex_, waferType_, partialType_, true, false);
0139         std::pair<double, double> xy1 = wafer2.cellUV2XY2(uv5.first, uv5.second, placeIndex_, waferType_);
0140         std::string cellType = (HGCalWaferMask::goodCell(uv5.first, uv5.second, partialType_)) ? "Goodcell" : "Badcell";
0141         std::string pointType = goodPoint ? "GoodPoint" : "BadPoint";
0142         std::string comment =
0143             ((uv1.first != uv5.first) || (uv1.second != uv5.second) || (xy1.first - xi > 6) || (xy1.second - yi > 6))
0144                 ? " ***** ERROR *****"
0145                 : "";
0146         edm::LogVerbatim("HGCalGeom") << pointType << " " << cellType << " x = " << xi << ":" << xy1.first
0147                                       << " y = " << yi << ":" << xy1.second << " type = " << waferType_
0148                                       << " placement index " << placeIndex_ << " u " << uv1.first << ":" << uv5.first
0149                                       << ":" << uv5.first << " v " << uv1.second << ":" << uv5.second << ":"
0150                                       << uv5.second << ":" << comment;
0151       }
0152     }
0153   } else {
0154     for (int i = 0; i < nTrials_; i++) {
0155       int ui = std::floor(2 * nCells * rand() / RAND_MAX);
0156       int vi = std::floor(2 * nCells * rand() / RAND_MAX);
0157       if ((ui < 2 * nCells) && (vi < 2 * nCells) && ((vi - ui) < nCells) && ((ui - vi) <= nCells) &&
0158           HGCalWaferMask::goodCell(ui, vi, partialType_)) {
0159         //Only allowing (U, V) inside a wafer
0160         std::pair<double, double> xy1 = wafer2.cellUV2XY2(ui, vi, placeIndex_, waferType_);
0161         std::pair<int32_t, int32_t> uv1 =
0162             wafer.cellUVFromXY1(xy1.first, xy1.second, placeIndex_, waferType_, true, false);
0163         std::pair<int32_t, int32_t> uv5 =
0164             wafer.cellUVFromXY1(xy1.first, xy1.second, placeIndex_, waferType_, partialType_, true, false);
0165         std::string comment = ((uv1.first != ui) || (uv1.second != vi) || (uv5.first != ui) || (uv5.second != vi))
0166                                   ? " ***** ERROR *****"
0167                                   : "";
0168         edm::LogVerbatim("HGCalGeom") << "u = " << ui << " v = " << vi << " type = " << waferType_
0169                                       << " placement index " << placeIndex_ << " x " << xy1.first << " ,y "
0170                                       << xy1.second << " u " << uv5.first << " v " << uv5.second << comment;
0171       }
0172     }
0173   }
0174   auto end_t = std::chrono::high_resolution_clock::now();
0175   auto diff_t = end_t - start_t;
0176   edm::LogVerbatim("HGCalGeom") << "Execution time for " << nTrials_
0177                                 << " events = " << std::chrono::duration<double, std::milli>(diff_t).count() << " ms";
0178 }
0179 
0180 // define this as a plug-in
0181 DEFINE_FWK_MODULE(HGCalPartialCellTester);