File indexing completed on 2025-06-03 00:12:44
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include <string>
0015
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020
0021 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0022 #include "DQMServices/Core/interface/DQMStore.h"
0023
0024 #include "DataFormats/ForwardDetId/interface/BTLDetId.h"
0025
0026 class BtlSimHitsHarvester : public DQMEDHarvester {
0027 public:
0028 explicit BtlSimHitsHarvester(const edm::ParameterSet& iConfig);
0029 ~BtlSimHitsHarvester() override;
0030
0031 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0032
0033 protected:
0034 void dqmEndJob(DQMStore::IBooker&, DQMStore::IGetter&) override;
0035
0036 private:
0037 const std::string folder_;
0038
0039
0040 MonitorElement* meHitOccupancy_;
0041 MonitorElement* meHitOccupancyCell_;
0042 MonitorElement* meHitOccupancySM_;
0043 static constexpr int nRU_ = BTLDetId::kRUPerRod;
0044 MonitorElement* meHitOccupancyRUSlice_[nRU_];
0045 MonitorElement* meHitOccupancyCellRUSlice_[nRU_];
0046 MonitorElement* meHitOccupancySMRUSlice_[nRU_];
0047 };
0048
0049
0050 BtlSimHitsHarvester::BtlSimHitsHarvester(const edm::ParameterSet& iConfig)
0051 : folder_(iConfig.getParameter<std::string>("folder")) {}
0052
0053 BtlSimHitsHarvester::~BtlSimHitsHarvester() {}
0054
0055
0056 void BtlSimHitsHarvester::dqmEndJob(DQMStore::IBooker& ibook, DQMStore::IGetter& igetter) {
0057
0058 MonitorElement* meBtlHitLogEnergy = igetter.get(folder_ + "BtlHitLogEnergy");
0059 MonitorElement* meBtlHitLogEnergyRUSlice[nRU_];
0060 MonitorElement* meNevents = igetter.get(folder_ + "BtlNevents");
0061 MonitorElement* meBtlHitMultCell = igetter.get(folder_ + "BtlHitMultCell");
0062 MonitorElement* meBtlHitMultCellRUSlice[nRU_];
0063 MonitorElement* meBtlHitMultSM = igetter.get(folder_ + "BtlHitMultSM");
0064 MonitorElement* meBtlHitMultSMRUSlice[nRU_];
0065 bool missing_ru_slice = false;
0066 for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0067 meBtlHitLogEnergyRUSlice[ihistoRU] = igetter.get(folder_ + "BtlHitLogEnergyRUSlice" + std::to_string(ihistoRU));
0068 meBtlHitMultCellRUSlice[ihistoRU] = igetter.get(folder_ + "BtlHitMultCellRUSlice" + std::to_string(ihistoRU));
0069 meBtlHitMultSMRUSlice[ihistoRU] = igetter.get(folder_ + "BtlHitMultSMRUSlice" + std::to_string(ihistoRU));
0070 if (!meBtlHitLogEnergyRUSlice[ihistoRU] || !meBtlHitMultCellRUSlice[ihistoRU] || !meBtlHitMultSMRUSlice[ihistoRU]) {
0071 missing_ru_slice = true;
0072 edm::LogInfo("BtlSimHitsHarvester") << "Monitoring histograms per RUSlice in optional plots!" << std::endl;
0073 }
0074 }
0075 if (!meBtlHitLogEnergy || !meNevents || !meBtlHitMultCell || !meBtlHitMultSM) {
0076 edm::LogError("BtlSimHitsHarvester") << "Monitoring histograms not found!" << std::endl;
0077 return;
0078 }
0079
0080
0081 const float NBtlCrystals = BTLDetId::kCrystalsBTL;
0082 const float NBtlSMs = NBtlCrystals / BTLDetId::kCrystalsPerModuleV2;
0083 const float Nevents = meNevents->getEntries();
0084 const float scale_Crystals = (Nevents > 0 ? 1. / (Nevents * NBtlCrystals) : 1.);
0085 const float scale_Crystals_RU = (Nevents > 0 ? 1. / (Nevents * NBtlCrystals / nRU_) : 1.);
0086 const float scale_SMs = (Nevents > 0 ? 1. / (Nevents * NBtlSMs) : 1.);
0087 const float scale_SMs_RU = (Nevents > 0 ? 1. / (Nevents * NBtlSMs / nRU_) : 1.);
0088
0089
0090 ibook.cd(folder_);
0091 meHitOccupancy_ = ibook.book1D("BtlHitOccupancy",
0092 "BTL cell occupancy vs hit energy;log_{10}(E_{SIM} [MeV]); Occupancy per event",
0093 meBtlHitLogEnergy->getNbinsX(),
0094 meBtlHitLogEnergy->getTH1()->GetXaxis()->GetXmin(),
0095 meBtlHitLogEnergy->getTH1()->GetXaxis()->GetXmax());
0096 meHitOccupancyCell_ =
0097 ibook.book1D("BtlHitOccupancyCell",
0098 "BTL cell occupancy vs energy threshold;log_{10}(E_{th} [MeV]); Occupancy per event",
0099 meBtlHitMultCell->getNbinsX(),
0100 meBtlHitMultCell->getTH1()->GetXaxis()->GetXmin(),
0101 meBtlHitMultCell->getTH1()->GetXaxis()->GetXmax());
0102 meHitOccupancySM_ = ibook.book1D("BtlHitOccupancySM",
0103 "BTL SM occupancy vs energy threshold;log_{10}(E_{th} [MeV]); Occupancy per event",
0104 meBtlHitMultSM->getNbinsX(),
0105 meBtlHitMultSM->getTH1()->GetXaxis()->GetXmin(),
0106 meBtlHitMultSM->getTH1()->GetXaxis()->GetXmax());
0107 if (!missing_ru_slice) {
0108 for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0109 std::string name_LogEnergy = "BtlHitOccupancyRUSlice" + std::to_string(ihistoRU);
0110 std::string title_LogEnergy = "BTL cell occupancy vs hit energy (RU " + std::to_string(ihistoRU) +
0111 ");log_{10}(E_{SIM} [MeV]); Occupancy per event";
0112 meHitOccupancyRUSlice_[ihistoRU] =
0113 ibook.book1D(name_LogEnergy,
0114 title_LogEnergy,
0115 meBtlHitLogEnergyRUSlice[ihistoRU]->getNbinsX(),
0116 meBtlHitLogEnergyRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmin(),
0117 meBtlHitLogEnergyRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmax());
0118 std::string name_cell = "BtlHitOccupancyCellRUSlice" + std::to_string(ihistoRU);
0119 std::string title_cell = "BTL cell occupancy vs energy threshold (RU " + std::to_string(ihistoRU) +
0120 ");log_{10}(E_{th} [MeV]); Occupancy per event";
0121 meHitOccupancyCellRUSlice_[ihistoRU] =
0122 ibook.book1D(name_cell,
0123 title_cell,
0124 meBtlHitMultCellRUSlice[ihistoRU]->getNbinsX(),
0125 meBtlHitMultCellRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmin(),
0126 meBtlHitMultCellRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmax());
0127 std::string name_SM = "BtlHitOccupancySMRUSlice" + std::to_string(ihistoRU);
0128 std::string title_SM = "BTL SM occupancy vs energy threshold (RU " + std::to_string(ihistoRU) +
0129 ");log_{10}(E_{th} [MeV]); Occupancy per event";
0130 meHitOccupancySMRUSlice_[ihistoRU] =
0131 ibook.book1D(name_SM,
0132 title_SM,
0133 meBtlHitMultSMRUSlice[ihistoRU]->getNbinsX(),
0134 meBtlHitMultSMRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmin(),
0135 meBtlHitMultSMRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmax());
0136 }
0137 }
0138
0139
0140 double bin_sum = meBtlHitLogEnergy->getBinContent(meBtlHitLogEnergy->getNbinsX() + 1);
0141 for (int ibin = meBtlHitLogEnergy->getNbinsX(); ibin >= 1; --ibin) {
0142 bin_sum += meBtlHitLogEnergy->getBinContent(ibin);
0143 meHitOccupancy_->setBinContent(ibin, scale_Crystals * bin_sum);
0144 }
0145 if (!missing_ru_slice) {
0146 for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0147 double bin_sum_RUSlice =
0148 meBtlHitLogEnergyRUSlice[ihistoRU]->getBinContent(meBtlHitLogEnergyRUSlice[ihistoRU]->getNbinsX() + 1);
0149 for (int ibin = meBtlHitLogEnergyRUSlice[ihistoRU]->getNbinsX(); ibin >= 1; --ibin) {
0150 bin_sum_RUSlice += meBtlHitLogEnergyRUSlice[ihistoRU]->getBinContent(ibin);
0151 meHitOccupancyRUSlice_[ihistoRU]->setBinContent(ibin, scale_Crystals_RU * bin_sum_RUSlice);
0152 }
0153 }
0154 }
0155
0156
0157 for (int ibin = 0; ibin < meBtlHitMultCell->getNbinsX(); ibin++) {
0158 double bin_content = meBtlHitMultCell->getBinContent(ibin);
0159 meHitOccupancyCell_->setBinContent(ibin, bin_content * scale_Crystals);
0160 }
0161 for (int ibin = 0; ibin < meBtlHitMultSM->getNbinsX(); ibin++) {
0162 double bin_content = meBtlHitMultSM->getBinContent(ibin);
0163 meHitOccupancySM_->setBinContent(ibin, bin_content * scale_SMs);
0164 }
0165 if (!missing_ru_slice) {
0166 for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0167 for (int ibin = 0; ibin < meBtlHitMultCellRUSlice[ihistoRU]->getNbinsX(); ibin++) {
0168 double bin_content_RUSlice = meBtlHitMultCellRUSlice[ihistoRU]->getBinContent(ibin);
0169 meHitOccupancyCellRUSlice_[ihistoRU]->setBinContent(ibin, bin_content_RUSlice * scale_Crystals_RU);
0170 }
0171 for (int ibin = 0; ibin < meBtlHitMultSMRUSlice[ihistoRU]->getNbinsX(); ibin++) {
0172 double bin_content_RUSlice = meBtlHitMultSMRUSlice[ihistoRU]->getBinContent(ibin);
0173 meHitOccupancySMRUSlice_[ihistoRU]->setBinContent(ibin, bin_content_RUSlice * scale_SMs_RU);
0174 }
0175 }
0176 }
0177 }
0178
0179
0180 void BtlSimHitsHarvester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0181 edm::ParameterSetDescription desc;
0182
0183 desc.add<std::string>("folder", "MTD/BTL/SimHits/");
0184
0185 descriptions.add("btlSimHitsPostProcessor", desc);
0186 }
0187
0188 DEFINE_FWK_MODULE(BtlSimHitsHarvester);