Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    DQMOffline/ECALMultifitAnalyzer_HI
0004 // Class:      ECALMultifitAnalyzer_HI
0005 //
0006 /**\class ECALMultifitAnalyzer_HI ECALMultifitAnalyzer_HI.cc
0007    DQMOffline/ECALMultifitAnalyzer_HI/plugins/ECALMultifitAnalyzer_HI.cc
0008    Description: [one line class summary]
0009    Implementation:
0010    [Notes on implementation]
0011 */
0012 //
0013 // Original Author:  R. Alex Barbieri
0014 //         Created:  Tue, 17 Nov 2015 17:18:50 GMT
0015 // Modified for DQM: Raghav Kunnawalkam Elayavalli
0016 //                   Wednesday 18 Nov 2015 9:45:AM GVA time
0017 //
0018 //
0019 
0020 // system include files
0021 #include <memory>
0022 #include <vector>
0023 
0024 // user include files
0025 #include "FWCore/Framework/interface/ESHandle.h"
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0028 
0029 #include "FWCore/Framework/interface/Event.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 
0034 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0035 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0036 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0037 
0038 #include "Geometry/Records/interface/CaloTopologyRecord.h"
0039 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
0040 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0041 #include "Geometry/CaloTopology/interface/CaloTopology.h"
0042 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0043 
0044 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0045 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0046 #include "DQMServices/Core/interface/DQMStore.h"
0047 #include "FWCore/ServiceRegistry/interface/Service.h"
0048 
0049 #include "DataFormats/Math/interface/deltaR.h"
0050 
0051 #include "TH2F.h"
0052 
0053 //
0054 // class declaration
0055 //
0056 
0057 class ECALMultifitAnalyzer_HI : public DQMEDAnalyzer {
0058 public:
0059   explicit ECALMultifitAnalyzer_HI(const edm::ParameterSet &);
0060   ~ECALMultifitAnalyzer_HI() override {}
0061 
0062 private:
0063   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0064   void analyze(const edm::Event &, const edm::EventSetup &) override;
0065 
0066   edm::EDGetTokenT<std::vector<reco::Photon>> recoPhotonsCollection_;
0067   edm::EDGetTokenT<reco::CaloJetCollection> caloJetToken_;
0068   edm::EDGetTokenT<EcalRecHitCollection> RecHitCollection_EB_;
0069   edm::EDGetTokenT<EcalRecHitCollection> RecHitCollection_EE_;
0070   edm::ESGetToken<CaloGeometry, CaloGeometryRecord> geomH;
0071 
0072   double mRechitEnergyThreshold;
0073   double mRecoPhotonPtThreshold;
0074   double mRecoJetPtThreshold;
0075   double mDeltaRPhotonThreshold;
0076   double mDeltaRJetThreshold;
0077 
0078   MonitorElement *eb_chi2;
0079   MonitorElement *eb_chi2_eta;
0080   MonitorElement *eb_chi2_e5;
0081   MonitorElement *eb_chi2_e5_eta;
0082   MonitorElement *eb_errors;
0083   MonitorElement *eb_errors_eta;
0084   MonitorElement *eb_errors_e5;
0085   MonitorElement *eb_errors_e5_eta;
0086   MonitorElement *eb_chi2_photon15;
0087   MonitorElement *eb_errors_photon15;
0088   MonitorElement *eb_chi2_jet30;
0089   MonitorElement *eb_errors_jet30;
0090 
0091   MonitorElement *ee_chi2;
0092   MonitorElement *ee_chi2_eta;
0093   MonitorElement *ee_chi2_e5;
0094   MonitorElement *ee_chi2_e5_eta;
0095   MonitorElement *ee_errors;
0096   MonitorElement *ee_errors_eta;
0097   MonitorElement *ee_errors_e5;
0098   MonitorElement *ee_errors_e5_eta;
0099   MonitorElement *ee_chi2_photon15;
0100   MonitorElement *ee_errors_photon15;
0101   MonitorElement *ee_chi2_jet30;
0102   MonitorElement *ee_errors_jet30;
0103 };
0104 
0105 //
0106 // constructors and destructor
0107 //
0108 ECALMultifitAnalyzer_HI::ECALMultifitAnalyzer_HI(const edm::ParameterSet &iConfig)
0109     : recoPhotonsCollection_(consumes<std::vector<reco::Photon>>(iConfig.getParameter<edm::InputTag>("recoPhotonSrc"))),
0110       caloJetToken_(consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("recoJetSrc"))),
0111       RecHitCollection_EB_(consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("RecHitCollection_EB"))),
0112       RecHitCollection_EE_(consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("RecHitCollection_EE"))),
0113       geomH(esConsumes()),
0114       mRechitEnergyThreshold(iConfig.getParameter<double>("rechitEnergyThreshold")),
0115       mRecoPhotonPtThreshold(iConfig.getParameter<double>("recoPhotonPtThreshold")),
0116       mRecoJetPtThreshold(iConfig.getParameter<double>("recoJetPtThreshold")),
0117       mDeltaRPhotonThreshold(iConfig.getParameter<double>("deltaRPhotonThreshold")),
0118       mDeltaRJetThreshold(iConfig.getParameter<double>("deltaRJetThreshold")) {}
0119 
0120 //
0121 // member functions
0122 //
0123 
0124 // ------------ method called for each event  ------------
0125 void ECALMultifitAnalyzer_HI::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0126   using namespace edm;
0127 
0128   const CaloGeometry *geom = &iSetup.getData(geomH);
0129 
0130   Handle<std::vector<reco::Photon>> recoPhotonsHandle;
0131   iEvent.getByToken(recoPhotonsCollection_, recoPhotonsHandle);
0132 
0133   Handle<reco::CaloJetCollection> recoJetHandle;
0134   iEvent.getByToken(caloJetToken_, recoJetHandle);
0135 
0136   Handle<EcalRecHitCollection> ebHandle;
0137   iEvent.getByToken(RecHitCollection_EB_, ebHandle);
0138 
0139   Handle<EcalRecHitCollection> eeHandle;
0140   iEvent.getByToken(RecHitCollection_EE_, eeHandle);
0141 
0142   for (EcalRecHitCollection::const_iterator hit = ebHandle->begin(); hit != ebHandle->end(); ++hit) {
0143     eb_chi2->Fill(hit->chi2());
0144     eb_errors->Fill(hit->energyError());
0145     double eta = geom->getGeometry(hit->detid())->getPosition().eta();
0146     double phi = geom->getGeometry(hit->detid())->getPosition().phi();
0147     eb_chi2_eta->Fill(eta, hit->chi2());
0148     eb_errors_eta->Fill(eta, hit->energyError());
0149     if (hit->energy() > mRechitEnergyThreshold) {
0150       eb_chi2_e5->Fill(hit->chi2());
0151       eb_errors_e5->Fill(hit->energyError());
0152       eb_chi2_e5_eta->Fill(eta, hit->chi2());
0153       eb_errors_e5_eta->Fill(eta, hit->energyError());
0154     }
0155 
0156     for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end();
0157          ++pho) {
0158       if (pho->et() < mRecoPhotonPtThreshold)
0159         continue;
0160       double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
0161       if (dr < mDeltaRPhotonThreshold) {
0162         eb_chi2_photon15->Fill(hit->chi2());
0163         eb_errors_photon15->Fill(hit->energyError());
0164       }
0165     }
0166     for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
0167       if (jet->pt() < mRecoJetPtThreshold)
0168         continue;
0169       double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
0170       if (dr < mDeltaRJetThreshold) {
0171         eb_chi2_jet30->Fill(hit->chi2());
0172         eb_errors_jet30->Fill(hit->energyError());
0173       }
0174     }
0175   }
0176 
0177   for (EcalRecHitCollection::const_iterator hit = eeHandle->begin(); hit != eeHandle->end(); ++hit) {
0178     ee_chi2->Fill(hit->chi2());
0179     ee_errors->Fill(hit->energyError());
0180     double eta = geom->getGeometry(hit->detid())->getPosition().eta();
0181     double phi = geom->getGeometry(hit->detid())->getPosition().phi();
0182     ee_chi2_eta->Fill(eta, hit->chi2());
0183     ee_errors_eta->Fill(eta, hit->energyError());
0184     if (hit->energy() > mRechitEnergyThreshold) {
0185       ee_chi2_e5->Fill(hit->chi2());
0186       ee_errors_e5->Fill(hit->energyError());
0187       ee_chi2_e5_eta->Fill(eta, hit->chi2());
0188       ee_errors_e5_eta->Fill(eta, hit->energyError());
0189     }
0190 
0191     for (std::vector<reco::Photon>::const_iterator pho = recoPhotonsHandle->begin(); pho != recoPhotonsHandle->end();
0192          ++pho) {
0193       if (pho->et() < mRecoPhotonPtThreshold)
0194         continue;
0195       double dr = reco::deltaR(eta, phi, pho->eta(), pho->phi());
0196       if (dr < mDeltaRPhotonThreshold) {
0197         ee_chi2_photon15->Fill(hit->chi2());
0198         ee_errors_photon15->Fill(hit->energyError());
0199       }
0200     }
0201     for (std::vector<reco::CaloJet>::const_iterator jet = recoJetHandle->begin(); jet != recoJetHandle->end(); ++jet) {
0202       if (jet->pt() < mRecoJetPtThreshold)
0203         continue;
0204       double dr = reco::deltaR(eta, phi, jet->eta(), jet->phi());
0205       if (dr < mDeltaRJetThreshold) {
0206         ee_chi2_jet30->Fill(hit->chi2());
0207         ee_errors_jet30->Fill(hit->energyError());
0208       }
0209     }
0210   }
0211 }
0212 
0213 // ------------ method called once each job just before starting event loop
0214 // ------------
0215 void ECALMultifitAnalyzer_HI::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0216   iBooker.setCurrentFolder("EcalCalibration/HIN_MultiFitAnalyzer");
0217 
0218   eb_chi2 = nullptr;
0219   eb_chi2_eta = nullptr;
0220   eb_chi2_e5 = nullptr;
0221   eb_chi2_e5_eta = nullptr;
0222   eb_errors = nullptr;
0223   eb_errors_eta = nullptr;
0224   eb_errors_e5 = nullptr;
0225   eb_errors_e5_eta = nullptr;
0226   eb_chi2_photon15 = nullptr;
0227   eb_errors_photon15 = nullptr;
0228   eb_chi2_jet30 = nullptr;
0229   eb_errors_jet30 = nullptr;
0230 
0231   ee_chi2 = nullptr;
0232   ee_chi2_eta = nullptr;
0233   ee_chi2_e5 = nullptr;
0234   ee_chi2_e5_eta = nullptr;
0235   ee_errors = nullptr;
0236   ee_errors_eta = nullptr;
0237   ee_errors_e5 = nullptr;
0238   ee_errors_e5_eta = nullptr;
0239   ee_chi2_photon15 = nullptr;
0240   ee_errors_photon15 = nullptr;
0241   ee_chi2_jet30 = nullptr;
0242   ee_errors_jet30 = nullptr;
0243 
0244   const int nBins = 500;
0245   const float maxChi2 = 70;
0246   const float maxError = 0.5;
0247 
0248   TH2F *hProfile_Chi2 = new TH2F("hProfile_Chi2", "", nBins, -5, 5, nBins, 0, maxChi2);
0249   TH2F *hProfile_Err = new TH2F("hProfile_Err", "", nBins, -5, 5, nBins, 0, maxError);
0250 
0251   eb_chi2 = iBooker.book1D("rechit_eb_chi2", "Rechit eb_chi2;chi2 fit value;", nBins, 0, maxChi2);
0252   eb_chi2_eta = iBooker.book2D("rechit_eb_eta_Vs_mean_Chi2", hProfile_Chi2);
0253   eb_chi2_e5 = iBooker.book1D(Form("rechit_eb_chi2_e%2.0f", mRechitEnergyThreshold),
0254                               Form("Rechit eb_chi2, e>%2.0fGeV;chi2 fit value;", mRechitEnergyThreshold),
0255                               nBins,
0256                               0,
0257                               maxChi2);
0258   eb_chi2_e5_eta = iBooker.book2D(Form("rechit_eb_chi2_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Chi2);
0259   eb_errors = iBooker.book1D("rechit_eb_errors", "Rechit eb_errors;error on the energy;", nBins, 0, maxError);
0260   eb_errors_eta = iBooker.book2D("rechit_eb_errors_eta", hProfile_Err);
0261   eb_errors_e5 = iBooker.book1D(Form("rechit_eb_errors_e%2.0f", mRechitEnergyThreshold),
0262                                 "Rechit eb_errors, e>5GeV;error on the energy;",
0263                                 nBins,
0264                                 0,
0265                                 maxError);
0266   eb_errors_e5_eta = iBooker.book2D(Form("rechit_eb_errors_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Err);
0267   eb_chi2_photon15 = iBooker.book1D(Form("rechit_eb_chi2_photon%2.0f", mRecoPhotonPtThreshold),
0268                                     "Rechit eb_chi2 near photons;chi2 fit value;",
0269                                     nBins,
0270                                     0,
0271                                     maxChi2);
0272   eb_errors_photon15 = iBooker.book1D(Form("rechit_eb_errors_photon%2.0f", mRecoPhotonPtThreshold),
0273                                       "Rechit eb_errors near photons;error on the energy;",
0274                                       nBins,
0275                                       0,
0276                                       maxError);
0277   eb_chi2_jet30 = iBooker.book1D(Form("rechit_eb_chi2_jet%2.0f", mRecoJetPtThreshold),
0278                                  "Rechit eb_chi2 near jets;chi2 fit value;",
0279                                  nBins,
0280                                  0,
0281                                  maxChi2);
0282   eb_errors_jet30 = iBooker.book1D(Form("rechit_eb_errors_jet%2.0f", mRecoJetPtThreshold),
0283                                    "Rechit eb_errors near jets;error on the energy;",
0284                                    nBins,
0285                                    0,
0286                                    maxError);
0287 
0288   ee_chi2 = iBooker.book1D("rechit_ee_chi2", "Rechit ee_chi2;chi2 fit value;", nBins, 0, maxChi2);
0289   ee_chi2_eta = iBooker.book2D("rechit_ee_chi2_eta", hProfile_Chi2);
0290   ee_chi2_e5 = iBooker.book1D(Form("rechit_ee_chi2_e%2.0f", mRechitEnergyThreshold),
0291                               "Rechit ee_chi2, e>5GeV;chi2 fit value;",
0292                               nBins,
0293                               0,
0294                               maxChi2);
0295   ee_chi2_e5_eta = iBooker.book2D(Form("rechit_ee_chi2_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Chi2);
0296   ee_errors = iBooker.book1D("rechit_ee_errors", "Rechit ee_errors;error on the energy;", nBins, 0, maxError);
0297   ee_errors_eta = iBooker.book2D("rechit_ee_errors_eta", hProfile_Err);
0298   ee_errors_e5 = iBooker.book1D(Form("rechit_ee_errors_e%2.0f", mRechitEnergyThreshold),
0299                                 "Rechit ee_errors, e>5GeV;error on the energy;",
0300                                 nBins,
0301                                 0,
0302                                 maxError);
0303   ee_errors_e5_eta = iBooker.book2D(Form("rechit_ee_errors_e%2.0f_eta", mRechitEnergyThreshold), hProfile_Err);
0304   ee_chi2_photon15 = iBooker.book1D(Form("rechit_ee_chi2_photon%2.0f", mRecoPhotonPtThreshold),
0305                                     "Rechit ee_chi2 near photons;chi2 fit value;",
0306                                     nBins,
0307                                     0,
0308                                     maxChi2);
0309   ee_errors_photon15 = iBooker.book1D(Form("rechit_ee_errors_photon%2.0f", mRecoPhotonPtThreshold),
0310                                       "Rechit ee_errors near photons;error on the energy;",
0311                                       nBins,
0312                                       0,
0313                                       maxError);
0314   ee_chi2_jet30 = iBooker.book1D(Form("rechit_ee_chi2_jet%2.0f", mRecoJetPtThreshold),
0315                                  "Rechit ee_chi2 near jets;chi2 fit value;",
0316                                  nBins,
0317                                  0,
0318                                  maxChi2);
0319   ee_errors_jet30 = iBooker.book1D(Form("rechit_ee_errors_jet%2.0f", mRecoJetPtThreshold),
0320                                    "Rechit ee_errors near jets;error on the energy;",
0321                                    nBins,
0322                                    0,
0323                                    maxError);
0324 
0325   delete hProfile_Chi2;
0326   delete hProfile_Err;
0327 }
0328 
0329 // define this as a plug-in
0330 DEFINE_FWK_MODULE(ECALMultifitAnalyzer_HI);