1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
|
// Z->ee Filter
// user includes
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/EgammaCandidates/interface/Electron.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
#include "DataFormats/EgammaReco/interface/BasicCluster.h"
#include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
#include "DataFormats/TrackReco/interface/HitPattern.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/stream/EDFilter.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
// ROOT includes
#include "TLorentzVector.h"
class ZtoEEEventSelector : public edm::stream::EDFilter<> {
public:
explicit ZtoEEEventSelector(const edm::ParameterSet&);
bool filter(edm::Event&, edm::EventSetup const&) override;
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
// module config parameters
const edm::InputTag electronTag_;
const edm::InputTag bsTag_;
const edm::EDGetTokenT<reco::GsfElectronCollection> electronToken_;
const edm::EDGetTokenT<reco::BeamSpot> bsToken_;
const double maxEta_;
const double minPt_;
const double maxDeltaPhiInEB_;
const double maxDeltaEtaInEB_;
const double maxHOEEB_;
const double maxSigmaiEiEEB_;
const double maxDeltaPhiInEE_;
const double maxDeltaEtaInEE_;
const double maxHOEEE_;
const double maxSigmaiEiEEE_;
const double maxNormChi2_;
const double maxD0_;
const double maxDz_;
const int minPixelHits_;
const int minStripHits_;
const double maxIso_;
const double minPtHighest_;
const double minInvMass_;
const double maxInvMass_;
};
using namespace std;
using namespace edm;
void ZtoEEEventSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.addUntracked<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"));
desc.addUntracked<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"));
desc.addUntracked<double>("maxEta", 2.4);
desc.addUntracked<double>("minPt", 5);
desc.addUntracked<double>("maxDeltaPhiInEB", .15);
desc.addUntracked<double>("maxDeltaEtaInEB", .007);
desc.addUntracked<double>("maxHOEEB", .12);
desc.addUntracked<double>("maxSigmaiEiEEB", .01);
desc.addUntracked<double>("maxDeltaPhiInEE", .1);
desc.addUntracked<double>("maxDeltaEtaInEE", .009);
desc.addUntracked<double>("maxHOEEB_", .10);
desc.addUntracked<double>("maxSigmaiEiEEE", .03);
desc.addUntracked<double>("maxNormChi2", 1000);
desc.addUntracked<double>("maxD0", 0.02);
desc.addUntracked<double>("maxDz", 20.);
desc.addUntracked<uint32_t>("minPixelHits", 1);
desc.addUntracked<uint32_t>("minStripHits", 8);
desc.addUntracked<double>("maxIso", 0.3);
desc.addUntracked<double>("minPtHighest", 24);
desc.addUntracked<double>("minInvMass", 75);
desc.addUntracked<double>("maxInvMass", 105);
descriptions.addWithDefaultLabel(desc);
}
ZtoEEEventSelector::ZtoEEEventSelector(const edm::ParameterSet& ps)
: electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
bsToken_(consumes<reco::BeamSpot>(bsTag_)),
maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.4)),
minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
maxDeltaPhiInEB_(ps.getUntrackedParameter<double>("maxDeltaPhiInEB", .15)),
maxDeltaEtaInEB_(ps.getUntrackedParameter<double>("maxDeltaEtaInEB", .007)),
maxHOEEB_(ps.getUntrackedParameter<double>("maxHOEEB", .12)),
maxSigmaiEiEEB_(ps.getUntrackedParameter<double>("maxSigmaiEiEEB", .01)),
maxDeltaPhiInEE_(ps.getUntrackedParameter<double>("maxDeltaPhiInEE", .1)),
maxDeltaEtaInEE_(ps.getUntrackedParameter<double>("maxDeltaEtaInEE", .009)),
maxHOEEE_(ps.getUntrackedParameter<double>("maxHOEEB_", .10)),
maxSigmaiEiEEE_(ps.getUntrackedParameter<double>("maxSigmaiEiEEE", .03)),
maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 1000)),
maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 75)),
maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 105)) {}
bool ZtoEEEventSelector::filter(edm::Event& iEvent, edm::EventSetup const& iSetup) {
// Read Electron Collection
edm::Handle<reco::GsfElectronCollection> electronColl;
iEvent.getByToken(electronToken_, electronColl);
edm::Handle<reco::BeamSpot> beamSpot;
iEvent.getByToken(bsToken_, beamSpot);
std::vector<TLorentzVector> list;
std::vector<int> chrgeList;
if (electronColl.isValid()) {
for (auto const& ele : *electronColl) {
if (!ele.ecalDriven())
continue;
if (ele.pt() < minPt_)
continue;
// set a max Eta cut
if (!(ele.isEB() || ele.isEE()))
continue;
double hOverE = ele.hadronicOverEm();
double sigmaee = ele.sigmaIetaIeta();
double deltaPhiIn = ele.deltaPhiSuperClusterTrackAtVtx();
double deltaEtaIn = ele.deltaEtaSuperClusterTrackAtVtx();
// separate cut for barrel and endcap
if (ele.isEB()) {
if (fabs(deltaPhiIn) >= maxDeltaPhiInEB_ && fabs(deltaEtaIn) >= maxDeltaEtaInEB_ && hOverE >= maxHOEEB_ &&
sigmaee >= maxSigmaiEiEEB_)
continue;
} else if (ele.isEE()) {
if (fabs(deltaPhiIn) >= maxDeltaPhiInEE_ && fabs(deltaEtaIn) >= maxDeltaEtaInEE_ && hOverE >= maxHOEEE_ &&
sigmaee >= maxSigmaiEiEEE_)
continue;
}
reco::GsfTrackRef trk = ele.gsfTrack();
if (!trk.isNonnull())
continue; // only electrons with tracks
double chi2 = trk->chi2();
double ndof = trk->ndof();
double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
if (chbyndof >= maxNormChi2_)
continue;
double trkd0 = trk->d0();
if (beamSpot.isValid()) {
trkd0 = -(trk->dxy(beamSpot->position()));
} else {
edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get BeamSpot for label: " << bsTag_;
}
if (std::fabs(trkd0) >= maxD0_)
continue;
const reco::HitPattern& hitp = trk->hitPattern();
int nPixelHits = hitp.numberOfValidPixelHits();
if (nPixelHits < minPixelHits_)
continue;
int nStripHits = hitp.numberOfValidStripHits();
if (nStripHits < minStripHits_)
continue;
// PF Isolation
reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
float absiso =
pfIso.sumChargedHadronPt + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
float eiso = absiso / (ele.pt());
if (eiso > maxIso_)
continue;
TLorentzVector lv;
lv.SetPtEtaPhiE(ele.pt(), ele.eta(), ele.phi(), ele.energy());
list.push_back(lv);
chrgeList.push_back(ele.charge());
}
} else {
edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get ElectronCollection for label: " << electronTag_;
}
if (list.size() < 2)
return false;
if (chrgeList[0] + chrgeList[1] != 0)
return false;
if (list[0].Pt() < minPtHighest_)
return false;
TLorentzVector zv = list[0] + list[1];
if (zv.M() < minInvMass_ || zv.M() > maxInvMass_)
return false;
return true;
}
// Define this as a plug-in
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(ZtoEEEventSelector);
|