File indexing completed on 2024-04-06 12:32:35
0001 #include <memory>
0002 #include <unistd.h>
0003 #include <iostream>
0004 #include <fstream>
0005 #include <vector>
0006
0007 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Framework/interface/Run.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "DQMServices/Core/interface/DQMStore.h"
0015 #include "DataFormats/Common/interface/Handle.h"
0016
0017 class HGCalGeometryClient : public DQMEDHarvester {
0018 public:
0019 explicit HGCalGeometryClient(const edm::ParameterSet&);
0020 ~HGCalGeometryClient() override = default;
0021
0022 void beginRun(const edm::Run& run, const edm::EventSetup& c) override {}
0023 void dqmEndJob(DQMStore::IBooker& ib, DQMStore::IGetter& ig) override;
0024
0025 private:
0026 int geometryEndjob(const std::vector<MonitorElement*>& hcalMEs);
0027
0028 const std::string subDirectory_;
0029 };
0030
0031 HGCalGeometryClient::HGCalGeometryClient(const edm::ParameterSet& iConfig)
0032 : subDirectory_(iConfig.getParameter<std::string>("DirectoryName")) {}
0033
0034 void HGCalGeometryClient::dqmEndJob(DQMStore::IBooker& ib, DQMStore::IGetter& ig) {
0035 ig.setCurrentFolder("/");
0036 #ifdef EDM_ML_DEBUG
0037 edm::LogVerbatim("HGCalValid") << "HGCalGeometry :: runClient";
0038 #endif
0039 std::vector<MonitorElement*> hgcalMEs;
0040 std::vector<std::string> fullDirPath = ig.getSubdirs();
0041
0042 for (unsigned int i = 0; i < fullDirPath.size(); i++) {
0043 #ifdef EDM_ML_DEBUG
0044 edm::LogVerbatim("HGCalValid") << "HGCalGeometry::fullPath: " << fullDirPath.at(i);
0045 #endif
0046 ig.setCurrentFolder(fullDirPath.at(i));
0047 std::vector<std::string> fullSubDirPath = ig.getSubdirs();
0048
0049 for (unsigned int j = 0; j < fullSubDirPath.size(); j++) {
0050 #ifdef EDM_ML_DEBUG
0051 edm::LogVerbatim("HGCalValid") << "HGCalGeometry:: fullSubPath: " << fullSubDirPath.at(j);
0052 #endif
0053 if (strcmp(fullSubDirPath.at(j).c_str(), subDirectory_.c_str()) == 0) {
0054 hgcalMEs = ig.getContents(fullSubDirPath.at(j));
0055 #ifdef EDM_ML_DEBUG
0056 edm::LogVerbatim("HGCalValid") << "HGCalGeometry:: hgcalMES size : " << hgcalMEs.size();
0057 #endif
0058 if (!geometryEndjob(hgcalMEs))
0059 edm::LogWarning("HGCalValid") << "\nError in GeometryEndjob!";
0060 }
0061 }
0062 }
0063 }
0064
0065 int HGCalGeometryClient::geometryEndjob(const std::vector<MonitorElement*>& hgcalMEs) {
0066 std::string dets[3] = {"hee", "hef", "heb"};
0067 std::string hist1[4] = {"TotEdepStep", "dX", "dY", "dZ"};
0068 std::string hist2[10] = {"LayerVsEnStep",
0069 "XG4VsId",
0070 "YG4VsId",
0071 "ZG4VsId",
0072 "dxVsX",
0073 "dyVsY",
0074 "dzVsZ",
0075 "dxVsLayer",
0076 "dyVsLayer",
0077 "dzVsLayer"};
0078 std::vector<MonitorElement*> hist1_;
0079 std::vector<MonitorElement*> hist2_;
0080
0081
0082 for (unsigned int idet = 0; idet < 3; ++idet) {
0083 char name[100];
0084 for (unsigned int kh = 0; kh < 4; ++kh) {
0085 sprintf(name, "%s%s", dets[idet].c_str(), hist1[kh].c_str());
0086 for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
0087 if (strcmp(hgcalMEs[ih]->getName().c_str(), name) == 0) {
0088 hist1_.push_back(hgcalMEs[ih]);
0089 double nevent = hist1_.back()->getEntries();
0090 int nbinsx = hist1_.back()->getNbinsX();
0091 for (int i = 1; i <= nbinsx; ++i) {
0092 double binValue = hist1_.back()->getBinContent(i) / nevent;
0093 hist1_.back()->setBinContent(i, binValue);
0094 }
0095 }
0096 }
0097 }
0098 for (unsigned int kh = 0; kh < 10; ++kh) {
0099 sprintf(name, "%s%s", dets[idet].c_str(), hist2[kh].c_str());
0100 for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
0101 if (strcmp(hgcalMEs[ih]->getName().c_str(), name) == 0) {
0102 hist2_.push_back(hgcalMEs[ih]);
0103 double nevent = hist2_.back()->getEntries();
0104 int nbinsx = hist2_.back()->getNbinsX();
0105 int nbinsy = hist2_.back()->getNbinsY();
0106 for (int i = 1; i <= nbinsx; ++i) {
0107 for (int j = 1; j <= nbinsy; ++j) {
0108 double binValue = hist2_.back()->getBinContent(i, j) / nevent;
0109 hist2_.back()->setBinContent(i, j, binValue);
0110 }
0111 }
0112 }
0113 }
0114 }
0115 }
0116
0117 return 1;
0118 }
0119
0120 #include "FWCore/Framework/interface/MakerMacros.h"
0121
0122 DEFINE_FWK_MODULE(HGCalGeometryClient);