Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:01

0001 #include <map>
0002 #include <string>
0003 
0004 #include "TH1.h"
0005 
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0008 #include "FWCore/Utilities/interface/InputTag.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ServiceRegistry/interface/Service.h"
0011 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0012 
0013 #include "DataFormats/PatCandidates/interface/Electron.h"
0014 
0015 class PatZjetsElectronAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0016 public:
0017   explicit PatZjetsElectronAnalyzer(const edm::ParameterSet&);
0018   ~PatZjetsElectronAnalyzer() override;
0019 
0020 private:
0021   void beginJob() override;
0022   void analyze(const edm::Event&, const edm::EventSetup&) override;
0023   void endJob() override;
0024 
0025   // simple map to contain all histograms;
0026   // histograms are booked in the beginJob()
0027   // method
0028   std::map<std::string, TH1F*> histContainer_;
0029 
0030   // input tags
0031   edm::EDGetTokenT<edm::View<pat::Electron> > srcToken_;
0032 };
0033 
0034 PatZjetsElectronAnalyzer::PatZjetsElectronAnalyzer(const edm::ParameterSet& iConfig)
0035     : histContainer_(),
0036       srcToken_(consumes<edm::View<pat::Electron> >(iConfig.getUntrackedParameter<edm::InputTag>("src"))) {
0037   usesResource(TFileService::kSharedResource);
0038 }
0039 
0040 PatZjetsElectronAnalyzer::~PatZjetsElectronAnalyzer() {}
0041 
0042 void PatZjetsElectronAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0043   // get electron collection
0044   edm::Handle<edm::View<pat::Electron> > elecs;
0045   iEvent.getByToken(srcToken_, elecs);
0046 
0047   // loop electrons
0048   for (edm::View<pat::Electron>::const_iterator elec = elecs->begin(); elec != elecs->end(); ++elec) {
0049     // fill simple histograms
0050     histContainer_["pt"]->Fill(elec->pt());
0051     histContainer_["eta"]->Fill(elec->eta());
0052     histContainer_["phi"]->Fill(elec->phi());
0053     histContainer_["iso"]->Fill((elec->trackIso() + elec->caloIso()) / elec->pt());
0054     histContainer_["eop"]->Fill(elec->eSeedClusterOverP());
0055     histContainer_["clus"]->Fill(elec->e1x5() / elec->e5x5());
0056     // fill enegry flow histogram for isolation
0057     for (int bin = 1; bin <= histContainer_["dr"]->GetNbinsX(); ++bin) {
0058       double lowerEdge = histContainer_["dr"]->GetBinLowEdge(bin);
0059       double upperEdge = histContainer_["dr"]->GetBinLowEdge(bin) + histContainer_["dr"]->GetBinWidth(bin);
0060       histContainer_["dr"]->Fill(
0061           histContainer_["dr"]->GetBinCenter(bin),
0062           elec->trackIsoDeposit()->depositWithin(upperEdge) - elec->trackIsoDeposit()->depositWithin(lowerEdge));
0063     }
0064     // fill electron id histograms
0065     if (elec->electronID("eidRobustLoose") > 0.5)
0066       histContainer_["eIDs"]->Fill(0);
0067     if (elec->electronID("eidRobustTight") > 0.5)
0068       histContainer_["eIDs"]->Fill(1);
0069     if (elec->electronID("eidLoose") > 0.5)
0070       histContainer_["eIDs"]->Fill(2);
0071     if (elec->electronID("eidTight") > 0.5)
0072       histContainer_["eIDs"]->Fill(3);
0073     if (elec->electronID("eidRobustHighEnergy") > 0.5)
0074       histContainer_["eIDs"]->Fill(4);
0075   }
0076 }
0077 
0078 void PatZjetsElectronAnalyzer::beginJob() {
0079   // register to the TFileService
0080   edm::Service<TFileService> fs;
0081 
0082   // book histograms:
0083   histContainer_["pt"] = fs->make<TH1F>("pt", "pt", 150, 0., 150.);
0084   histContainer_["eta"] = fs->make<TH1F>("eta", "eta", 50, 0., 5.);
0085   histContainer_["phi"] = fs->make<TH1F>("phi", "phi", 60, 3.14, 3.14);
0086   histContainer_["iso"] = fs->make<TH1F>("iso", "iso", 30, 0., 10.);
0087   histContainer_["dr"] = fs->make<TH1F>("dr", "dr", 40, 0., 1.);
0088   histContainer_["eop"] = fs->make<TH1F>("eop", "eop", 40, 0., 1.);
0089   histContainer_["clus"] = fs->make<TH1F>("clus", "clus", 40, 0., 1.);
0090   histContainer_["eIDs"] = fs->make<TH1F>("eIDs", "eIDS", 5, 0., 5.);
0091 }
0092 
0093 void PatZjetsElectronAnalyzer::endJob() {}
0094 
0095 #include "FWCore/Framework/interface/MakerMacros.h"
0096 DEFINE_FWK_MODULE(PatZjetsElectronAnalyzer);