CMS 3D CMS Logo

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