Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:57:32

0001 /*
0002  * \file DQMHcalIterativePhiSymAlCaReco.cc
0003  *
0004  * \author Sunanda Banerjee
0005  *
0006  *
0007  *
0008  * Description: Monitoring of Iterative Phi Symmetry Calibration Stream
0009  */
0010 
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 
0021 // DQM include files
0022 
0023 #include "DQMServices/Core/interface/DQMStore.h"
0024 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0025 
0026 // work on collections
0027 
0028 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0029 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0030 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0031 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0032 
0033 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
0034 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0035 
0036 #include <string>
0037 
0038 class DQMHcalIterativePhiSymAlCaReco : public DQMOneEDAnalyzer<> {
0039 public:
0040   DQMHcalIterativePhiSymAlCaReco(const edm::ParameterSet &);
0041   ~DQMHcalIterativePhiSymAlCaReco() override;
0042 
0043   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0044 
0045 protected:
0046   //  void beginRun(const edm::Run& r, const edm::EventSetup& c);
0047   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0048   void analyze(const edm::Event &e, const edm::EventSetup &c) override;
0049 
0050   void dqmEndRun(const edm::Run &r, const edm::EventSetup &c) override {}
0051 
0052 private:
0053   static constexpr int maxDepth_ = 7;
0054   //
0055   // Monitor elements
0056   //
0057   MonitorElement *hiDistr2D_[maxDepth_];
0058 
0059   MonitorElement *hiDistrHBHEsize1D_;
0060   MonitorElement *hiDistrHFsize1D_;
0061   MonitorElement *hiDistrHOsize1D_;
0062 
0063   std::string folderName_;
0064   int hiDistr_y_nbin_;
0065   int hiDistr_x_nbin_;
0066   double hiDistr_y_min_;
0067   double hiDistr_y_max_;
0068   double hiDistr_x_min_;
0069   double hiDistr_x_max_;
0070   /// object to monitor
0071 
0072   edm::EDGetTokenT<HBHERecHitCollection> tok_hbhe_;
0073   edm::EDGetTokenT<HORecHitCollection> tok_ho_;
0074   edm::EDGetTokenT<HFRecHitCollection> tok_hf_;
0075 };
0076 
0077 // ******************************************
0078 // constructors
0079 // *****************************************
0080 
0081 DQMHcalIterativePhiSymAlCaReco::DQMHcalIterativePhiSymAlCaReco(const edm::ParameterSet &ps) {
0082   //
0083   // Input from configurator file
0084   //
0085   folderName_ = ps.getParameter<std::string>("folderName");
0086 
0087   // histogram parameters
0088   tok_ho_ = consumes<HORecHitCollection>(ps.getParameter<edm::InputTag>("hoInput"));
0089   tok_hf_ = consumes<HFRecHitCollection>(ps.getParameter<edm::InputTag>("hfInput"));
0090   tok_hbhe_ = consumes<HBHERecHitCollection>(ps.getParameter<edm::InputTag>("hbheInput"));
0091 
0092   // Distribution of rechits in iPhi, iEta
0093   hiDistr_y_nbin_ = ps.getUntrackedParameter<int>("hiDistr_y_nbin", 72);
0094   hiDistr_y_min_ = ps.getUntrackedParameter<double>("hiDistr_y_min", 0.5);
0095   hiDistr_y_max_ = ps.getUntrackedParameter<double>("hiDistr_y_max", 72.5);
0096   hiDistr_x_nbin_ = ps.getUntrackedParameter<int>("hiDistr_x_nbin", 83);
0097   hiDistr_x_min_ = ps.getUntrackedParameter<double>("hiDistr_x_min", -41.5);
0098   hiDistr_x_max_ = ps.getUntrackedParameter<double>("hiDistr_x_max", 41.5);
0099 }
0100 
0101 DQMHcalIterativePhiSymAlCaReco::~DQMHcalIterativePhiSymAlCaReco() {}
0102 
0103 void DQMHcalIterativePhiSymAlCaReco::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0104   edm::ParameterSetDescription desc;
0105   desc.add<std::string>("folderName", "ALCAStreamHcalIterativePhiSym");
0106   desc.add<edm::InputTag>("hbheInput", edm::InputTag("hbhereco"));
0107   desc.add<edm::InputTag>("hfInput", edm::InputTag("hfreco"));
0108   desc.add<edm::InputTag>("hoInput", edm::InputTag("horeco"));
0109   desc.addUntracked<int>("hiDistr_y_nbin", 72);
0110   desc.addUntracked<double>("hiDistr_y_min", 0.5);
0111   desc.addUntracked<double>("hiDistr_y_max", 72.5);
0112   desc.addUntracked<int>("hiDistr_x_nbin", 83);
0113   desc.addUntracked<double>("hiDistr_x_min", -41.5);
0114   desc.addUntracked<double>("hiDistr_x_max", 41.5);
0115   descriptions.add("dqmHcalIterativePhiSymAlCaReco", desc);
0116 }
0117 
0118 //--------------------------------------------------------
0119 void DQMHcalIterativePhiSymAlCaReco::bookHistograms(DQMStore::IBooker &ibooker,
0120                                                     edm::Run const &irun,
0121                                                     edm::EventSetup const &isetup) {
0122   // create and cd into new folder
0123   ibooker.setCurrentFolder(folderName_);
0124 
0125   // book some histograms 1D
0126   hiDistrHBHEsize1D_ = ibooker.book1D("DistrHBHEsize", "Size of HBHE Collection", 100, 0.0, 10000.0);
0127   hiDistrHFsize1D_ = ibooker.book1D("DistrHFsize", "Size of HF Collection", 100, 0.0, 10000.0);
0128   hiDistrHOsize1D_ = ibooker.book1D("DistrHOsize", "Size of HO Collection", 100, 0.0, 3000.0);
0129 
0130   // Eta-phi occupancy
0131   for (int k = 0; k < maxDepth_; ++k) {
0132     char name[20], title[20];
0133     sprintf(name, "MBdepth%d", (k + 1));
0134     sprintf(title, "Depth %d", (k + 1));
0135     hiDistr2D_[k] = ibooker.book2D(
0136         name, title, hiDistr_x_nbin_, hiDistr_x_min_, hiDistr_x_max_, hiDistr_y_nbin_, hiDistr_y_min_, hiDistr_y_max_);
0137     hiDistr2D_[k]->setAxisTitle("i#phi ", 2);
0138     hiDistr2D_[k]->setAxisTitle("i#eta ", 1);
0139   }
0140 }
0141 
0142 //--------------------------------------------------------
0143 
0144 //-------------------------------------------------------------
0145 
0146 void DQMHcalIterativePhiSymAlCaReco::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0147   // First HBHE RecHits
0148   edm::Handle<HBHERecHitCollection> hbheMB;
0149   iEvent.getByToken(tok_hbhe_, hbheMB);
0150   if (!hbheMB.isValid()) {
0151     edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HBHE RecHit Collection!";
0152   } else {
0153     const HBHERecHitCollection Hithbhe = *(hbheMB.product());
0154     hiDistrHBHEsize1D_->Fill(Hithbhe.size());
0155     for (HBHERecHitCollection::const_iterator hbheItr = Hithbhe.begin(); hbheItr != Hithbhe.end(); hbheItr++) {
0156       HcalDetId hid = HcalDetId(hbheItr->detid());
0157       int id = hid.depth() - 1;
0158       if ((id >= 0) && (id < maxDepth_))
0159         hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hbheItr->energy());
0160     }
0161   }
0162 
0163   // Then HF RecHits
0164   edm::Handle<HFRecHitCollection> hfMB;
0165   iEvent.getByToken(tok_hf_, hfMB);
0166   if (!hfMB.isValid()) {
0167     edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HF RecHit Collection!";
0168   } else {
0169     const HFRecHitCollection Hithf = *(hfMB.product());
0170     hiDistrHFsize1D_->Fill(Hithf.size());
0171     for (HFRecHitCollection::const_iterator hfItr = Hithf.begin(); hfItr != Hithf.end(); hfItr++) {
0172       HcalDetId hid = HcalDetId(hfItr->detid());
0173       int id = hid.depth() - 1;
0174       if ((id >= 0) && (id < maxDepth_))
0175         hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hfItr->energy());
0176     }
0177   }
0178 
0179   // And finally HO RecHits
0180   edm::Handle<HORecHitCollection> hoMB;
0181   iEvent.getByToken(tok_ho_, hoMB);
0182   if (!hoMB.isValid()) {
0183     edm::LogWarning("DQMHcal") << "DQMHcalIterativePhiSymAlCaReco: Error! can't get HO RecHit Collection!";
0184   } else {
0185     const HORecHitCollection Hitho = *(hoMB.product());
0186     hiDistrHOsize1D_->Fill(Hitho.size());
0187     for (HORecHitCollection::const_iterator hoItr = Hitho.begin(); hoItr != Hitho.end(); hoItr++) {
0188       HcalDetId hid = HcalDetId(hoItr->detid());
0189       int id = hid.depth() - 1;
0190       if ((id >= 0) && (id < maxDepth_))
0191         hiDistr2D_[id]->Fill(hid.ieta(), hid.iphi(), hoItr->energy());
0192     }
0193   }
0194 
0195 }  // analyze
0196 
0197 DEFINE_FWK_MODULE(DQMHcalIterativePhiSymAlCaReco);