CMS 3D CMS Logo

PatZjetsElectronAnalyzer.cc
Go to the documentation of this file.
1 #include <map>
2 #include <string>
3 
4 #include "TH1.h"
5 
12 
14 
16 public:
18  ~PatZjetsElectronAnalyzer() override;
19 
20 private:
21  void beginJob() override;
22  void analyze(const edm::Event&, const edm::EventSetup&) override;
23  void endJob() override;
24 
25  // simple map to contain all histograms;
26  // histograms are booked in the beginJob()
27  // method
28  std::map<std::string, TH1F*> histContainer_;
29 
30  // input tags
32 };
33 
35  : histContainer_(),
36  srcToken_(consumes<edm::View<pat::Electron> >(iConfig.getUntrackedParameter<edm::InputTag>("src"))) {}
37 
39 
41  // get electron collection
43  iEvent.getByToken(srcToken_, elecs);
44 
45  // loop electrons
46  for (edm::View<pat::Electron>::const_iterator elec = elecs->begin(); elec != elecs->end(); ++elec) {
47  // fill simple histograms
48  histContainer_["pt"]->Fill(elec->pt());
49  histContainer_["eta"]->Fill(elec->eta());
50  histContainer_["phi"]->Fill(elec->phi());
51  histContainer_["iso"]->Fill((elec->trackIso() + elec->caloIso()) / elec->pt());
52  histContainer_["eop"]->Fill(elec->eSeedClusterOverP());
53  histContainer_["clus"]->Fill(elec->e1x5() / elec->e5x5());
54  // fill enegry flow histogram for isolation
55  for (int bin = 1; bin <= histContainer_["dr"]->GetNbinsX(); ++bin) {
56  double lowerEdge = histContainer_["dr"]->GetBinLowEdge(bin);
57  double upperEdge = histContainer_["dr"]->GetBinLowEdge(bin) + histContainer_["dr"]->GetBinWidth(bin);
58  histContainer_["dr"]->Fill(
59  histContainer_["dr"]->GetBinCenter(bin),
60  elec->trackIsoDeposit()->depositWithin(upperEdge) - elec->trackIsoDeposit()->depositWithin(lowerEdge));
61  }
62  // fill electron id histograms
63  if (elec->electronID("eidRobustLoose") > 0.5)
64  histContainer_["eIDs"]->Fill(0);
65  if (elec->electronID("eidRobustTight") > 0.5)
66  histContainer_["eIDs"]->Fill(1);
67  if (elec->electronID("eidLoose") > 0.5)
68  histContainer_["eIDs"]->Fill(2);
69  if (elec->electronID("eidTight") > 0.5)
70  histContainer_["eIDs"]->Fill(3);
71  if (elec->electronID("eidRobustHighEnergy") > 0.5)
72  histContainer_["eIDs"]->Fill(4);
73  }
74 }
75 
77  // register to the TFileService
79 
80  // book histograms:
81  histContainer_["pt"] = fs->make<TH1F>("pt", "pt", 150, 0., 150.);
82  histContainer_["eta"] = fs->make<TH1F>("eta", "eta", 50, 0., 5.);
83  histContainer_["phi"] = fs->make<TH1F>("phi", "phi", 60, 3.14, 3.14);
84  histContainer_["iso"] = fs->make<TH1F>("iso", "iso", 30, 0., 10.);
85  histContainer_["dr"] = fs->make<TH1F>("dr", "dr", 40, 0., 1.);
86  histContainer_["eop"] = fs->make<TH1F>("eop", "eop", 40, 0., 1.);
87  histContainer_["clus"] = fs->make<TH1F>("clus", "clus", 40, 0., 1.);
88  histContainer_["eIDs"] = fs->make<TH1F>("eIDs", "eIDS", 5, 0., 5.);
89 }
90 
92 
sistrip::View
View
Definition: ConstantsForView.h:26
PatZjetsElectronAnalyzer
Definition: PatZjetsElectronAnalyzer.cc:15
edm::EDGetTokenT
Definition: EDGetToken.h:33
Electron
Definition: Electron.py:1
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
singleTopDQM_cfi.upperEdge
upperEdge
Definition: singleTopDQM_cfi.py:72
EDAnalyzer.h
PatZjetsElectronAnalyzer::endJob
void endJob() override
Definition: PatZjetsElectronAnalyzer.cc:91
edm::Handle
Definition: AssociativeIterator.h:50
edm::EDAnalyzer
Definition: EDAnalyzer.h:28
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Service.h
singleTopDQM_cfi.elecs
elecs
Definition: singleTopDQM_cfi.py:41
TFileService.h
singleTopDQM_cfi.lowerEdge
lowerEdge
Definition: singleTopDQM_cfi.py:71
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
edm::Service< TFileService >
iEvent
int iEvent
Definition: GenABIO.cc:224
PatZjetsElectronAnalyzer::srcToken_
edm::EDGetTokenT< edm::View< pat::Electron > > srcToken_
Definition: PatZjetsElectronAnalyzer.cc:31
edm::EventSetup
Definition: EventSetup.h:58
pat
Definition: HeavyIon.h:7
PatZjetsElectronAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: PatZjetsElectronAnalyzer.cc:40
TFileService::make
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
InputTag.h
newFWLiteAna.bin
bin
Definition: newFWLiteAna.py:161
PatZjetsElectronAnalyzer::~PatZjetsElectronAnalyzer
~PatZjetsElectronAnalyzer() override
Definition: PatZjetsElectronAnalyzer.cc:38
PatZjetsElectronAnalyzer::histContainer_
std::map< std::string, TH1F * > histContainer_
Definition: PatZjetsElectronAnalyzer.cc:28
PatZjetsElectronAnalyzer::beginJob
void beginJob() override
Definition: PatZjetsElectronAnalyzer.cc:76
PatZjetsElectronAnalyzer::PatZjetsElectronAnalyzer
PatZjetsElectronAnalyzer(const edm::ParameterSet &)
Definition: PatZjetsElectronAnalyzer.cc:34
Electron.h
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
ParameterSet.h
edm::Event
Definition: Event.h:73