CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTBCalo/HcalPlotter/src/HcalQLPlotAnal.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HcalQLPlotAnal
00004 // Class:      HcalQLPlotAnal
00005 // 
00013 //
00014 // Original Author:  Phillip R. Dudero
00015 //         Created:  Tue Jan 16 21:11:37 CST 2007
00016 // $Id: HcalQLPlotAnal.cc,v 1.5 2007/12/31 18:43:18 ratnik Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDAnalyzer.h"
00027 #include "FWCore/Framework/interface/Event.h"
00028 #include "DataFormats/Common/interface/Handle.h"
00029 #include "FWCore/Framework/interface/ESHandle.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00033 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00034 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00035 #include "TBDataFormats/HcalTBObjects/interface/HcalTBTriggerData.h"
00036 #include "RecoTBCalo/HcalPlotter/src/HcalQLPlotAnalAlgos.h"
00037 #include <string>
00038 //
00039 // class declaration
00040 //
00041 
00042 class HcalQLPlotAnal : public edm::EDAnalyzer {
00043    public:
00044       explicit HcalQLPlotAnal(const edm::ParameterSet&);
00045       ~HcalQLPlotAnal();
00046 
00047 
00048    private:
00049       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00050       virtual void endJob() ;
00051 
00052       // ----------member data ---------------------------
00053   edm::InputTag hbheRHLabel_,hoRHLabel_,hfRHLabel_;
00054   edm::InputTag hcalDigiLabel_, hcalTrigLabel_;
00055   bool doCalib_;
00056   double calibFC2GeV_;
00057   HcalQLPlotAnalAlgos * algo_;
00058 
00059 };
00060 
00061 //
00062 // constants, enums and typedefs
00063 //
00064 
00065 //
00066 // static data member definitions
00067 //
00068 
00069 //
00070 // constructors and destructor
00071 //
00072 HcalQLPlotAnal::HcalQLPlotAnal(const edm::ParameterSet& iConfig) :
00073   hbheRHLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hbheRHtag")),
00074   hoRHLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hoRHtag")),
00075   hfRHLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hfRHtag")),
00076   hcalDigiLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hcalDigiTag")),
00077   hcalTrigLabel_(iConfig.getUntrackedParameter<edm::InputTag>("hcalTrigTag")),
00078   doCalib_(iConfig.getUntrackedParameter<bool>("doCalib",false)),
00079   calibFC2GeV_(iConfig.getUntrackedParameter<double>("calibFC2GeV",0.2))
00080 {
00081   algo_ = new
00082     HcalQLPlotAnalAlgos(iConfig.getUntrackedParameter<std::string>("outputFilename").c_str(),
00083                         iConfig.getParameter<edm::ParameterSet>("HistoParameters"));
00084 }
00085 
00086 
00087 HcalQLPlotAnal::~HcalQLPlotAnal()
00088 {
00089    // do anything here that needs to be done at destruction time
00090    // (e.g. close files, deallocate resources etc.)
00091 }
00092 
00093 
00094 //
00095 // member functions
00096 //
00097 
00098 // ------------ method called to for each event  ------------
00099 void
00100 HcalQLPlotAnal::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00101 {
00102   // Step A/C: Get Inputs and process (repeatedly)
00103   edm::Handle<HcalTBTriggerData> trig;
00104   iEvent.getByLabel(hcalTrigLabel_,trig);
00105   if (!trig.isValid()) {
00106     edm::LogError("HcalQLPlotAnal::analyze") << "No Trigger Data found, skip event";
00107     return;
00108   } else {
00109     algo_->SetEventType(*trig);
00110   }
00111   edm::Handle<HBHEDigiCollection> hbhedg;
00112   iEvent.getByLabel(hcalDigiLabel_,hbhedg);
00113   if (!hbhedg.isValid()) {
00114     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HBHE Digis/RecHits not found";
00115   } else {
00116     algo_->processDigi(*hbhedg);
00117   }
00118   edm::Handle<HBHERecHitCollection> hbherh;  
00119   iEvent.getByLabel(hbheRHLabel_,hbherh);
00120   if (!hbherh.isValid()) {
00121     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HBHE Digis/RecHits not found";
00122   } else {
00123     algo_->processRH(*hbherh,*hbhedg);
00124   }
00125 
00126   edm::Handle<HODigiCollection> hodg; 
00127   iEvent.getByLabel(hcalDigiLabel_,hodg);
00128   if (!hodg.isValid()) {
00129     // can't find it!
00130     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HO Digis/RecHits not found";
00131   } else {
00132     algo_->processDigi(*hodg);
00133   }
00134   edm::Handle<HORecHitCollection> horh;
00135   iEvent.getByLabel(hoRHLabel_,horh);
00136   if (!horh.isValid()) {
00137     // can't find it!
00138     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HO Digis/RecHits not found";
00139   } else {
00140     algo_->processRH(*horh,*hodg);
00141   }
00142   
00143   edm::Handle<HFDigiCollection> hfdg;
00144   iEvent.getByLabel(hcalDigiLabel_,hfdg);
00145 
00146   if (!hfdg.isValid()) {
00147     // can't find it!
00148     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HF Digis/RecHits not found";
00149   } else {
00150     algo_->processDigi(*hfdg);
00151   }
00152 
00153   edm::Handle<HFRecHitCollection> hfrh;
00154   iEvent.getByLabel(hfRHLabel_,hfrh);
00155   if (!hfrh.isValid()) {
00156     // can't find it!
00157     edm::LogWarning("HcalQLPlotAnal::analyze") << "One of HF Digis/RecHits not found";
00158   } else {
00159     algo_->processRH(*hfrh,*hfdg);
00160   }
00161 
00162   if (doCalib_) {
00163     // No rechits as of yet...
00164     edm::Handle<HcalCalibDigiCollection> calibdg;
00165     iEvent.getByLabel(hcalDigiLabel_,calibdg);
00166     if (!calibdg.isValid()) {
00167       edm::LogWarning("HcalQLPlotAnal::analyze") << "Hcal Calib Digis not found";
00168     } else {
00169       algo_->processDigi(*calibdg,calibFC2GeV_);
00170     }
00171   }
00172 
00173 }
00174 
00175 
00176 // ------------ method called once each job just after ending the event loop  ------------
00177 void 
00178 HcalQLPlotAnal::endJob()
00179 {
00180   algo_->end();
00181 }
00182 
00183 //define this as a plug-in
00184 DEFINE_FWK_MODULE(HcalQLPlotAnal);