Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <array>
0002 #include <iostream>
0003 #include <sstream>
0004 #include <map>
0005 
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 
0014 #include "Geometry/HGCalTBCommonData/interface/HGCalTBParameters.h"
0015 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0016 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0017 
0018 class HGCalTBParameterTester : public edm::one::EDAnalyzer<> {
0019 public:
0020   explicit HGCalTBParameterTester(const edm::ParameterSet&);
0021   ~HGCalTBParameterTester() override = default;
0022   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0023 
0024   void beginJob() override {}
0025   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0026   void endJob() override {}
0027 
0028 private:
0029   template <typename T>
0030   void myPrint(std::string const& s, std::vector<T> const& obj, int n) const;
0031   void myPrint(std::string const& s, std::vector<double> const& obj1, std::vector<double> const& obj2, int n) const;
0032   void printTrform(HGCalTBParameters const*) const;
0033   void printWaferType(HGCalTBParameters const* phgp) const;
0034 
0035   const std::string name_;
0036   edm::ESGetToken<HGCalTBParameters, IdealGeometryRecord> token_;
0037   const int mode_;
0038 };
0039 
0040 HGCalTBParameterTester::HGCalTBParameterTester(const edm::ParameterSet& ic)
0041     : name_(ic.getParameter<std::string>("name")),
0042       token_(esConsumes<HGCalTBParameters, IdealGeometryRecord>(edm::ESInputTag{"", name_})),
0043       mode_(ic.getParameter<int>("mode")) {}
0044 
0045 void HGCalTBParameterTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046   edm::ParameterSetDescription desc;
0047   desc.add<std::string>("name", "HGCalEESensitive");
0048   desc.add<int>("mode", 0);
0049   descriptions.add("hgcTBParameterTesterEE", desc);
0050 }
0051 
0052 void HGCalTBParameterTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0053   edm::LogVerbatim("HGCalGeomr") << "HGCalTBParameter::Here I am";
0054   auto start = std::chrono::high_resolution_clock::now();
0055 
0056   const auto& hgp = iSetup.getData(token_);
0057   const auto* phgp = &hgp;
0058 
0059   edm::LogVerbatim("HGCalGeom") << phgp->name_;
0060 
0061   // Wafers of 6-inch format
0062   edm::LogVerbatim("HGCalGeom") << "DetectorType: " << phgp->detectorType_;
0063   edm::LogVerbatim("HGCalGeom") << "WaferR_: " << phgp->waferR_;
0064   edm::LogVerbatim("HGCalGeom") << "nCells_: " << phgp->nCells_;
0065   edm::LogVerbatim("HGCalGeom") << "nSectors_: " << phgp->nSectors_;
0066   edm::LogVerbatim("HGCalGeom") << "FirstLayer: " << phgp->firstLayer_;
0067   edm::LogVerbatim("HGCalGeom") << "mode_: " << phgp->mode_;
0068 
0069   myPrint("CellSize", phgp->cellSize_, 10);
0070   myPrint("slopeMin", phgp->slopeMin_, 10);
0071   myPrint("slopeTop", phgp->slopeTop_, 10);
0072   myPrint("zFrontTop", phgp->zFrontTop_, 10);
0073   myPrint("rMaxFront", phgp->rMaxFront_, 10);
0074   myPrint("zRanges", phgp->zRanges_, 10);
0075   myPrint("moduleBlS", phgp->moduleBlS_, 10);
0076   myPrint("moduleTlS", phgp->moduleTlS_, 10);
0077   myPrint("moduleHS", phgp->moduleHS_, 10);
0078   myPrint("moduleDzS", phgp->moduleDzS_, 10);
0079   myPrint("moduleAlphaS", phgp->moduleAlphaS_, 10);
0080   myPrint("moduleCellS", phgp->moduleCellS_, 10);
0081   myPrint("moduleBlR", phgp->moduleBlR_, 10);
0082   myPrint("moduleTlR", phgp->moduleTlR_, 10);
0083   myPrint("moduleHR", phgp->moduleHR_, 10);
0084   myPrint("moduleDzR", phgp->moduleDzR_, 10);
0085   myPrint("moduleAlphaR", phgp->moduleAlphaR_, 10);
0086   myPrint("moduleCellR", phgp->moduleCellR_, 10);
0087   myPrint("trformTranX", phgp->trformTranX_, 10);
0088   myPrint("trformTranY", phgp->trformTranY_, 10);
0089   myPrint("trformTranZ", phgp->trformTranZ_, 10);
0090   myPrint("trformRotXX", phgp->trformRotXX_, 10);
0091   myPrint("trformRotYX", phgp->trformRotYX_, 10);
0092   myPrint("trformRotZX", phgp->trformRotZX_, 10);
0093   myPrint("trformRotXY", phgp->trformRotXY_, 10);
0094   myPrint("trformRotYY", phgp->trformRotYY_, 10);
0095   myPrint("trformRotZY", phgp->trformRotZY_, 10);
0096   myPrint("trformRotXZ", phgp->trformRotXZ_, 10);
0097   myPrint("trformRotYZ", phgp->trformRotYZ_, 10);
0098   myPrint("trformRotZZ", phgp->trformRotZZ_, 10);
0099   myPrint("zLayerHex", phgp->zLayerHex_, 10);
0100   myPrint("rMinLayHex", phgp->rMinLayHex_, 10);
0101   myPrint("rMaxLayHex", phgp->rMaxLayHex_, 10);
0102   myPrint("waferPos", phgp->waferPosX_, phgp->waferPosY_, 4);
0103   myPrint("cellFine", phgp->cellFineX_, phgp->cellFineY_, 4);
0104   myPrint("cellFineHalf", phgp->cellFineHalf_, 10);
0105   myPrint("cellCoarse", phgp->cellCoarseX_, phgp->cellCoarseY_, 4);
0106   myPrint("cellCoarseHalf", phgp->cellCoarseHalf_, 10);
0107   myPrint("boundR", phgp->boundR_, 10);
0108   myPrint("moduleLayS", phgp->moduleLayS_, 10);
0109   myPrint("moduleLayR", phgp->moduleLayR_, 10);
0110   myPrint("layer", phgp->layer_, 18);
0111   myPrint("layerIndex", phgp->layerIndex_, 18);
0112   myPrint("layerGroup", phgp->layerGroup_, 18);
0113   myPrint("cellFactor", phgp->cellFactor_, 10);
0114   myPrint("depth", phgp->depth_, 18);
0115   myPrint("depthIndex", phgp->depthIndex_, 18);
0116   myPrint("depthLayerF", phgp->depthLayerF_, 18);
0117   myPrint("waferCopy", phgp->waferCopy_, 10);
0118   myPrint("waferTypeL", phgp->waferTypeL_, 25);
0119   myPrint("waferTypeT", phgp->waferTypeT_, 25);
0120   myPrint("layerGroupM", phgp->layerGroupM_, 18);
0121   myPrint("layerGroupO", phgp->layerGroupO_, 18);
0122   printTrform(phgp);
0123   myPrint("levelTop", phgp->levelT_, 10);
0124   printWaferType(phgp);
0125 
0126   auto finish = std::chrono::high_resolution_clock::now();
0127   std::chrono::duration<double> elapsed = finish - start;
0128   edm::LogVerbatim("HGCalGeom") << "Elapsed time: " << elapsed.count() << " s";
0129 }
0130 
0131 template <typename T>
0132 void HGCalTBParameterTester::myPrint(std::string const& s, std::vector<T> const& obj, int n) const {
0133   int k(0), kk(0);
0134   edm::LogVerbatim("HGCalGeom") << s << " with " << obj.size() << " elements with n " << n << ": 1000";
0135   std::ostringstream st1[1000];
0136   for (auto const& it : obj) {
0137     st1[kk] << it << ", ";
0138     ++k;
0139     if (k == n) {
0140       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0141       ++kk;
0142       k = 0;
0143     }
0144   }
0145   if (k > 0)
0146     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0147 }
0148 
0149 void HGCalTBParameterTester::myPrint(std::string const& s,
0150                                      std::vector<double> const& obj1,
0151                                      std::vector<double> const& obj2,
0152                                      int n) const {
0153   int k(0), kk(0);
0154   std::ostringstream st1[250];
0155   edm::LogVerbatim("HGCalGeom") << s << " with " << obj1.size() << " elements with n " << n << ": 250";
0156   for (unsigned int k1 = 0; k1 < obj1.size(); ++k1) {
0157     st1[kk] << "(" << obj1[k1] << ", " << obj2[k1] << ") ";
0158     ++k;
0159     if (k == n) {
0160       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0161       ++kk;
0162       k = 0;
0163     }
0164   }
0165   if (k > 0)
0166     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0167 }
0168 
0169 void HGCalTBParameterTester::printTrform(HGCalTBParameters const* phgp) const {
0170   int k(0), kk(0);
0171   std::ostringstream st1[20];
0172   edm::LogVerbatim("HGCalGeom") << "TrformIndex with " << phgp->trformIndex_.size() << " elements with n 7:20";
0173   for (unsigned int i = 0; i < phgp->trformIndex_.size(); ++i) {
0174     std::array<int, 4> id = phgp->getID(i);
0175     st1[kk] << id[0] << ":" << id[1] << ":" << id[2] << ":" << id[3] << ", ";
0176     ++k;
0177     if (k == 7) {
0178       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0179       ++kk;
0180       k = 0;
0181     }
0182   }
0183   if (k > 0)
0184     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0185 }
0186 
0187 void HGCalTBParameterTester::printWaferType(HGCalTBParameters const* phgp) const {
0188   int k(0);
0189   edm::LogVerbatim("HGCalGeom") << "waferTypes with " << phgp->waferTypes_.size() << " elements";
0190   std::map<std::pair<int, int>, int> kounts;
0191   std::map<std::pair<int, int>, int>::iterator itr;
0192   for (auto const& it : phgp->waferTypes_) {
0193     std::ostringstream st1;
0194     st1 << " [" << k << "] " << HGCalWaferIndex::waferLayer(it.first);
0195     if (HGCalWaferIndex::waferFormat(it.first)) {
0196       st1 << ":" << HGCalWaferIndex::waferU(it.first) << ":" << HGCalWaferIndex::waferV(it.first);
0197     } else {
0198       st1 << ":" << HGCalWaferIndex::waferCopy(it.first);
0199     }
0200     edm::LogVerbatim("HGCalGeom") << st1.str() << " ==> (" << (it.second).first << ":" << (it.second).second << ")";
0201     itr = kounts.find(it.second);
0202     if (itr == kounts.end())
0203       kounts[it.second] = 1;
0204     else
0205       ++(itr->second);
0206     ++k;
0207   }
0208   if (!kounts.empty()) {
0209     edm::LogVerbatim("HGCalGeom") << "Summary of waferTypes ==========================";
0210     for (itr = kounts.begin(); itr != kounts.end(); ++itr)
0211       edm::LogVerbatim("HGCalGeom") << "Type (" << (itr->first).first << ":" << (itr->first).second << ") Kount "
0212                                     << itr->second;
0213   }
0214 }
0215 
0216 DEFINE_FWK_MODULE(HGCalTBParameterTester);