Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    HGCalWaferTester
0004 // Class:      HGCalWaferTester
0005 //
0006 /**\class HGCalWaferTester HGCalWaferTester.cc
0007  test/HGCalWaferTester.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 2019/07/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/HGCalDDDConstants.h"
0039 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0040 
0041 class HGCalWaferTester : public edm::one::EDAnalyzer<> {
0042 public:
0043   explicit HGCalWaferTester(const edm::ParameterSet&);
0044   ~HGCalWaferTester() override = default;
0045   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0046 
0047   void beginJob() override {}
0048   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0049   void endJob() override {}
0050 
0051 private:
0052   const std::string nameSense_, nameDetector_;
0053   const bool reco_;
0054   edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord> dddToken_;
0055 };
0056 
0057 HGCalWaferTester::HGCalWaferTester(const edm::ParameterSet& iC)
0058     : nameSense_(iC.getParameter<std::string>("NameSense")),
0059       nameDetector_(iC.getParameter<std::string>("NameDevice")),
0060       reco_(iC.getParameter<bool>("Reco")) {
0061   dddToken_ = esConsumes<HGCalDDDConstants, IdealGeometryRecord>(edm::ESInputTag{"", nameSense_});
0062 
0063   edm::LogVerbatim("HGCalGeom") << "Test numbering for " << nameDetector_ << " using constants of " << nameSense_
0064                                 << " for  RecoFlag " << reco_;
0065 }
0066 
0067 void HGCalWaferTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0068   edm::ParameterSetDescription desc;
0069   desc.add<std::string>("NameSense", "HGCalEESensitive");
0070   desc.add<std::string>("NameDevice", "HGCal EE");
0071   desc.add<bool>("Reco", false);
0072   descriptions.add("hgcalWaferTesterEE", desc);
0073 }
0074 
0075 // ------------ method called to produce the data  ------------
0076 void HGCalWaferTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0077   const HGCalDDDConstants& hgdc = iSetup.getData(dddToken_);
0078   edm::LogVerbatim("HGCalGeom") << nameDetector_ << " Layers = " << hgdc.layers(reco_)
0079                                 << " Sectors = " << hgdc.sectors() << std::endl;
0080   if (hgdc.waferHexagon8()) {
0081     int layer = hgdc.firstLayer();
0082     for (int u = -12; u <= 12; ++u) {
0083       std::pair<double, double> xy = hgdc.waferPosition(layer, u, 0, reco_, false);
0084       edm::LogVerbatim("HGCalGeom") << " iz = +, u = " << u << ", v = 0: x = " << xy.first << " y = " << xy.second
0085                                     << "\n"
0086                                     << " iz = -, u = " << u << ", v = 0: x = " << -xy.first << " y = " << xy.second;
0087     }
0088     for (int v = -12; v <= 12; ++v) {
0089       std::pair<double, double> xy = hgdc.waferPosition(layer, 0, v, reco_, false);
0090       edm::LogVerbatim("HGCalGeom") << " iz = +, u = 0, v = " << v << ": x = " << xy.first << " y = " << xy.second
0091                                     << "\n"
0092                                     << " iz = -, u = 0, v = " << v << ": x = " << -xy.first << " y = " << xy.second;
0093     }
0094   }
0095 }
0096 
0097 // define this as a plug-in
0098 DEFINE_FWK_MODULE(HGCalWaferTester);