Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:25

0001 // -*- C++ -*-//
0002 // Package:    Hcal
0003 // Class:      HcalCollapseAnalyzer
0004 //
0005 /**\class HcalCollapseAnalyzer HcalCollapseAnalyzer.cc
0006  DQMOffline/Hcal/src/HcalCollapseAnalyzer.cc
0007 
0008  Description: Studies the collapsing of HcalRecHits
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Sunanda Banerjee
0015 //         Created:  Thu Dec 26 18:52:02 CST 2017
0016 //
0017 //
0018 
0019 // system include files
0020 #include <string>
0021 #include <vector>
0022 
0023 // Root objects
0024 #include "TH1.h"
0025 
0026 // user include files
0027 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0028 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0029 
0030 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0031 #include "DQMServices/Core/interface/DQMStore.h"
0032 
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/Frameworkfwd.h"
0035 #include "FWCore/Framework/interface/MakerMacros.h"
0036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0038 #include "FWCore/Utilities/interface/EDGetToken.h"
0039 #include "FWCore/Utilities/interface/InputTag.h"
0040 
0041 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0042 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0043 
0044 class HcalCollapseAnalyzer : public DQMEDAnalyzer {
0045 public:
0046   explicit HcalCollapseAnalyzer(const edm::ParameterSet &);
0047   ~HcalCollapseAnalyzer() override {}
0048 
0049   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0050 
0051 private:
0052   void analyze(edm::Event const &, edm::EventSetup const &) override;
0053   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0054 
0055   // ----------member data ---------------------------
0056   const std::string topFolderName_;
0057   const int verbosity_;
0058   const edm::InputTag recHitHBHE_, preRecHitHBHE_;
0059   const bool doHE_, doHB_;
0060   const HcalTopology *theHBHETopology = nullptr;
0061 
0062   edm::EDGetTokenT<HBHERecHitCollection> tok_hbhe_;
0063   edm::EDGetTokenT<HBHERecHitCollection> tok_prehbhe_;
0064   edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> hcalTopologyToken_;
0065 
0066   MonitorElement *h_merge, *h_size, *h_depth;
0067   MonitorElement *h_sfrac, *h_frac, *h_balance;
0068 };
0069 
0070 HcalCollapseAnalyzer::HcalCollapseAnalyzer(const edm::ParameterSet &iConfig)
0071     : topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
0072       verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
0073       recHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"))),
0074       preRecHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"))),
0075       doHE_(iConfig.getUntrackedParameter<bool>("doHE", true)),
0076       doHB_(iConfig.getUntrackedParameter<bool>("doHB", false)),
0077       hcalTopologyToken_{esConsumes<HcalTopology, HcalRecNumberingRecord>()} {
0078   // define tokens for access
0079   tok_hbhe_ = consumes<HBHERecHitCollection>(recHitHBHE_);
0080   tok_prehbhe_ = consumes<HBHERecHitCollection>(preRecHitHBHE_);
0081 
0082   edm::LogVerbatim("Collapse") << "Verbosity " << verbosity_ << " with tags " << recHitHBHE_ << " and "
0083                                << preRecHitHBHE_ << " and Do " << doHB_ << ":" << doHE_;
0084 }
0085 
0086 void HcalCollapseAnalyzer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0087   edm::ParameterSetDescription desc;
0088   desc.add<std::string>("topFolderName", "HcalCollapse");
0089   desc.addUntracked<int>("verbosity", 0);
0090   desc.addUntracked<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"));
0091   desc.addUntracked<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"));
0092   desc.addUntracked<bool>("doHE", true);
0093   desc.addUntracked<bool>("doHB", false);
0094   descriptions.add("hcalCollapseAnalyzer", desc);
0095 }
0096 
0097 void HcalCollapseAnalyzer::analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) {
0098   if (verbosity_ > 0)
0099     edm::LogVerbatim("Collapse") << "Run " << iEvent.id().run() << " Event " << iEvent.id().event()
0100                                  << " starts ==============";
0101 
0102   theHBHETopology = &iSetup.getData(hcalTopologyToken_);
0103 
0104   edm::Handle<HBHERecHitCollection> hbhereco;
0105   iEvent.getByToken(tok_hbhe_, hbhereco);
0106   edm::Handle<HBHERecHitCollection> hbheprereco;
0107   iEvent.getByToken(tok_prehbhe_, hbheprereco);
0108   if (verbosity_ > 0) {
0109     edm::LogVerbatim("Collapse") << "Handle Reco " << hbhereco << " Size " << hbhereco->size();
0110     edm::LogVerbatim("Collapse") << "Handle PreReco " << hbheprereco << " Size " << hbheprereco->size();
0111   }
0112   if (hbhereco.isValid() && hbheprereco.isValid()) {
0113     const HBHERecHitCollection *recohbhe = hbhereco.product();
0114     const HBHERecHitCollection *prerecohbhe = hbheprereco.product();
0115     if (verbosity_ > 0)
0116       edm::LogVerbatim("Collapse") << "Size of hbhereco: " << recohbhe->size()
0117                                    << " and hbheprereco: " << prerecohbhe->size();
0118     double sfrac = (prerecohbhe->empty()) ? 1 : ((double)(recohbhe->size())) / ((double)(prerecohbhe->size()));
0119     h_sfrac->Fill(sfrac);
0120     h_size->Fill(recohbhe->size());
0121     for (const auto &hit : (*recohbhe)) {
0122       HcalDetId id = hit.id();
0123       if (((id.subdet() == HcalEndcap) && doHE_) || ((id.subdet() == HcalBarrel) && doHB_)) {
0124         h_depth->Fill(id.depth());
0125         std::vector<HcalDetId> ids;
0126         theHBHETopology->unmergeDepthDetId(id, ids);
0127         if (verbosity_ > 0) {
0128           edm::LogVerbatim("Collapse") << id << " is derived from " << ids.size() << " components";
0129           for (unsigned int k = 0; k < ids.size(); ++k)
0130             edm::LogVerbatim("Collapse") << "[" << k << "] " << ids[k];
0131         }
0132         h_merge->Fill(ids.size());
0133         double energy = hit.energy();
0134         double etot(0);
0135         unsigned int k(0);
0136         for (const auto phit : (*prerecohbhe)) {
0137           if (std::find(ids.begin(), ids.end(), phit.id()) != ids.end()) {
0138             etot += phit.energy();
0139             double frac = (energy > 0) ? phit.energy() / energy : 1.0;
0140             h_frac->Fill(frac);
0141             if (verbosity_ > 0)
0142               edm::LogVerbatim("Collapse") << "Frac[" << k << "] " << frac << " for " << phit.id();
0143             ++k;
0144           }
0145         }
0146         double frac = (energy > 0) ? etot / energy : 1.0;
0147         h_balance->Fill(frac);
0148         if (verbosity_ > 0)
0149           edm::LogVerbatim("Collapse") << "Total energy " << energy << ":" << etot << ":" << frac;
0150       }
0151     }
0152   }
0153 }
0154 
0155 void HcalCollapseAnalyzer::bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &, edm::EventSetup const &) {
0156   // Book histograms
0157   ibooker.setCurrentFolder(topFolderName_);
0158   h_merge = ibooker.book1D("h_merge", "Number of hits merged", 10, 0.0, 10.0);
0159   h_size = ibooker.book1D("h_size", "Size of the RecHit collection", 100, 500.0, 1500.0);
0160   h_depth = ibooker.book1D("h_depth", "Depth of the Id's used", 10, 0.0, 10.0);
0161   h_sfrac = ibooker.book1D("h_sfrac", "Ratio of sizes of preRecHit and RecHit collections", 150, 0.0, 1.5);
0162   h_frac = ibooker.book1D("h_frac", "Fraction of energy before collapse", 150, 0.0, 1.5);
0163   h_balance = ibooker.book1D("h_balance", "Balance of energy between pre- and post-collapse", 100, 0.5, 1.5);
0164 }
0165 
0166 DEFINE_FWK_MODULE(HcalCollapseAnalyzer);