Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-02 02:40:41

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/HGCalCommonData/interface/HGCalParameters.h"
0015 #include "Geometry/HGCalCommonData/interface/HGCalTileIndex.h"
0016 #include "Geometry/HGCalCommonData/interface/HGCalWaferIndex.h"
0017 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0018 
0019 class HGCalParameterTester : public edm::one::EDAnalyzer<> {
0020 public:
0021   explicit HGCalParameterTester(const edm::ParameterSet&);
0022   ~HGCalParameterTester() override = default;
0023   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024 
0025   void beginJob() override {}
0026   void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0027   void endJob() override {}
0028 
0029 private:
0030   template <typename T>
0031   void myPrint(std::string const& s, std::vector<T> const& obj, int n) const;
0032   template <typename T>
0033   void myPrint(std::string const& s, std::vector<std::pair<T, T> > const& obj, int n) const;
0034   void myPrint(std::string const& s, std::vector<double> const& obj1, std::vector<double> const& obj2, int n) const;
0035   void myPrint(std::string const& s, HGCalParameters::wafer_map const& obj, int n) const;
0036   void printTrform(HGCalParameters const*) const;
0037   void printWaferType(HGCalParameters const* phgp) const;
0038 
0039   const std::string name_;
0040   edm::ESGetToken<HGCalParameters, IdealGeometryRecord> token_;
0041   const int mode_;
0042 };
0043 
0044 HGCalParameterTester::HGCalParameterTester(const edm::ParameterSet& ic)
0045     : name_(ic.getParameter<std::string>("Name")),
0046       token_(esConsumes<HGCalParameters, IdealGeometryRecord>(edm::ESInputTag{"", name_})),
0047       mode_(ic.getParameter<int>("Mode")) {}
0048 
0049 void HGCalParameterTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050   edm::ParameterSetDescription desc;
0051   desc.add<std::string>("Name", "HGCalEESensitive");
0052   desc.add<int>("Mode", 1);
0053   descriptions.add("hgcParameterTesterEE", desc);
0054 }
0055 
0056 void HGCalParameterTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0057   edm::LogVerbatim("HGCalGeomr") << "HGCalParameter::Here I am";
0058   auto start = std::chrono::high_resolution_clock::now();
0059 
0060   const auto& hgp = iSetup.getData(token_);
0061   const auto* phgp = &hgp;
0062 
0063   edm::LogVerbatim("HGCalGeom") << phgp->name_;
0064   if (mode_ == 0) {
0065     // Wafers of 6-inch format
0066     edm::LogVerbatim("HGCalGeom") << "DetectorType: " << phgp->detectorType_;
0067     edm::LogVerbatim("HGCalGeom") << "WaferR_: " << phgp->waferR_;
0068     edm::LogVerbatim("HGCalGeom") << "nCells_: " << phgp->nCells_;
0069     edm::LogVerbatim("HGCalGeom") << "nSectors_: " << phgp->nSectors_;
0070     edm::LogVerbatim("HGCalGeom") << "FirstLayer: " << phgp->firstLayer_;
0071     edm::LogVerbatim("HGCalGeom") << "FirstMixedLayer: " << phgp->firstMixedLayer_;
0072     edm::LogVerbatim("HGCalGeom") << "mode_: " << phgp->mode_;
0073 
0074     myPrint("CellSize", phgp->cellSize_, 10);
0075     myPrint("slopeMin", phgp->slopeMin_, 10);
0076     myPrint("slopeTop", phgp->slopeTop_, 10);
0077     myPrint("zFrontTop", phgp->zFrontTop_, 10);
0078     myPrint("rMaxFront", phgp->rMaxFront_, 10);
0079     myPrint("zRanges", phgp->zRanges_, 10);
0080     myPrint("moduleBlS", phgp->moduleBlS_, 10);
0081     myPrint("moduleTlS", phgp->moduleTlS_, 10);
0082     myPrint("moduleHS", phgp->moduleHS_, 10);
0083     myPrint("moduleDzS", phgp->moduleDzS_, 10);
0084     myPrint("moduleAlphaS", phgp->moduleAlphaS_, 10);
0085     myPrint("moduleCellS", phgp->moduleCellS_, 10);
0086     myPrint("moduleBlR", phgp->moduleBlR_, 10);
0087     myPrint("moduleTlR", phgp->moduleTlR_, 10);
0088     myPrint("moduleHR", phgp->moduleHR_, 10);
0089     myPrint("moduleDzR", phgp->moduleDzR_, 10);
0090     myPrint("moduleAlphaR", phgp->moduleAlphaR_, 10);
0091     myPrint("moduleCellR", phgp->moduleCellR_, 10);
0092     myPrint("trformTranX", phgp->trformTranX_, 10);
0093     myPrint("trformTranY", phgp->trformTranY_, 10);
0094     myPrint("trformTranZ", phgp->trformTranZ_, 10);
0095     myPrint("trformRotXX", phgp->trformRotXX_, 10);
0096     myPrint("trformRotYX", phgp->trformRotYX_, 10);
0097     myPrint("trformRotZX", phgp->trformRotZX_, 10);
0098     myPrint("trformRotXY", phgp->trformRotXY_, 10);
0099     myPrint("trformRotYY", phgp->trformRotYY_, 10);
0100     myPrint("trformRotZY", phgp->trformRotZY_, 10);
0101     myPrint("trformRotXZ", phgp->trformRotXZ_, 10);
0102     myPrint("trformRotYZ", phgp->trformRotYZ_, 10);
0103     myPrint("trformRotZZ", phgp->trformRotZZ_, 10);
0104     myPrint("zLayerHex", phgp->zLayerHex_, 10);
0105     myPrint("rMinLayHex", phgp->rMinLayHex_, 10);
0106     myPrint("rMaxLayHex", phgp->rMaxLayHex_, 10);
0107     myPrint("waferPos", phgp->waferPosX_, phgp->waferPosY_, 4);
0108     myPrint("cellFine", phgp->cellFineX_, phgp->cellFineY_, 4);
0109     myPrint("cellFineHalf", phgp->cellFineHalf_, 10);
0110     myPrint("cellCoarse", phgp->cellCoarseX_, phgp->cellCoarseY_, 4);
0111     myPrint("cellCoarseHalf", phgp->cellCoarseHalf_, 10);
0112     myPrint("boundR", phgp->boundR_, 10);
0113     myPrint("moduleLayS", phgp->moduleLayS_, 10);
0114     myPrint("moduleLayR", phgp->moduleLayR_, 10);
0115     myPrint("layer", phgp->layer_, 18);
0116     myPrint("layerIndex", phgp->layerIndex_, 18);
0117     myPrint("layerGroup", phgp->layerGroup_, 18);
0118     myPrint("cellFactor", phgp->cellFactor_, 10);
0119     myPrint("depth", phgp->depth_, 18);
0120     myPrint("depthIndex", phgp->depthIndex_, 18);
0121     myPrint("depthLayerF", phgp->depthLayerF_, 18);
0122     myPrint("waferCopy", phgp->waferCopy_, 10);
0123     myPrint("waferTypeL", phgp->waferTypeL_, 25);
0124     myPrint("waferTypeT", phgp->waferTypeT_, 25);
0125     myPrint("layerGroupM", phgp->layerGroupM_, 18);
0126     myPrint("layerGroupO", phgp->layerGroupO_, 18);
0127     printTrform(phgp);
0128     myPrint("levelTop", phgp->levelT_, 10);
0129     printWaferType(phgp);
0130 
0131   } else if (mode_ == 1) {
0132     // Wafers of 8-inch format
0133     edm::LogVerbatim("HGCalGeom") << "DetectorType: " << phgp->detectorType_;
0134     edm::LogVerbatim("HGCalGeom") << "UseSimWt: " << phgp->useSimWt_;
0135     edm::LogVerbatim("HGCalGeom") << "Wafer Parameters: " << phgp->waferSize_ << ":" << phgp->waferR_ << ":"
0136                                   << phgp->waferThick_ << ":" << phgp->sensorSeparation_ << ":"
0137                                   << phgp->sensorSizeOffset_ << ":" << phgp->guardRingOffset_ << ":" << phgp->mouseBite_
0138                                   << ":" << phgp->useOffset_;
0139     myPrint("waferThickness", phgp->waferThickness_, 10);
0140     edm::LogVerbatim("HGCalGeom") << "nCells_: " << phgp->nCellsFine_ << ":" << phgp->nCellsCoarse_;
0141     edm::LogVerbatim("HGCalGeom") << "nSectors_: " << phgp->nSectors_;
0142     edm::LogVerbatim("HGCalGeom") << "FirstLayer: " << phgp->firstLayer_;
0143     edm::LogVerbatim("HGCalGeom") << "FirstMixedLayer: " << phgp->firstMixedLayer_;
0144     edm::LogVerbatim("HGCalGeom") << "LayerOffset: " << phgp->layerOffset_;
0145     edm::LogVerbatim("HGCalGeom") << "mode_: " << phgp->mode_;
0146     edm::LogVerbatim("HGCalGeom") << "cassettes_: " << phgp->cassettes_;
0147 
0148     edm::LogVerbatim("HGCalGeom") << "waferUVMax: " << phgp->waferUVMax_;
0149     myPrint("waferUVMaxLayer", phgp->waferUVMaxLayer_, 20);
0150     myPrint("CellThickness", phgp->cellThickness_, 10);
0151     myPrint("radius100to200", phgp->radius100to200_, 10);
0152     myPrint("radius200to300", phgp->radius200to300_, 10);
0153     edm::LogVerbatim("HGCalGeom") << "choiceType " << phgp->choiceType_ << "   nCornerCut " << phgp->nCornerCut_
0154                                   << "  fracAreaMin " << phgp->fracAreaMin_ << "  zMinForRad " << phgp->zMinForRad_;
0155 
0156     myPrint("CellSize", phgp->cellSize_, 10);
0157     myPrint("radiusMixBoundary", phgp->radiusMixBoundary_, 10);
0158     myPrint("LayerCenter", phgp->layerCenter_, 20);
0159     myPrint("LayerType", phgp->layerType_, 20);
0160     edm::LogVerbatim("HGCalGeom") << "Layer Rotation " << phgp->layerRotation_ << "   with " << phgp->layerRotV_.size()
0161                                   << "  parameters";
0162     for (unsigned int k = 0; k < phgp->layerRotV_.size(); ++k)
0163       edm::LogVerbatim("HGCalGeom") << "Element[" << k << "] " << phgp->layerRotV_[k].first << ":"
0164                                     << phgp->layerRotV_[k].second;
0165     edm::LogVerbatim("HGCalGeom") << "CalibCellRadiusHD " << phgp->calibCellRHD_;
0166     myPrint("CalibCellFullHD", phgp->calibCellFullHD_, 12);
0167     myPrint("CalibCellPartHD", phgp->calibCellPartHD_, 12);
0168     edm::LogVerbatim("HGCalGeom") << "CalibCellRadiusHD " << phgp->calibCellRLD_;
0169     myPrint("CalibCellFullLD", phgp->calibCellFullLD_, 12);
0170     myPrint("CalibCellPartLD", phgp->calibCellPartLD_, 12);
0171     myPrint("cassetteShift", phgp->cassetteShift_, 8);
0172     myPrint("cassetteShiftTile", phgp->cassetteShiftTile_, 8);
0173     myPrint("cassetteRetractTile", phgp->cassetteRetractTile_, 8);
0174 
0175     edm::LogVerbatim("HGCalGeom") << "MaskMode: " << phgp->waferMaskMode_;
0176     if (phgp->waferMaskMode_ > 1) {
0177       edm::LogVerbatim("HGCalGeom") << "WaferInfo with " << phgp->waferInfoMap_.size() << " elements";
0178       unsigned int kk(0);
0179       std::unordered_map<int32_t, HGCalParameters::waferInfo>::const_iterator itr = phgp->waferInfoMap_.begin();
0180       for (; itr != phgp->waferInfoMap_.end(); ++itr, ++kk)
0181         edm::LogVerbatim("HGCalGeom") << "[" << kk << "] " << itr->first << "["
0182                                       << HGCalWaferIndex::waferLayer(itr->first) << ", "
0183                                       << HGCalWaferIndex::waferU(itr->first) << ", "
0184                                       << HGCalWaferIndex::waferV(itr->first) << "] (" << (itr->second).type << ", "
0185                                       << (itr->second).part << ", " << (itr->second).orient << ")";
0186     }
0187 
0188     myPrint("slopeMin", phgp->slopeMin_, 10);
0189     myPrint("zFrontMin", phgp->zFrontMin_, 10);
0190     myPrint("rMinFront", phgp->rMinFront_, 10);
0191     myPrint("slopeTop", phgp->slopeTop_, 10);
0192     myPrint("zFrontTop", phgp->zFrontTop_, 10);
0193     myPrint("rMaxFront", phgp->rMaxFront_, 10);
0194     myPrint("zRanges", phgp->zRanges_, 10);
0195     myPrint("moduleBlS", phgp->moduleBlS_, 10);
0196     myPrint("moduleTlS", phgp->moduleTlS_, 10);
0197     myPrint("moduleHS", phgp->moduleHS_, 10);
0198     myPrint("moduleDzS", phgp->moduleDzS_, 10);
0199     myPrint("moduleAlphaS", phgp->moduleAlphaS_, 10);
0200     myPrint("moduleCellS", phgp->moduleCellS_, 10);
0201     myPrint("moduleBlR", phgp->moduleBlR_, 10);
0202     myPrint("moduleTlR", phgp->moduleTlR_, 10);
0203     myPrint("moduleHR", phgp->moduleHR_, 10);
0204     myPrint("moduleDzR", phgp->moduleDzR_, 10);
0205     myPrint("moduleAlphaR", phgp->moduleAlphaR_, 10);
0206     myPrint("moduleCellR", phgp->moduleCellR_, 10);
0207     myPrint("trformTranX", phgp->trformTranX_, 8);
0208     myPrint("trformTranY", phgp->trformTranY_, 8);
0209     myPrint("trformTranZ", phgp->trformTranZ_, 8);
0210     myPrint("trformRotXX", phgp->trformRotXX_, 10);
0211     myPrint("trformRotYX", phgp->trformRotYX_, 10);
0212     myPrint("trformRotZX", phgp->trformRotZX_, 10);
0213     myPrint("trformRotXY", phgp->trformRotXY_, 10);
0214     myPrint("trformRotYY", phgp->trformRotYY_, 10);
0215     myPrint("trformRotZY", phgp->trformRotZY_, 10);
0216     myPrint("trformRotXZ", phgp->trformRotXZ_, 10);
0217     myPrint("trformRotYZ", phgp->trformRotYZ_, 10);
0218     myPrint("trformRotZZ", phgp->trformRotZZ_, 10);
0219     myPrint("xLayerHex", phgp->xLayerHex_, 8);
0220     myPrint("yLayerHex", phgp->yLayerHex_, 8);
0221     myPrint("zLayerHex", phgp->zLayerHex_, 8);
0222     myPrint("rMinLayHex", phgp->rMinLayHex_, 8);
0223     myPrint("rMaxLayHex", phgp->rMaxLayHex_, 8);
0224     myPrint("waferPos", phgp->waferPosX_, phgp->waferPosY_, 4);
0225     myPrint("cellFineIndex", phgp->cellFineIndex_, 8);
0226     myPrint("cellFine", phgp->cellFineX_, phgp->cellFineY_, 4);
0227     myPrint("cellCoarseIndex", phgp->cellCoarseIndex_, 8);
0228     myPrint("cellCoarse", phgp->cellCoarseX_, phgp->cellCoarseY_, 4);
0229     myPrint("layer", phgp->layer_, 18);
0230     myPrint("layerIndex", phgp->layerIndex_, 18);
0231     myPrint("depth", phgp->depth_, 18);
0232     myPrint("depthIndex", phgp->depthIndex_, 18);
0233     myPrint("depthLayerF", phgp->depthLayerF_, 18);
0234     myPrint("waferCopy", phgp->waferCopy_, 10);
0235     myPrint("waferTypeL", phgp->waferTypeL_, 25);
0236     printTrform(phgp);
0237     myPrint("levelTop", phgp->levelT_, 10);
0238     printWaferType(phgp);
0239   } else {
0240     // Tpaezoid (scintillator) type
0241     edm::LogVerbatim("HGCalGeom") << "DetectorType: " << phgp->detectorType_;
0242     edm::LogVerbatim("HGCalGeom") << "UseSimWt: " << phgp->useSimWt_;
0243     edm::LogVerbatim("HGCalGeom") << "nCells_: " << phgp->nCellsFine_ << ":" << phgp->nCellsCoarse_;
0244     edm::LogVerbatim("HGCalGeom") << "MinTileZize: " << phgp->minTileSize_;
0245     edm::LogVerbatim("HGCalGeom") << "FirstLayer: " << phgp->firstLayer_;
0246     edm::LogVerbatim("HGCalGeom") << "FirstMixedLayer: " << phgp->firstMixedLayer_;
0247     edm::LogVerbatim("HGCalGeom") << "LayerOffset: " << phgp->layerOffset_;
0248     edm::LogVerbatim("HGCalGeom") << "mode_: " << phgp->mode_;
0249     edm::LogVerbatim("HGCalGeom") << "waferUVMax: " << phgp->waferUVMax_;
0250     edm::LogVerbatim("HGCalGeom") << "nSectors_: " << phgp->nSectors_;
0251     edm::LogVerbatim("HGCalGeom") << "nCells_: " << phgp->nCellsFine_ << ":" << phgp->nCellsCoarse_;
0252 
0253     myPrint("CellSize", phgp->cellSize_, 10);
0254     myPrint("radiusMixBoundary", phgp->radiusMixBoundary_, 10);
0255     myPrint("nPhiBinBH", phgp->nPhiBinBH_, 18);
0256     myPrint("layerFrontBH", phgp->layerFrontBH_, 10);
0257     myPrint("LayerCenter", phgp->layerCenter_, 20);
0258     myPrint("rMinLayerBH", phgp->rMinLayerBH_, 10);
0259     myPrint("slopeMin", phgp->slopeMin_, 10);
0260     myPrint("zFrontMin", phgp->zFrontMin_, 10);
0261     myPrint("rMinFront", phgp->rMinFront_, 10);
0262     myPrint("radiusLayer[0]", phgp->radiusLayer_[0], 10);
0263     myPrint("radiusLayer[1]", phgp->radiusLayer_[1], 10);
0264     myPrint("iradMinBH", phgp->iradMinBH_, 20);
0265     myPrint("iradMaxBH", phgp->iradMaxBH_, 20);
0266     edm::LogVerbatim("HGCalGeom") << "MaskMode: " << phgp->waferMaskMode_;
0267     if (phgp->waferMaskMode_ > 1) {
0268       myPrint("tileRingR", phgp->tileRingR_, 4);
0269       myPrint("tileRingRange", phgp->tileRingRange_, 8);
0270       edm::LogVerbatim("HGCalGeom") << "TileInfo with " << phgp->tileInfoMap_.size() << " elements";
0271       unsigned int kk(0);
0272       std::unordered_map<int32_t, HGCalParameters::tileInfo>::const_iterator itr = phgp->tileInfoMap_.begin();
0273       for (; itr != phgp->tileInfoMap_.end(); ++itr, ++kk)
0274         edm::LogVerbatim("HGCalGeom") << "[" << kk << "] " << itr->first << "[" << HGCalTileIndex::tileLayer(itr->first)
0275                                       << ", " << HGCalTileIndex::tileRing(itr->first) << ", "
0276                                       << HGCalTileIndex::tilePhi(itr->first) << "] (" << (itr->second).type << ", "
0277                                       << (itr->second).sipm << std::hex << ", " << (itr->second).hex[0] << ", "
0278                                       << (itr->second).hex[1] << ", " << (itr->second).hex[2] << ", "
0279                                       << (itr->second).hex[3] << ", " << (itr->second).hex[4] << ", "
0280                                       << (itr->second).hex[5] << ")" << std::dec;
0281     }
0282 
0283     myPrint("slopeTop", phgp->slopeTop_, 10);
0284     myPrint("zFrontTop", phgp->zFrontTop_, 10);
0285     myPrint("rMaxFront", phgp->rMaxFront_, 10);
0286     myPrint("zRanges", phgp->zRanges_, 10);
0287     myPrint("firstModule", phgp->firstModule_, 10);
0288     myPrint("lastModule", phgp->lastModule_, 10);
0289     myPrint("moduleBlS", phgp->moduleBlS_, 10);
0290     myPrint("moduleTlS", phgp->moduleTlS_, 10);
0291     myPrint("moduleHS", phgp->moduleHS_, 10);
0292     myPrint("moduleDzS", phgp->moduleDzS_, 10);
0293     myPrint("moduleAlphaS", phgp->moduleAlphaS_, 10);
0294     myPrint("moduleCellS", phgp->moduleCellS_, 10);
0295     myPrint("moduleBlR", phgp->moduleBlR_, 10);
0296     myPrint("moduleTlR", phgp->moduleTlR_, 10);
0297     myPrint("moduleHR", phgp->moduleHR_, 10);
0298     myPrint("moduleDzR", phgp->moduleDzR_, 10);
0299     myPrint("moduleAlphaR", phgp->moduleAlphaR_, 10);
0300     myPrint("moduleCellR", phgp->moduleCellR_, 9);
0301     myPrint("trformTranX", phgp->trformTranY_, 9);
0302     myPrint("trformTranY", phgp->trformTranY_, 9);
0303     myPrint("trformTranZ", phgp->trformTranZ_, 9);
0304     myPrint("trformRotXX", phgp->trformRotXX_, 10);
0305     myPrint("trformRotYX", phgp->trformRotYX_, 10);
0306     myPrint("trformRotZX", phgp->trformRotZX_, 10);
0307     myPrint("trformRotXY", phgp->trformRotXY_, 10);
0308     myPrint("trformRotYY", phgp->trformRotYY_, 10);
0309     myPrint("trformRotZY", phgp->trformRotZY_, 10);
0310     myPrint("trformRotXZ", phgp->trformRotXZ_, 10);
0311     myPrint("trformRotYZ", phgp->trformRotYZ_, 10);
0312     myPrint("trformRotZZ", phgp->trformRotZZ_, 10);
0313     myPrint("xLayerHex", phgp->xLayerHex_, 10);
0314     myPrint("yLayerHex", phgp->yLayerHex_, 10);
0315     myPrint("zLayerHex", phgp->zLayerHex_, 10);
0316     myPrint("rMinLayHex", phgp->rMinLayHex_, 9);
0317     myPrint("rMaxLayHex", phgp->rMaxLayHex_, 9);
0318     myPrint("layer", phgp->layer_, 18);
0319     myPrint("layerIndex", phgp->layerIndex_, 18);
0320     myPrint("depth", phgp->depth_, 18);
0321     myPrint("depthIndex", phgp->depthIndex_, 18);
0322     myPrint("depthLayerF", phgp->depthLayerF_, 18);
0323     printTrform(phgp);
0324     myPrint("levelTop", phgp->levelT_, 10);
0325     printWaferType(phgp);
0326   }
0327 
0328   auto finish = std::chrono::high_resolution_clock::now();
0329   std::chrono::duration<double> elapsed = finish - start;
0330   edm::LogVerbatim("HGCalGeom") << "Elapsed time: " << elapsed.count() << " s";
0331 }
0332 
0333 template <typename T>
0334 void HGCalParameterTester::myPrint(std::string const& s, std::vector<T> const& obj, int n) const {
0335   int k(0), kk(0);
0336   edm::LogVerbatim("HGCalGeom") << s << " with " << obj.size() << " elements with n " << n << ": 1000";
0337   std::ostringstream st1[1000];
0338   for (auto const& it : obj) {
0339     st1[kk] << it << ", ";
0340     ++k;
0341     if (k == n) {
0342       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0343       ++kk;
0344       k = 0;
0345     }
0346   }
0347   if (k > 0)
0348     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0349 }
0350 
0351 template <typename T>
0352 void HGCalParameterTester::myPrint(std::string const& s, std::vector<std::pair<T, T> > const& obj, int n) const {
0353   int k(0), kk(0);
0354   edm::LogVerbatim("HGCalGeom") << s << " with " << obj.size() << " elements with n " << n << ":200";
0355   std::ostringstream st1[200];
0356   for (auto const& it : obj) {
0357     st1[kk] << "(" << it.first << ", " << it.second << ") ";
0358     ++k;
0359     if (k == n) {
0360       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0361       ++kk;
0362       k = 0;
0363     }
0364   }
0365   if (k > 0)
0366     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0367 }
0368 
0369 void HGCalParameterTester::myPrint(std::string const& s,
0370                                    std::vector<double> const& obj1,
0371                                    std::vector<double> const& obj2,
0372                                    int n) const {
0373   int k(0), kk(0);
0374   std::ostringstream st1[250];
0375   edm::LogVerbatim("HGCalGeom") << s << " with " << obj1.size() << " elements with n " << n << ": 250";
0376   for (unsigned int k1 = 0; k1 < obj1.size(); ++k1) {
0377     st1[kk] << "(" << obj1[k1] << ", " << obj2[k1] << ") ";
0378     ++k;
0379     if (k == n) {
0380       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0381       ++kk;
0382       k = 0;
0383     }
0384   }
0385   if (k > 0)
0386     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0387 }
0388 
0389 void HGCalParameterTester::myPrint(std::string const& s, HGCalParameters::wafer_map const& obj, int n) const {
0390   int k(0), kk(0);
0391   std::ostringstream st1[100];
0392   edm::LogVerbatim("HGCalGeom") << s << " with " << obj.size() << " elements with n " << n << ": 100";
0393   for (auto const& it : obj) {
0394     st1[kk] << it.first << ":" << it.second << ", ";
0395     ++k;
0396     if (k == n) {
0397       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0398       ++kk;
0399       k = 0;
0400     }
0401   }
0402   if (k > 0)
0403     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0404 }
0405 
0406 void HGCalParameterTester::printTrform(HGCalParameters const* phgp) const {
0407   int k(0), kk(0);
0408   std::ostringstream st1[20];
0409   edm::LogVerbatim("HGCalGeom") << "TrformIndex with " << phgp->trformIndex_.size() << " elements with n 7:20";
0410   for (unsigned int i = 0; i < phgp->trformIndex_.size(); ++i) {
0411     std::array<int, 4> id = phgp->getID(i);
0412     st1[kk] << id[0] << ":" << id[1] << ":" << id[2] << ":" << id[3] << ", ";
0413     ++k;
0414     if (k == 7) {
0415       edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0416       ++kk;
0417       k = 0;
0418     }
0419   }
0420   if (k > 0)
0421     edm::LogVerbatim("HGCalGeom") << st1[kk].str();
0422 }
0423 
0424 void HGCalParameterTester::printWaferType(HGCalParameters const* phgp) const {
0425   int k(0);
0426   edm::LogVerbatim("HGCalGeom") << "waferTypes with " << phgp->waferTypes_.size() << " elements";
0427   std::map<std::pair<int, int>, int> kounts;
0428   std::map<std::pair<int, int>, int>::iterator itr;
0429   for (auto const& it : phgp->waferTypes_) {
0430     std::ostringstream st1;
0431     st1 << " [" << k << "] " << HGCalWaferIndex::waferLayer(it.first);
0432     if (HGCalWaferIndex::waferFormat(it.first)) {
0433       st1 << ":" << HGCalWaferIndex::waferU(it.first) << ":" << HGCalWaferIndex::waferV(it.first);
0434     } else {
0435       st1 << ":" << HGCalWaferIndex::waferCopy(it.first);
0436     }
0437     edm::LogVerbatim("HGCalGeom") << st1.str() << " ==> (" << (it.second).first << ":" << (it.second).second << ")";
0438     itr = kounts.find(it.second);
0439     if (itr == kounts.end())
0440       kounts[it.second] = 1;
0441     else
0442       ++(itr->second);
0443     ++k;
0444   }
0445   if (!kounts.empty()) {
0446     edm::LogVerbatim("HGCalGeom") << "Summary of waferTypes ==========================";
0447     for (itr = kounts.begin(); itr != kounts.end(); ++itr)
0448       edm::LogVerbatim("HGCalGeom") << "Type (" << (itr->first).first << ":" << (itr->first).second << ") Kount "
0449                                     << itr->second;
0450   }
0451 }
0452 
0453 DEFINE_FWK_MODULE(HGCalParameterTester);