Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:22:03

0001 // -*- C++ -*-
0002 //
0003 // Package:    HcalQLPlotAnal
0004 // Class:      HcalQLPlotAnal
0005 //
0006 /**\class HcalQLPlotAnal HcalQLPlotAnal.cc RecoTBCalo/HcalQLPlotAnal/src/HcalQLPlotAnal.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Phillip R. Dudero
0015 //         Created:  Tue Jan 16 21:11:37 CST 2007
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "DataFormats/Common/interface/Handle.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0030 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0031 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0032 #include "TBDataFormats/HcalTBObjects/interface/HcalTBTriggerData.h"
0033 #include "RecoTBCalo/HcalPlotter/src/HcalQLPlotAnalAlgos.h"
0034 #include <string>
0035 //
0036 // class declaration
0037 //
0038 
0039 class HcalQLPlotAnal : public edm::one::EDAnalyzer<> {
0040 public:
0041   explicit HcalQLPlotAnal(const edm::ParameterSet&);
0042   ~HcalQLPlotAnal() override = default;
0043 
0044 private:
0045   void analyze(const edm::Event&, const edm::EventSetup&) override;
0046   void endJob() override;
0047 
0048   // ----------member data ---------------------------
0049   edm::InputTag hcalDigiLabel_;
0050   edm::EDGetTokenT<HBHERecHitCollection> tok_hbherec_;
0051   edm::EDGetTokenT<HORecHitCollection> tok_horec_;
0052   edm::EDGetTokenT<HFRecHitCollection> tok_hfrec_;
0053   edm::EDGetTokenT<HBHEDigiCollection> tok_hbhe_;
0054   edm::EDGetTokenT<HODigiCollection> tok_ho_;
0055   edm::EDGetTokenT<HFDigiCollection> tok_hf_;
0056   edm::EDGetTokenT<HcalCalibDigiCollection> tok_calib_;
0057   edm::EDGetTokenT<HcalTBTriggerData> tok_tb_;
0058   bool doCalib_;
0059   double calibFC2GeV_;
0060   HcalQLPlotAnalAlgos* algo_;
0061 };
0062 
0063 //
0064 // constants, enums and typedefs
0065 //
0066 
0067 //
0068 // static data member definitions
0069 //
0070 
0071 //
0072 // constructors and destructor
0073 //
0074 HcalQLPlotAnal::HcalQLPlotAnal(const edm::ParameterSet& iConfig)
0075     : hcalDigiLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hcalDigiTag")),
0076       doCalib_(iConfig.getUntrackedParameter<bool>("doCalib", false)),
0077       calibFC2GeV_(iConfig.getUntrackedParameter<double>("calibFC2GeV", 0.2)) {
0078   algo_ = new HcalQLPlotAnalAlgos(iConfig.getUntrackedParameter<std::string>("outputFilename").c_str(),
0079                                   iConfig.getParameter<edm::ParameterSet>("HistoParameters"));
0080 
0081   tok_hbherec_ = consumes<HBHERecHitCollection>(iConfig.getUntrackedParameter<edm::InputTag>("hbheRHtag"));
0082   tok_horec_ = consumes<HORecHitCollection>(iConfig.getUntrackedParameter<edm::InputTag>("hoRHtag"));
0083   tok_hfrec_ = consumes<HFRecHitCollection>(iConfig.getUntrackedParameter<edm::InputTag>("hfRHtag"));
0084   tok_hbhe_ = consumes<HBHEDigiCollection>(hcalDigiLabel_);
0085   tok_ho_ = consumes<HODigiCollection>(hcalDigiLabel_);
0086   tok_hf_ = consumes<HFDigiCollection>(hcalDigiLabel_);
0087   tok_calib_ = consumes<HcalCalibDigiCollection>(hcalDigiLabel_);
0088   tok_tb_ = consumes<HcalTBTriggerData>(iConfig.getUntrackedParameter<edm::InputTag>("hcalTrigTag"));
0089 }
0090 
0091 //
0092 // member functions
0093 //
0094 
0095 // ------------ method called to for each event  ------------
0096 void HcalQLPlotAnal::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0097   // Step A/C: Get Inputs and process (repeatedly)
0098   edm::Handle<HcalTBTriggerData> trig;
0099   iEvent.getByToken(tok_tb_, trig);
0100   if (!trig.isValid()) {
0101     edm::LogError("HcalQLPlotAnal::analyze") << "No Trigger Data found, skip event";
0102     return;
0103   } else {
0104     algo_->SetEventType(*trig);
0105   }
0106   edm::Handle<HBHEDigiCollection> hbhedg;
0107   iEvent.getByToken(tok_hbhe_, hbhedg);
0108   if (!hbhedg.isValid()) {
0109     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HBHE Digis/RecHits not found";
0110   } else {
0111     algo_->processDigi(*hbhedg);
0112   }
0113   edm::Handle<HBHERecHitCollection> hbherh;
0114   iEvent.getByToken(tok_hbherec_, hbherh);
0115   if (!hbherh.isValid()) {
0116     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HBHE Digis/RecHits not found";
0117   } else {
0118     algo_->processRH(*hbherh, *hbhedg);
0119   }
0120 
0121   edm::Handle<HODigiCollection> hodg;
0122   iEvent.getByToken(tok_ho_, hodg);
0123   if (!hodg.isValid()) {
0124     // can't find it!
0125     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HO Digis/RecHits not found";
0126   } else {
0127     algo_->processDigi(*hodg);
0128   }
0129   edm::Handle<HORecHitCollection> horh;
0130   iEvent.getByToken(tok_horec_, horh);
0131   if (!horh.isValid()) {
0132     // can't find it!
0133     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HO Digis/RecHits not found";
0134   } else {
0135     algo_->processRH(*horh, *hodg);
0136   }
0137 
0138   edm::Handle<HFDigiCollection> hfdg;
0139   iEvent.getByToken(tok_hf_, hfdg);
0140 
0141   if (!hfdg.isValid()) {
0142     // can't find it!
0143     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HF Digis/RecHits not found";
0144   } else {
0145     algo_->processDigi(*hfdg);
0146   }
0147 
0148   edm::Handle<HFRecHitCollection> hfrh;
0149   iEvent.getByToken(tok_hfrec_, hfrh);
0150   if (!hfrh.isValid()) {
0151     // can't find it!
0152     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HF Digis/RecHits not found";
0153   } else {
0154     algo_->processRH(*hfrh, *hfdg);
0155   }
0156 
0157   if (doCalib_) {
0158     // No rechits as of yet...
0159     edm::Handle<HcalCalibDigiCollection> calibdg;
0160     iEvent.getByToken(tok_calib_, calibdg);
0161     if (!calibdg.isValid()) {
0162       edm::LogWarning("HcalQLPlotAnal::analyze") << "Hcal Calib Digis not found";
0163     } else {
0164       algo_->processDigi(*calibdg, calibFC2GeV_);
0165     }
0166   }
0167 }
0168 
0169 // ------------ method called once each job just after ending the event loop  ------------
0170 void HcalQLPlotAnal::endJob() { algo_->end(); }
0171 
0172 //define this as a plug-in
0173 DEFINE_FWK_MODULE(HcalQLPlotAnal);