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 
15 class PatZjetsElectronAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
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  usesResource(TFileService::kSharedResource);
38 }
39 
41 
43  // get electron collection
45  iEvent.getByToken(srcToken_, elecs);
46 
47  // loop electrons
48  for (edm::View<pat::Electron>::const_iterator elec = elecs->begin(); elec != elecs->end(); ++elec) {
49  // fill simple histograms
50  histContainer_["pt"]->Fill(elec->pt());
51  histContainer_["eta"]->Fill(elec->eta());
52  histContainer_["phi"]->Fill(elec->phi());
53  histContainer_["iso"]->Fill((elec->trackIso() + elec->caloIso()) / elec->pt());
54  histContainer_["eop"]->Fill(elec->eSeedClusterOverP());
55  histContainer_["clus"]->Fill(elec->e1x5() / elec->e5x5());
56  // fill enegry flow histogram for isolation
57  for (int bin = 1; bin <= histContainer_["dr"]->GetNbinsX(); ++bin) {
58  double lowerEdge = histContainer_["dr"]->GetBinLowEdge(bin);
59  double upperEdge = histContainer_["dr"]->GetBinLowEdge(bin) + histContainer_["dr"]->GetBinWidth(bin);
60  histContainer_["dr"]->Fill(
61  histContainer_["dr"]->GetBinCenter(bin),
62  elec->trackIsoDeposit()->depositWithin(upperEdge) - elec->trackIsoDeposit()->depositWithin(lowerEdge));
63  }
64  // fill electron id histograms
65  if (elec->electronID("eidRobustLoose") > 0.5)
66  histContainer_["eIDs"]->Fill(0);
67  if (elec->electronID("eidRobustTight") > 0.5)
68  histContainer_["eIDs"]->Fill(1);
69  if (elec->electronID("eidLoose") > 0.5)
70  histContainer_["eIDs"]->Fill(2);
71  if (elec->electronID("eidTight") > 0.5)
72  histContainer_["eIDs"]->Fill(3);
73  if (elec->electronID("eidRobustHighEnergy") > 0.5)
74  histContainer_["eIDs"]->Fill(4);
75  }
76 }
77 
79  // register to the TFileService
81 
82  // book histograms:
83  histContainer_["pt"] = fs->make<TH1F>("pt", "pt", 150, 0., 150.);
84  histContainer_["eta"] = fs->make<TH1F>("eta", "eta", 50, 0., 5.);
85  histContainer_["phi"] = fs->make<TH1F>("phi", "phi", 60, 3.14, 3.14);
86  histContainer_["iso"] = fs->make<TH1F>("iso", "iso", 30, 0., 10.);
87  histContainer_["dr"] = fs->make<TH1F>("dr", "dr", 40, 0., 1.);
88  histContainer_["eop"] = fs->make<TH1F>("eop", "eop", 40, 0., 1.);
89  histContainer_["clus"] = fs->make<TH1F>("clus", "clus", 40, 0., 1.);
90  histContainer_["eIDs"] = fs->make<TH1F>("eIDs", "eIDS", 5, 0., 5.);
91 }
92 
94 
static const std::string kSharedResource
Definition: TFileService.h:76
void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< std::string, TH1F * > histContainer_
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< pat::Electron > > srcToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
PatZjetsElectronAnalyzer(const edm::ParameterSet &)