Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // system include files
0002 #include <memory>
0003 
0004 // user include files
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0007 
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0014 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0015 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0016 
0017 #include "EgammaAnalysis/ElectronTools/interface/EcalIsolationCorrector.h"
0018 
0019 #include "TH1F.h"
0020 
0021 #include <iostream>
0022 
0023 class CorrectECALIsolation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0024 public:
0025   explicit CorrectECALIsolation(const edm::ParameterSet&);
0026   ~CorrectECALIsolation() override;
0027 
0028 private:
0029   void analyze(const edm::Event&, const edm::EventSetup&) override;
0030 
0031   edm::ParameterSet conf_;
0032 
0033   bool isData_;
0034   edm::EDGetTokenT<reco::GsfElectronCollection> tokenGsfElectrons_;
0035 
0036   TH1F* uncorrectedIsolationEB_;
0037   TH1F* correctedIsolationEB_;
0038   TH1F* uncorrectedIsolationEE_;
0039   TH1F* correctedIsolationEE_;
0040 };
0041 
0042 CorrectECALIsolation::CorrectECALIsolation(const edm::ParameterSet& iConfig) : conf_(iConfig) {
0043   usesResource(TFileService::kSharedResource);
0044 
0045   isData_ = iConfig.getUntrackedParameter<bool>("isData", false);
0046   tokenGsfElectrons_ = consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("Electrons"));
0047 
0048   edm::Service<TFileService> fs;
0049   uncorrectedIsolationEB_ = fs->make<TH1F>("uncorrectedIsolationEB", "uncorrected IsolationEB", 50, 0, 10);
0050   correctedIsolationEB_ = fs->make<TH1F>("correctedIsolationEB", "corrected IsolationEB", 50, 0, 10);
0051   uncorrectedIsolationEE_ = fs->make<TH1F>("uncorrectedIsolationEE", "uncorrected IsolationEE", 50, 0, 10);
0052   correctedIsolationEE_ = fs->make<TH1F>("correctedIsolationEE", "corrected IsolationEE", 50, 0, 10);
0053 }
0054 
0055 CorrectECALIsolation::~CorrectECALIsolation() {}
0056 
0057 void CorrectECALIsolation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0058   edm::Handle<reco::GsfElectronCollection> theEGammaCollection;
0059   iEvent.getByToken(tokenGsfElectrons_, theEGammaCollection);
0060   const reco::GsfElectronCollection theEGamma = *(theEGammaCollection.product());
0061 
0062   // Setup a corrector for electrons
0063   EcalIsolationCorrector ecalIsoCorr(true);
0064 
0065   unsigned nele = theEGammaCollection->size();
0066 
0067   for (unsigned iele = 0; iele < nele; ++iele) {
0068     reco::GsfElectronRef myElectronRef(theEGammaCollection, iele);
0069 
0070     float uncorrIso = myElectronRef->dr03EcalRecHitSumEt();
0071     float corrIso = ecalIsoCorr.correctForHLTDefinition(*myElectronRef, iEvent.id().run(), isData_);
0072     std::cout << "Uncorrected Isolation Sum: " << uncorrIso << " - Corrected: " << corrIso << std::endl;
0073 
0074     if (myElectronRef->isEB()) {
0075       uncorrectedIsolationEB_->Fill(uncorrIso);
0076       correctedIsolationEB_->Fill(corrIso);
0077     } else {
0078       uncorrectedIsolationEE_->Fill(uncorrIso);
0079       correctedIsolationEE_->Fill(corrIso);
0080     }
0081   }
0082 }
0083 
0084 //define this as a plug-in
0085 DEFINE_FWK_MODULE(CorrectECALIsolation);