File indexing completed on 2023-03-17 11:27:49
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
0018
0019 class HGCalHitClient : public DQMEDHarvester {
0020 public:
0021 explicit HGCalHitClient(const edm::ParameterSet&);
0022 ~HGCalHitClient() override = default;
0023
0024 void beginRun(const edm::Run& run, const edm::EventSetup& c) override {}
0025 void dqmEndJob(DQMStore::IBooker& ib, DQMStore::IGetter& ig) override;
0026
0027 private:
0028 const std::string subDirectory_;
0029
0030 int geometryEndjob(const std::vector<MonitorElement*>& hcalMEs);
0031 };
0032
0033 HGCalHitClient::HGCalHitClient(const edm::ParameterSet& iConfig)
0034 : subDirectory_(iConfig.getParameter<std::string>("DirectoryName")) {}
0035
0036 void HGCalHitClient::dqmEndJob(DQMStore::IBooker& ib, DQMStore::IGetter& ig) {
0037 ig.setCurrentFolder("/");
0038 #ifdef EDM_ML_DEBUG
0039 edm::LogVerbatim("HGCalValid") << "HGCalHitValidation :: runClient";
0040 #endif
0041 std::vector<MonitorElement*> hgcalMEs;
0042 std::vector<std::string> fullDirPath = ig.getSubdirs();
0043
0044 for (unsigned int i = 0; i < fullDirPath.size(); i++) {
0045 #ifdef EDM_ML_DEBUG
0046 edm::LogVerbatim("HGCalValid") << "HGCalHitValidation::fullPath: " << fullDirPath.at(i);
0047 #endif
0048 ig.setCurrentFolder(fullDirPath.at(i));
0049 std::vector<std::string> fullSubDirPath = ig.getSubdirs();
0050
0051 for (unsigned int j = 0; j < fullSubDirPath.size(); j++) {
0052 #ifdef EDM_ML_DEBUG
0053 edm::LogVerbatim("HGCalValid") << "HGCalHitValidation:: fullSubPath: " << fullSubDirPath.at(j);
0054 #endif
0055 if (strcmp(fullSubDirPath.at(j).c_str(), subDirectory_.c_str()) == 0) {
0056 hgcalMEs = ig.getContents(fullSubDirPath.at(j));
0057 #ifdef EDM_ML_DEBUG
0058 edm::LogVerbatim("HGCalValid") << "HGCalHitValidation:: hgcalMES size : " << hgcalMEs.size();
0059 #endif
0060 if (!geometryEndjob(hgcalMEs))
0061 edm::LogWarning("HGCalValid") << "\nError in HGcalHitEndjob!";
0062 }
0063 }
0064 }
0065 }
0066
0067 int HGCalHitClient::geometryEndjob(const std::vector<MonitorElement*>& hgcalMEs) {
0068 std::string dets[3] = {"hee", "hef", "heb"};
0069 std::string hist1[2] = {"EnRec", "EnSim"};
0070 std::string hist2[7] = {"dzVsZ", "dyVsY", "dxVsX", "RecVsSimZ", "RecVsSimY", "RecVsSimX", "EnSimRec"};
0071
0072 std::vector<MonitorElement*> hist1_;
0073 std::vector<MonitorElement*> hist2_;
0074
0075
0076 for (unsigned int idet = 0; idet < 3; ++idet) {
0077 char name[100];
0078 for (unsigned int kh = 0; kh < 2; ++kh) {
0079 sprintf(name, "%s%s", dets[idet].c_str(), hist1[kh].c_str());
0080 for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
0081 if (strcmp(hgcalMEs[ih]->getName().c_str(), name) == 0) {
0082 hist1_.push_back(hgcalMEs[ih]);
0083 double nevent = hist1_.back()->getEntries();
0084 int nbinsx = hist1_.back()->getNbinsX();
0085 for (int i = 1; i <= nbinsx; ++i) {
0086 double binValue = hist1_.back()->getBinContent(i) / nevent;
0087 hist1_.back()->setBinContent(i, binValue);
0088 }
0089 }
0090 }
0091 }
0092 for (unsigned int kh = 0; kh < 7; ++kh) {
0093 sprintf(name, "%s%s", dets[idet].c_str(), hist2[kh].c_str());
0094 for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
0095 if (strcmp(hgcalMEs[ih]->getName().c_str(), name) == 0) {
0096 hist2_.push_back(hgcalMEs[ih]);
0097 double nevent = hist2_.back()->getEntries();
0098 int nbinsx = hist2_.back()->getNbinsX();
0099 int nbinsy = hist2_.back()->getNbinsY();
0100 for (int i = 1; i <= nbinsx; ++i) {
0101 for (int j = 1; j <= nbinsy; ++j) {
0102 double binValue = hist2_.back()->getBinContent(i, j) / nevent;
0103 hist2_.back()->setBinContent(i, j, binValue);
0104 }
0105 }
0106 }
0107 }
0108 }
0109 }
0110
0111 return 1;
0112 }
0113
0114 #include "FWCore/Framework/interface/MakerMacros.h"
0115
0116 DEFINE_FWK_MODULE(HGCalHitClient);