CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ZtoEEEventSelector.cc
Go to the documentation of this file.
1 // Z->ee Filter
10 #include "TLorentzVector.h"
13 
16 #include "TH1.h"
18 
19 using namespace std;
20 using namespace edm;
21 
23  : electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
24  bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
25  electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
26  bsToken_(consumes<reco::BeamSpot>(bsTag_)),
27  maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.4)),
28  minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
29  maxDeltaPhiInEB_(ps.getUntrackedParameter<double>("maxDeltaPhiInEB", .15)),
30  maxDeltaEtaInEB_(ps.getUntrackedParameter<double>("maxDeltaEtaInEB", .007)),
31  maxHOEEB_(ps.getUntrackedParameter<double>("maxHOEEB", .12)),
32  maxSigmaiEiEEB_(ps.getUntrackedParameter<double>("maxSigmaiEiEEB", .01)),
33  maxDeltaPhiInEE_(ps.getUntrackedParameter<double>("maxDeltaPhiInEE", .1)),
34  maxDeltaEtaInEE_(ps.getUntrackedParameter<double>("maxDeltaEtaInEE", .009)),
35  maxHOEEE_(ps.getUntrackedParameter<double>("maxHOEEB_", .10)),
36  maxSigmaiEiEEE_(ps.getUntrackedParameter<double>("maxSigmaiEiEEE", .03)),
37  maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 10)),
38  maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
39  maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
40  minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
41  minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
42  maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
43  minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
44  minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 60)),
45  maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 120)) {}
46 
48  // Read Electron Collection
50  iEvent.getByToken(electronToken_, electronColl);
51 
53  iEvent.getByToken(bsToken_, beamSpot);
54 
55  std::vector<TLorentzVector> list;
56  std::vector<int> chrgeList;
57 
58  if (electronColl.isValid()) {
59  for (auto const& ele : *electronColl) {
60  if (!ele.ecalDriven())
61  continue;
62  if (ele.pt() < minPt_)
63  continue;
64  // set a max Eta cut
65  if (!(ele.isEB() || ele.isEE()))
66  continue;
67 
68  double hOverE = ele.hadronicOverEm();
69  double sigmaee = ele.sigmaIetaIeta();
70  double deltaPhiIn = ele.deltaPhiSuperClusterTrackAtVtx();
71  double deltaEtaIn = ele.deltaEtaSuperClusterTrackAtVtx();
72 
73  // separate cut for barrel and endcap
74  if (ele.isEB()) {
75  if (fabs(deltaPhiIn) >= maxDeltaPhiInEB_ && fabs(deltaEtaIn) >= maxDeltaEtaInEB_ && hOverE >= maxHOEEB_ &&
76  sigmaee >= maxSigmaiEiEEB_)
77  continue;
78  } else if (ele.isEE()) {
79  if (fabs(deltaPhiIn) >= maxDeltaPhiInEE_ && fabs(deltaEtaIn) >= maxDeltaEtaInEE_ && hOverE >= maxHOEEE_ &&
80  sigmaee >= maxSigmaiEiEEE_)
81  continue;
82  }
83 
84  reco::GsfTrackRef trk = ele.gsfTrack();
85  if (!trk.isNonnull())
86  continue; // only electrons with tracks
87  double chi2 = trk->chi2();
88  double ndof = trk->ndof();
89  double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
90  if (chbyndof >= maxNormChi2_)
91  continue;
92 
93  double trkd0 = trk->d0();
94  if (beamSpot.isValid()) {
95  trkd0 = -(trk->dxy(beamSpot->position()));
96  } else {
97  edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get BeamSpot for label: " << bsTag_;
98  }
99  if (std::fabs(trkd0) >= maxD0_)
100  continue;
101 
102  const reco::HitPattern& hitp = trk->hitPattern();
103  int nPixelHits = hitp.numberOfValidPixelHits();
104  if (nPixelHits < minPixelHits_)
105  continue;
106 
107  int nStripHits = hitp.numberOfValidStripHits();
108  if (nStripHits < minStripHits_)
109  continue;
110 
111  // PF Isolation
112  reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
113  float absiso =
114  pfIso.sumChargedHadronPt + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
115  float eiso = absiso / (ele.pt());
116  if (eiso > maxIso_)
117  continue;
118 
119  TLorentzVector lv;
120  lv.SetPtEtaPhiE(ele.pt(), ele.eta(), ele.phi(), ele.energy());
121  list.push_back(lv);
122  chrgeList.push_back(ele.charge());
123  }
124  } else {
125  edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get ElectronCollection for label: " << electronTag_;
126  }
127  if (list.size() < 2)
128  return false;
129  if (chrgeList[0] + chrgeList[1] != 0)
130  return false;
131 
132  if (list[0].Pt() < minPtHighest_)
133  return false;
134  TLorentzVector zv = list[0] + list[1];
135  if (zv.M() < minInvMass_ || zv.M() > maxInvMass_)
136  return false;
137 
138  return true;
139 }
140 // Define this as a plug-in
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
const edm::InputTag electronTag_
const double maxSigmaiEiEEE_
const double maxDeltaPhiInEB_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
const edm::EDGetTokenT< reco::GsfElectronCollection > electronToken_
const double minPtHighest_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const double maxDeltaPhiInEE_
float *__restrict__ zv
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
Definition: GsfElectron.h:670
int numberOfValidStripHits() const
Definition: HitPattern.h:837
const double minInvMass_
const double maxDeltaEtaInEE_
Log< level::Error, false > LogError
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
int iEvent
Definition: GenABIO.cc:224
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:665
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:664
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
const double maxDeltaEtaInEB_
bool isValid() const
Definition: HandleBase.h:70
const double maxInvMass_
bool filter(edm::Event &, edm::EventSetup const &) override
const edm::InputTag bsTag_
ZtoEEEventSelector(const edm::ParameterSet &)
int numberOfValidPixelHits() const
Definition: HitPattern.h:825
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:663
const double maxSigmaiEiEEB_
const double maxNormChi2_