CMS 3D CMS Logo

PatElectronAnalyzer.cc
Go to the documentation of this file.
1 #include "TH1.h"
2 #include "TH2.h"
3 #include "TMath.h"
4 
11 
14 
16 public:
17  explicit PatElectronAnalyzer(const edm::ParameterSet &);
18  ~PatElectronAnalyzer() override;
19 
20 private:
21  void beginJob() override;
22  void analyze(const edm::Event &, const edm::EventSetup &) override;
23  void endJob() override;
24 
25  // restrictions for the electron to be
26  // considered
27  double minPt_;
28  double maxEta_;
29 
30  // decide in which mode to run the analyzer
31  // 0 : genMatch, 1 : tagAndProbe are
32  // supported depending on the line comments
33  unsigned int mode_;
34 
35  // choose a given electronID for the electron
36  // in consideration; the following types are
37  // available:
38  // * eidRobustLoose
39  // * eidRobustTight
40  // * eidLoose
41  // * eidTight
42  // * eidRobustHighEnergy
44 
45  // source of electrons
47  // source of generator particles
49 
52 
53  // internal variables for genMatchMode and
54  // tagAndProbeMode
55  double maxDeltaR_;
56  double maxDeltaM_;
57  double maxTagIso_;
58 
59  // book histograms of interest
60  TH1I *nr_;
61  TH1F *pt_;
62  TH1F *eta_;
63  TH1F *phi_;
64  TH1F *genPt_;
65  TH1F *genEta_;
66  TH1F *genPhi_;
67  TH1F *deltaR_;
68  TH1F *isoTag_;
69  TH1F *invMass_;
70  TH1F *deltaPhi_;
71 };
72 
73 #include "Math/VectorUtil.h"
74 
76  : minPt_(cfg.getParameter<double>("minPt")),
77  maxEta_(cfg.getParameter<double>("maxEta")),
78  mode_(cfg.getParameter<unsigned int>("mode")),
79  electronID_(cfg.getParameter<std::string>("electronID")),
80  electronSrcToken_(consumes<std::vector<pat::Electron> >(cfg.getParameter<edm::InputTag>("electronSrc"))),
81  particleSrcToken_(consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("particleSrc"))),
82  genMatchMode_(cfg.getParameter<edm::ParameterSet>("genMatchMode")),
83  tagAndProbeMode_(cfg.getParameter<edm::ParameterSet>("tagAndProbeMode")) {
84  // complete the configuration of the analyzer
85  maxDeltaR_ = genMatchMode_.getParameter<double>("maxDeltaR");
86  maxDeltaM_ = tagAndProbeMode_.getParameter<double>("maxDeltaM");
87  maxTagIso_ = tagAndProbeMode_.getParameter<double>("maxTagIso");
88 
89  // register histograms to the TFileService
91  nr_ = fs->make<TH1I>("nr", "nr", 10, 0, 10);
92  pt_ = fs->make<TH1F>("pt", "pt", 20, 0., 100.);
93  eta_ = fs->make<TH1F>("eta", "eta", 30, -3., 3.);
94  phi_ = fs->make<TH1F>("phi", "phi", 35, -3.5, 3.5);
95  genPt_ = fs->make<TH1F>("genPt", "pt", 20, 0., 100.);
96  genEta_ = fs->make<TH1F>("genEta", "eta", 30, -3., 3.);
97  genPhi_ = fs->make<TH1F>("genPhi", "phi", 35, -3.5, 3.5);
98  deltaR_ = fs->make<TH1F>("deltaR", "log(dR)", 50, -5., 0.);
99  isoTag_ = fs->make<TH1F>("isoTag", "iso", 50, 0., 10.);
100  invMass_ = fs->make<TH1F>("invMass", "m", 100, 50., 150.);
101  deltaPhi_ = fs->make<TH1F>("deltaPhi", "deltaPhi", 100, -3.5, 3.5);
102 }
103 
105 
107  // get electron collection
110  // get generator particle collection
113 
114  nr_->Fill(electrons->size());
115 
116  // ----------------------------------------------------------------------
117  //
118  // First Part Mode 0: genMatch
119  //
120  // ----------------------------------------------------------------------
121  if (mode_ == 0) {
122  // loop generator particles
123  for (reco::GenParticleCollection::const_iterator part = particles->begin(); part != particles->end(); ++part) {
124  // only loop stable electrons
125  if (part->status() == 1 && abs(part->pdgId()) == 11) {
126  if (part->pt() > minPt_ && fabs(part->eta()) < maxEta_) {
127  genPt_->Fill(part->pt());
128  genEta_->Fill(part->eta());
129  genPhi_->Fill(part->phi());
130  }
131  }
132  }
133 
134  // loop electrons
135  for (std::vector<pat::Electron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
136  if (elec->genLepton()) {
137  float deltaR = ROOT::Math::VectorUtil::DeltaR(elec->genLepton()->p4(), elec->p4());
138  deltaR_->Fill(TMath::Log10(deltaR));
139  if (deltaR < maxDeltaR_) {
140  if (electronID_ != "none") {
141  if (elec->electronID(electronID_) < 0.5)
142  continue;
143  }
144  if (elec->pt() > minPt_ && fabs(elec->eta()) < maxEta_) {
145  pt_->Fill(elec->pt());
146  eta_->Fill(elec->eta());
147  phi_->Fill(elec->phi());
148  }
149  }
150  }
151  }
152  }
153 
154  // ----------------------------------------------------------------------
155  //
156  // Second Part Mode 1: tagAndProbe
157  //
158  // ----------------------------------------------------------------------
159  if (mode_ == 1) {
160  // loop tag electron
161  for (std::vector<pat::Electron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
162  isoTag_->Fill(elec->trackIso());
163  if (elec->trackIso() < maxTagIso_ && elec->electronID("eidTight") > 0.5) {
164  // loop probe electron
165  for (std::vector<pat::Electron>::const_iterator probe = electrons->begin(); probe != electrons->end();
166  ++probe) {
167  // skip the tag electron itself
168  if (probe == elec)
169  continue;
170 
171  float zMass = (probe->p4() + elec->p4()).mass();
172  invMass_->Fill(zMass);
173  float deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(elec->p4(), probe->p4());
174  deltaPhi_->Fill(deltaPhi);
175 
176  // check for the Z mass
177  if (fabs(zMass - 90.) < maxDeltaM_) {
178  if (electronID_ != "none") {
179  if (probe->electronID(electronID_) < 0.5)
180  continue;
181  }
182  if (probe->pt() > minPt_ && fabs(probe->eta()) < maxEta_) {
183  pt_->Fill(probe->pt());
184  eta_->Fill(probe->eta());
185  phi_->Fill(probe->phi());
186  }
187  }
188  }
189  }
190  }
191  }
192 }
193 
195 
197 
PatElectronAnalyzer::beginJob
void beginJob() override
Definition: PatElectronAnalyzer.cc:194
PatElectronAnalyzer::maxDeltaR_
double maxDeltaR_
Definition: PatElectronAnalyzer.cc:55
PatElectronAnalyzer::isoTag_
TH1F * isoTag_
Definition: PatElectronAnalyzer.cc:68
edm::EDGetTokenT
Definition: EDGetToken.h:33
Electron
Definition: Electron.py:1
edm
HLT enums.
Definition: AlignableModifier.h:19
PatElectronAnalyzer::invMass_
TH1F * invMass_
Definition: PatElectronAnalyzer.cc:69
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:85964
reco::GenParticleCollection
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
Definition: GenParticleFwd.h:13
PatElectronAnalyzer::maxTagIso_
double maxTagIso_
Definition: PatElectronAnalyzer.cc:57
PatElectronAnalyzer::maxEta_
double maxEta_
Definition: PatElectronAnalyzer.cc:28
EDAnalyzer.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
PatElectronAnalyzer::genEta_
TH1F * genEta_
Definition: PatElectronAnalyzer.cc:65
edm::Handle
Definition: AssociativeIterator.h:50
PatElectronAnalyzer::particleSrcToken_
edm::EDGetTokenT< reco::GenParticleCollection > particleSrcToken_
Definition: PatElectronAnalyzer.cc:48
ecalTrigSettings_cff.particles
particles
Definition: ecalTrigSettings_cff.py:11
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
HLT_FULL_cff.DeltaPhi
DeltaPhi
Definition: HLT_FULL_cff.py:465
edm::EDAnalyzer
Definition: EDAnalyzer.h:28
GenParticle.h
PatElectronAnalyzer::maxDeltaM_
double maxDeltaM_
Definition: PatElectronAnalyzer.cc:56
MakerMacros.h
PatElectronAnalyzer::tagAndProbeMode_
edm::ParameterSet tagAndProbeMode_
Definition: PatElectronAnalyzer.cc:51
part
part
Definition: HCALResponse.h:20
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
PatElectronAnalyzer::genPhi_
TH1F * genPhi_
Definition: PatElectronAnalyzer.cc:66
Service.h
PatElectronAnalyzer::genPt_
TH1F * genPt_
Definition: PatElectronAnalyzer.cc:64
PatElectronAnalyzer::eta_
TH1F * eta_
Definition: PatElectronAnalyzer.cc:62
PatElectronAnalyzer::electronSrcToken_
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
Definition: PatElectronAnalyzer.cc:46
PatElectronAnalyzer::PatElectronAnalyzer
PatElectronAnalyzer(const edm::ParameterSet &)
Definition: PatElectronAnalyzer.cc:75
PatElectronAnalyzer::deltaPhi_
TH1F * deltaPhi_
Definition: PatElectronAnalyzer.cc:70
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:531
PatElectronAnalyzer::endJob
void endJob() override
Definition: PatElectronAnalyzer.cc:196
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
TFileService.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
ParameterSet
Definition: Functions.h:16
PatElectronAnalyzer::nr_
TH1I * nr_
Definition: PatElectronAnalyzer.cc:60
edm::Service< TFileService >
createfilelist.int
int
Definition: createfilelist.py:10
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::EventSetup
Definition: EventSetup.h:57
pat
Definition: HeavyIon.h:7
electronAnalyzer_cfi.DeltaR
DeltaR
Definition: electronAnalyzer_cfi.py:33
PatElectronAnalyzer::deltaR_
TH1F * deltaR_
Definition: PatElectronAnalyzer.cc:67
InputTag.h
looper.cfg
cfg
Definition: looper.py:297
std
Definition: JetResolutionObject.h:76
PatElectronAnalyzer
Definition: PatElectronAnalyzer.cc:15
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
PatElectronAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: PatElectronAnalyzer.cc:106
PatElectronAnalyzer::pt_
TH1F * pt_
Definition: PatElectronAnalyzer.cc:61
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Electron.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
PatElectronAnalyzer::electronID_
std::string electronID_
Definition: PatElectronAnalyzer.cc:43
PatElectronAnalyzer::phi_
TH1F * phi_
Definition: PatElectronAnalyzer.cc:63
PatElectronAnalyzer::mode_
unsigned int mode_
Definition: PatElectronAnalyzer.cc:33
edm::Event
Definition: Event.h:73
PatElectronAnalyzer::genMatchMode_
edm::ParameterSet genMatchMode_
Definition: PatElectronAnalyzer.cc:50
PatElectronAnalyzer::~PatElectronAnalyzer
~PatElectronAnalyzer() override
Definition: PatElectronAnalyzer.cc:104
HLTMuonOfflineAnalyzer_cfi.zMass
zMass
Definition: HLTMuonOfflineAnalyzer_cfi.py:101
PatElectronAnalyzer::minPt_
double minPt_
Definition: PatElectronAnalyzer.cc:27