Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:31

0001 // system include files
0002 #include <string>
0003 #include <iostream>
0004 
0005 // user include files
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/Utilities/interface/InputTag.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 
0015 #include "DataFormats/DetId/interface/DetId.h"
0016 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0017 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0018 #include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
0019 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0020 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0021 
0022 #include "Geometry/HcalCommonData/interface/HcalHitRelabeller.h"
0023 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0024 
0025 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
0026 #include "SimDataFormats/CaloTest/interface/HGCalTestNumbering.h"
0027 #include "SimDataFormats/CaloHit/interface/PCaloHit.h"
0028 
0029 class HcalDumpHits : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0030 public:
0031   explicit HcalDumpHits(const edm::ParameterSet&);
0032   ~HcalDumpHits() override = default;
0033 
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036 private:
0037   void beginJob() override {}
0038   void endJob() override {}
0039   void beginRun(edm::Run const&, edm::EventSetup const&) override;
0040   void endRun(edm::Run const&, edm::EventSetup const&) override {}
0041   void analyze(edm::Event const&, edm::EventSetup const&) override;
0042 
0043   // ----------member data ---------------------------
0044   const edm::EDGetTokenT<edm::PCaloHitContainer> simHitSource_;
0045   const edm::EDGetTokenT<QIE11DigiCollection> hbheDigiSource_;
0046   const edm::EDGetTokenT<QIE10DigiCollection> hfDigiSource_;
0047   const edm::EDGetTokenT<HODigiCollection> hoDigiSource_;
0048   const edm::EDGetTokenT<HBHERecHitCollection> hbheRecHitSource_;
0049   const edm::EDGetTokenT<HFRecHitCollection> hfRecHitSource_;
0050   const edm::EDGetTokenT<HORecHitCollection> hoRecHitSource_;
0051   const edm::ESGetToken<HcalDDDRecConstants, HcalRecNumberingRecord> tok_HRNDC_;
0052   const HcalDDDRecConstants* hcons_;
0053 };
0054 
0055 HcalDumpHits::HcalDumpHits(const edm::ParameterSet& iConfig)
0056     : simHitSource_(consumes<edm::PCaloHitContainer>(iConfig.getParameter<edm::InputTag>("simHitSource"))),
0057       hbheDigiSource_(consumes<QIE11DigiCollection>(iConfig.getParameter<edm::InputTag>("hbheDigiSource"))),
0058       hfDigiSource_(consumes<QIE10DigiCollection>(iConfig.getParameter<edm::InputTag>("hfDigiSource"))),
0059       hoDigiSource_(consumes<HODigiCollection>(iConfig.getParameter<edm::InputTag>("hoDigiSource"))),
0060       hbheRecHitSource_(consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheRecHitSource"))),
0061       hfRecHitSource_(consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfR(ecHitSource"))),
0062       hoRecHitSource_(consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoRecHitSource"))),
0063       tok_HRNDC_(esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord, edm::Transition::BeginRun>()),
0064       hcons_(nullptr) {}
0065 
0066 void HcalDumpHits::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0067   edm::ParameterSetDescription desc;
0068   desc.add<edm::InputTag>("simHitSource", edm::InputTag("g4SimHits", "HcalHits"));
0069   desc.add<edm::InputTag>("hbheDigiSource", edm::InputTag("hcalDigis"));
0070   desc.add<edm::InputTag>("hfDigiSource", edm::InputTag("hcalDigis"));
0071   desc.add<edm::InputTag>("hoDigiSource", edm::InputTag("hcalDigis"));
0072   desc.add<edm::InputTag>("hbheRecHitSource", edm::InputTag("hbhereco"));
0073   desc.add<edm::InputTag>("hfRecHitSource", edm::InputTag("hfreco"));
0074   desc.add<edm::InputTag>("hoRecHitSource", edm::InputTag("horeco"));
0075   descriptions.add("hcalDumpHits", desc);
0076 }
0077 
0078 void HcalDumpHits::beginRun(const edm::Run&, const edm::EventSetup& iSetup) { hcons_ = &iSetup.getData(tok_HRNDC_); }
0079 
0080 void HcalDumpHits::analyze(const edm::Event& iEvent, const edm::EventSetup&) {
0081   // first SimHits
0082   const edm::Handle<edm::PCaloHitContainer>& theCaloHitContainer = iEvent.getHandle(simHitSource_);
0083   if (theCaloHitContainer.isValid()) {
0084     edm::LogVerbatim("HcalValidation") << theCaloHitContainer->size() << " SimHits in HCAL";
0085     unsigned int k(0);
0086     for (auto const& hit : *(theCaloHitContainer.product())) {
0087       unsigned int id = hit.id();
0088       HcalDetId hid = HcalHitRelabeller::relabel(id, hcons_);
0089       edm::LogVerbatim("HcalValidation") << "[" << k << "] " << hid << " E " << hit.energy() << " T " << hit.time();
0090       ++k;
0091     }
0092   }
0093 
0094   // Digis (HBHE/HF/HO)
0095   const edm::Handle<QIE11DigiCollection>& hbheDigiCollection = iEvent.getHandle(hbheDigiSource_);
0096   if (hbheDigiCollection.isValid()) {
0097     edm::LogVerbatim("HcalValidation") << hbheDigiCollection->size() << " Digis for HB/HE";
0098     unsigned int k(0);
0099     for (auto const& it : *(hbheDigiCollection.product())) {
0100       QIE11DataFrame hit(it);
0101       std::ostringstream ost;
0102       ost << "[" << k << "] " << HcalDetId(hit.detid()) << " with " << hit.size() << " words:";
0103       unsigned int k1(0);
0104       for (auto itr = hit.begin(); itr != hit.end(); ++itr, ++k1)
0105         ost << " [" << k1 << "] " << (*itr);
0106       edm::LogVerbatim("HcalValidation") << ost.str();
0107       ++k;
0108     }
0109   }
0110   const edm::Handle<QIE10DigiCollection>& hfDigiCollection = iEvent.getHandle(hfDigiSource_);
0111   if (hfDigiCollection.isValid()) {
0112     edm::LogVerbatim("HcalValidation") << hfDigiCollection->size() << " Digis for HF";
0113     unsigned int k(0);
0114     for (auto const& it : *(hfDigiCollection.product())) {
0115       QIE10DataFrame hit(it);
0116       std::ostringstream ost;
0117       ost << "[" << k << "] " << HcalDetId(hit.detid()) << " with " << hit.size() << " words ";
0118       unsigned int k1(0);
0119       for (auto itr = hit.begin(); itr != hit.end(); ++itr, ++k1)
0120         ost << " [" << k1 << "] " << (*itr);
0121       edm::LogVerbatim("HcalValidation") << ost.str();
0122       ++k;
0123     }
0124   }
0125   const edm::Handle<HODigiCollection>& hoDigiCollection = iEvent.getHandle(hoDigiSource_);
0126   if (hoDigiCollection.isValid()) {
0127     edm::LogVerbatim("HcalValidation") << hoDigiCollection->size() << " Digis for HO";
0128     unsigned int k(0);
0129     for (auto const& it : *(hoDigiCollection.product())) {
0130       HODataFrame hit(it);
0131       std::ostringstream ost;
0132       ost << "[" << k << "] " << HcalDetId(hit.id()) << " with " << hit.size() << " words ";
0133       for (int k1 = 0; k1 < hit.size(); ++k1)
0134         ost << " [" << k1 << "] " << hit.sample(k1);
0135       edm::LogVerbatim("HcalValidation") << ost.str();
0136       ++k;
0137     }
0138   }
0139 
0140   // RecHits (HBHE/HF/HO)
0141   const edm::Handle<HBHERecHitCollection>& hbhecoll = iEvent.getHandle(hbheRecHitSource_);
0142   if (hbhecoll.isValid()) {
0143     edm::LogVerbatim("HcalValidation") << hbhecoll->size() << " RecHits for HB/HE";
0144     unsigned int k(0);
0145     for (const auto& it : *(hbhecoll.product())) {
0146       HBHERecHit hit(it);
0147       edm::LogVerbatim("HcalValidation") << "[" << k << "] = " << hit;
0148       ++k;
0149     }
0150   }
0151   const edm::Handle<HFRecHitCollection>& hfcoll = iEvent.getHandle(hfRecHitSource_);
0152   if (hfcoll.isValid()) {
0153     edm::LogVerbatim("HcalValidation") << hfcoll->size() << " RecHits for HF";
0154     unsigned int k(0);
0155     for (const auto& it : *(hfcoll.product())) {
0156       HFRecHit hit(it);
0157       edm::LogVerbatim("HcalValidation") << "[" << k << "] = " << hit;
0158       ++k;
0159     }
0160   }
0161   const edm::Handle<HORecHitCollection>& hocoll = iEvent.getHandle(hoRecHitSource_);
0162   if (hocoll.isValid()) {
0163     edm::LogVerbatim("HcalValidation") << hocoll->size() << " RecHits for HO";
0164     unsigned int k(0);
0165     for (const auto& it : *(hocoll.product())) {
0166       HORecHit hit(it);
0167       edm::LogVerbatim("HcalValidation") << "[" << k << "] = " << hit;
0168       ++k;
0169     }
0170   }
0171 }
0172 
0173 //define this as a plug-in
0174 #include "FWCore/Framework/interface/MakerMacros.h"
0175 #include "FWCore/PluginManager/interface/ModuleDef.h"
0176 
0177 DEFINE_FWK_MODULE(HcalDumpHits);