Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalCellPositionTester
0004 // Class:      HGCalCellPositionTester
0005 //
0006 /**\class HGCalCellPositionTester HGCalCellPositionTester.cc
0007  test/HGCalCellPositionTester.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 2022/01/15
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 
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0038 #include "Geometry/HGCalCommonData/interface/HGCalCell.h"
0039 
0040 class HGCalCellPositionTester : public edm::one::EDAnalyzer<> {
0041 public:
0042   explicit HGCalCellPositionTester(const edm::ParameterSet&);
0043   ~HGCalCellPositionTester() override = default;
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   const double waferSize_;
0052   const int waferType_;
0053   const int placeIndex_;
0054 };
0055 
0056 HGCalCellPositionTester::HGCalCellPositionTester(const edm::ParameterSet& iC)
0057     : waferSize_(iC.getParameter<double>("waferSize")),
0058       waferType_(iC.getParameter<int>("waferType")),
0059       placeIndex_(iC.getParameter<int>("cellPlacementIndex")) {
0060   edm::LogVerbatim("HGCalGeom") << "Test positions for wafer of size " << waferSize_ << " Type " << waferType_
0061                                 << " Placement Index " << placeIndex_;
0062 }
0063 
0064 void HGCalCellPositionTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0065   edm::ParameterSetDescription desc;
0066   desc.add<double>("waferSize", 166.4408);
0067   desc.add<int>("waferType", 0);
0068   desc.add<int>("cellPlacementIndex", 7);
0069   descriptions.add("hgcalCellPositionTester", desc);
0070 }
0071 
0072 // ------------ method called to produce the data  ------------
0073 void HGCalCellPositionTester::analyze(const edm::Event&, const edm::EventSetup&) {
0074   const int nFine(12), nCoarse(8);
0075   const double tol(0.00001);
0076   HGCalCell wafer(waferSize_, nFine, nCoarse);
0077   int nCells = (waferType_ == 0) ? nFine : nCoarse;
0078   int indexMin = (placeIndex_ >= 0) ? placeIndex_ : 0;
0079   int indexMax = (placeIndex_ >= 0) ? placeIndex_ : 11;
0080   edm::LogVerbatim("HGCalGeom") << "\nHGCalCellPositionTester:: nCells " << nCells << " and placement index between "
0081                                 << indexMin << " and " << indexMax << "\n\n";
0082   for (int placeIndex = indexMin; placeIndex <= indexMax; ++placeIndex) {
0083     for (int iu = 0; iu < 2 * nCells; ++iu) {
0084       for (int iv = 0; iv < 2 * nCells; ++iv) {
0085         int u(iu), v(iv);
0086         if (placeIndex < HGCalCell::cellPlacementExtra) {
0087           u = iv;
0088           v = iu;
0089         }
0090         if (((v - u) < nCells) && ((u - v) <= nCells)) {
0091           std::pair<double, double> xy1 = wafer.cellUV2XY1(u, v, placeIndex, waferType_);
0092           std::pair<double, double> xy2 = wafer.cellUV2XY2(u, v, placeIndex, waferType_);
0093           double dx = xy1.first - xy2.first;
0094           double dy = xy1.second - xy2.second;
0095           std::string comment = ((std::abs(dx) > tol) || (std::abs(dy) > tol)) ? " ***** ERROR *****" : "";
0096           edm::LogVerbatim("HGCalGeom") << "u = " << u << " v = " << v << " type = " << waferType_
0097                                         << " placement index " << placeIndex << " x " << xy1.first << ":" << xy2.first
0098                                         << ":" << dx << " y " << xy1.second << ":" << xy2.second << ":" << dy
0099                                         << comment;
0100         }
0101       }
0102     }
0103   }
0104 }
0105 
0106 // define this as a plug-in
0107 DEFINE_FWK_MODULE(HGCalCellPositionTester);