Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <iostream>
0002 #include <sstream>
0003 #include <string>
0004 #include <vector>
0005 
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Framework/interface/ESTransientHandle.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 
0016 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0017 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0018 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0019 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0020 
0021 class HcalTopologyTester : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0022 public:
0023   explicit HcalTopologyTester(const edm::ParameterSet&);
0024 
0025   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0026 
0027 private:
0028   void analyze(edm::Event const&, edm::EventSetup const&) override;
0029   void beginJob() override {}
0030   void beginRun(edm::Run const&, edm::EventSetup const&) override {}
0031   void endRun(edm::Run const&, edm::EventSetup const&) override {}
0032   void doTest(const HcalTopology& topology);
0033 
0034   // ----------member data ---------------------------
0035   const edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tokTopo_;
0036 };
0037 
0038 HcalTopologyTester::HcalTopologyTester(const edm::ParameterSet&)
0039     : tokTopo_{esConsumes<HcalTopology, HcalRecNumberingRecord>(edm::ESInputTag{})} {}
0040 
0041 void HcalTopologyTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0042   edm::ParameterSetDescription desc;
0043   desc.setUnknown();
0044   descriptions.add("hcalTopologyTester", desc);
0045 }
0046 
0047 void HcalTopologyTester::analyze(edm::Event const&, edm::EventSetup const& iSetup) { doTest(iSetup.getData(tokTopo_)); }
0048 
0049 void HcalTopologyTester::doTest(const HcalTopology& topology) {
0050   // First test on movements along eta/phi directions
0051   edm::LogVerbatim("HCalGeom") << "\nTest on movements along eta/phi directions"
0052                                << "\n==========================================";
0053   for (int idet = 0; idet < 4; idet++) {
0054     HcalSubdetector subdet = HcalBarrel;
0055     if (idet == 1)
0056       subdet = HcalOuter;
0057     else if (idet == 2)
0058       subdet = HcalEndcap;
0059     else if (idet == 3)
0060       subdet = HcalForward;
0061     for (int depth = 1; depth < 4; ++depth) {
0062       for (int ieta = -41; ieta <= 41; ieta++) {
0063         for (int iphi = 1; iphi <= 72; iphi++) {
0064           const HcalDetId id(subdet, ieta, iphi, depth);
0065           if (topology.valid(id)) {
0066             std::vector<DetId> idE = topology.east(id);
0067             std::vector<DetId> idW = topology.west(id);
0068             std::vector<DetId> idN = topology.north(id);
0069             std::vector<DetId> idS = topology.south(id);
0070             std::vector<DetId> idU = topology.up(id);
0071             edm::LogVerbatim("HCalGeom") << "Neighbours for : Tower " << id;
0072             std::ostringstream st1;
0073             st1 << "          " << idE.size() << " sets along East:";
0074             for (auto& i : idE)
0075               st1 << " " << (HcalDetId)(i());
0076             edm::LogVerbatim("HCalGeom") << st1.str();
0077             std::ostringstream st2;
0078             st2 << "          " << idW.size() << " sets along West:";
0079             for (auto& i : idW)
0080               st2 << " " << (HcalDetId)(i());
0081             edm::LogVerbatim("HCalGeom") << st2.str();
0082             std::ostringstream st3;
0083             st3 << "          " << idN.size() << " sets along North:";
0084             for (auto& i : idN)
0085               st3 << " " << (HcalDetId)(i());
0086             edm::LogVerbatim("HCalGeom") << st3.str();
0087             std::ostringstream st4;
0088             st4 << "          " << idS.size() << " sets along South:";
0089             for (auto& i : idS)
0090               st4 << " " << (HcalDetId)(i());
0091             edm::LogVerbatim("HCalGeom") << st4.str();
0092             std::ostringstream st5;
0093             st5 << "          " << idU.size() << " sets up in depth:";
0094             for (auto& i : idU)
0095               st5 << " " << (HcalDetId)(i());
0096             edm::LogVerbatim("HCalGeom") << st5.str();
0097           }
0098         }
0099       }
0100     }
0101   }
0102 
0103   // Check on Dense Index
0104   edm::LogVerbatim("HCalGeom") << "\nCheck on Dense Index"
0105                                << "\n=====================";
0106   int maxDepthHB = topology.maxDepthHB();
0107   int maxDepthHE = topology.maxDepthHE();
0108   for (int det = 1; det <= HcalForward; det++) {
0109     for (int eta = -HcalDetId::kHcalEtaMask2; eta <= (int)(HcalDetId::kHcalEtaMask2); eta++) {
0110       for (unsigned int phi = 0; phi <= HcalDetId::kHcalPhiMask2; phi++) {
0111         for (int depth = 1; depth < maxDepthHB + maxDepthHE; depth++) {
0112           HcalDetId cell((HcalSubdetector)det, eta, phi, depth);
0113           if (topology.valid(cell)) {
0114             unsigned int dense = topology.detId2denseId(DetId(cell));
0115             DetId id = topology.denseId2detId(dense);
0116             std::string cherr = (cell != HcalDetId(id)) ? " **** ERROR *****" : "";
0117             edm::LogVerbatim("HCalGeom") << cell << " Dense " << std::hex << dense << std::dec << " o/p "
0118                                          << HcalDetId(id) << cherr;
0119           }
0120         }
0121       }
0122     }
0123   }
0124 
0125   // Check list of depths
0126   edm::LogVerbatim("HCalGeom") << "\nCheck list of Depths"
0127                                << "\n====================";
0128   for (int eta = topology.lastHERing() - 2; eta <= topology.lastHERing(); ++eta) {
0129     for (unsigned int phi = 0; phi <= HcalDetId::kHcalPhiMask2; phi++) {
0130       for (int depth = 1; depth <= maxDepthHE; depth++) {
0131         HcalDetId cell(HcalEndcap, eta, phi, depth);
0132         if (topology.valid(cell)) {
0133           std::vector<int> depths = topology.mergedDepthList29(cell);
0134           std::ostringstream st1;
0135           st1 << cell << " is with merge depth flag " << topology.mergedDepth29(cell) << " having " << depths.size()
0136               << " merged depths:";
0137           for (unsigned int k = 0; k < depths.size(); ++k)
0138             st1 << " [" << k << "]:" << depths[k];
0139           edm::LogVerbatim("HCalGeom") << st1.str();
0140         }
0141       }
0142     }
0143   }
0144 }
0145 
0146 //define this as a plug-in
0147 DEFINE_FWK_MODULE(HcalTopologyTester);