Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:43

0001 // -*- C++ -*-
0002 //
0003 // Package:    Validation/MtdValidation
0004 // Class:      BtlLocalRecoHarvester
0005 //
0006 /**\class BtlLocalRecoHarvester BtlLocalRecoHarvester.cc Validation/MtdValidation/plugins/BtlLocalRecoHarvester.cc
0007 
0008  Description: BTL SIM hits validation harvester
0009 
0010  Implementation:
0011      [Notes on implementation]
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 BtlLocalRecoHarvester : public DQMEDHarvester {
0027 public:
0028   explicit BtlLocalRecoHarvester(const edm::ParameterSet& iConfig);
0029   ~BtlLocalRecoHarvester() 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   // --- Histograms
0040   MonitorElement* meHitOccupancyLog_;
0041   MonitorElement* meHitOccupancy_;
0042   static constexpr int nRU_ = BTLDetId::kRUPerRod;
0043   MonitorElement* meHitOccupancyRUSlice_[nRU_];
0044 };
0045 
0046 // ------------ constructor and destructor --------------
0047 BtlLocalRecoHarvester::BtlLocalRecoHarvester(const edm::ParameterSet& iConfig)
0048     : folder_(iConfig.getParameter<std::string>("folder")) {}
0049 
0050 BtlLocalRecoHarvester::~BtlLocalRecoHarvester() {}
0051 
0052 // ------------ endjob tasks ----------------------------
0053 void BtlLocalRecoHarvester::dqmEndJob(DQMStore::IBooker& ibook, DQMStore::IGetter& igetter) {
0054   // --- Get the monitoring histograms
0055   MonitorElement* meBtlHitLogEnergy = igetter.get(folder_ + "BtlHitLogEnergy");
0056   MonitorElement* meBtlHitEnergy = igetter.get(folder_ + "BtlHitEnergy");
0057   MonitorElement* meBtlHitEnergyRUSlice[nRU_];
0058   MonitorElement* meNevents = igetter.get(folder_ + "BtlNevents");
0059   bool missing_ru_slice = false;
0060   for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0061     meBtlHitEnergyRUSlice[ihistoRU] = igetter.get(folder_ + "BtlHitEnergyRUSlice" + std::to_string(ihistoRU));
0062     if (!meBtlHitEnergyRUSlice[ihistoRU]) {
0063       missing_ru_slice = true;
0064     }
0065   }
0066 
0067   if (!meBtlHitLogEnergy || !meBtlHitEnergy || missing_ru_slice || !meNevents) {
0068     edm::LogError("BtlLocalRecoHarvester") << "Monitoring histograms not found!" << std::endl;
0069     return;
0070   }
0071 
0072   // --- Get the number of BTL crystals and the number of processed events
0073   const float NBtlCrystals = BTLDetId::kCrystalsBTL;
0074   const float Nevents = meNevents->getEntries();
0075   const float scale_Crystals = (Nevents > 0 ? 1. / (Nevents * NBtlCrystals) : 1.);
0076   const float scale_Crystals_RU = (Nevents > 0 ? 1. / (Nevents * NBtlCrystals / nRU_) : 1.);
0077 
0078   // --- Book the cumulative histograms
0079   ibook.cd(folder_);
0080   meHitOccupancyLog_ =
0081       ibook.book1D("BtlHitOccupancyLog",
0082                    "BTL cell occupancy vs RECO hit energy;log_{10}(E_{RECO} [MeV]); Occupancy per event",
0083                    meBtlHitLogEnergy->getNbinsX(),
0084                    meBtlHitLogEnergy->getTH1()->GetXaxis()->GetXmin(),
0085                    meBtlHitLogEnergy->getTH1()->GetXaxis()->GetXmax());
0086   meHitOccupancy_ = ibook.book1D("BtlHitOccupancy",
0087                                  "BTL cell occupancy vs RECO hit energy;E_{RECO} [MeV]; Occupancy per event",
0088                                  meBtlHitEnergy->getNbinsX(),
0089                                  meBtlHitEnergy->getTH1()->GetXaxis()->GetXmin(),
0090                                  meBtlHitEnergy->getTH1()->GetXaxis()->GetXmax());
0091   for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0092     std::string name_Energy = "BtlHitOccupancyRUSlice" + std::to_string(ihistoRU);
0093     std::string title_Energy = "BTL cell occupancy vs RECO hit energy (RU " + std::to_string(ihistoRU) +
0094                                ");E_{RECO} [MeV]; Occupancy per event";
0095     meHitOccupancyRUSlice_[ihistoRU] = ibook.book1D(name_Energy,
0096                                                     title_Energy,
0097                                                     meBtlHitEnergyRUSlice[ihistoRU]->getNbinsX(),
0098                                                     meBtlHitEnergyRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmin(),
0099                                                     meBtlHitEnergyRUSlice[ihistoRU]->getTH1()->GetXaxis()->GetXmax());
0100   }
0101 
0102   // --- Calculate the cumulative histograms
0103   double bin_sum_log = meBtlHitLogEnergy->getBinContent(meBtlHitLogEnergy->getNbinsX() + 1);
0104   for (int ibin = meBtlHitLogEnergy->getNbinsX(); ibin >= 1; --ibin) {
0105     bin_sum_log += meBtlHitLogEnergy->getBinContent(ibin);
0106     meHitOccupancyLog_->setBinContent(ibin, scale_Crystals * bin_sum_log);
0107   }
0108   double bin_sum = meBtlHitEnergy->getBinContent(meBtlHitEnergy->getNbinsX() + 1);
0109   for (int ibin = meBtlHitEnergy->getNbinsX(); ibin >= 1; --ibin) {
0110     bin_sum += meBtlHitEnergy->getBinContent(ibin);
0111     meHitOccupancy_->setBinContent(ibin, scale_Crystals * bin_sum);
0112   }
0113   for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0114     double bin_sum_RUSlice =
0115         meBtlHitEnergyRUSlice[ihistoRU]->getBinContent(meBtlHitEnergyRUSlice[ihistoRU]->getNbinsX() + 1);
0116     for (int ibin = meBtlHitEnergyRUSlice[ihistoRU]->getNbinsX(); ibin >= 1; --ibin) {
0117       bin_sum_RUSlice += meBtlHitEnergyRUSlice[ihistoRU]->getBinContent(ibin);
0118       meHitOccupancyRUSlice_[ihistoRU]->setBinContent(ibin, scale_Crystals_RU * bin_sum_RUSlice);
0119     }
0120   }
0121 }
0122 
0123 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0124 void BtlLocalRecoHarvester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0125   edm::ParameterSetDescription desc;
0126 
0127   desc.add<std::string>("folder", "MTD/BTL/LocalReco/");
0128 
0129   descriptions.add("btlLocalRecoPostProcessor", desc);
0130 }
0131 
0132 DEFINE_FWK_MODULE(BtlLocalRecoHarvester);