CMS 3D CMS Logo

ZtoEEEventSelector.cc
Go to the documentation of this file.
1 // Z->ee Filter
10 #include "TLorentzVector.h"
14 
17 #include "TH1.h"
19 
20 using namespace std;
21 using namespace edm;
22 
24  electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
25  bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
26  electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
27  bsToken_(consumes<reco::BeamSpot>(bsTag_)),
28  maxEta_(ps.getUntrackedParameter<double>("maxEta", 2.4)),
29  minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
30  maxDeltaPhiInEB_(ps.getUntrackedParameter<double>("maxDeltaPhiInEB", .15)),
31  maxDeltaEtaInEB_(ps.getUntrackedParameter<double>("maxDeltaEtaInEB", .007)),
32  maxHOEEB_(ps.getUntrackedParameter<double>("maxHOEEB", .12)),
33  maxSigmaiEiEEB_(ps.getUntrackedParameter<double>("maxSigmaiEiEEB", .01)),
34  maxDeltaPhiInEE_(ps.getUntrackedParameter<double>("maxDeltaPhiInEE", .1)),
35  maxDeltaEtaInEE_(ps.getUntrackedParameter<double>("maxDeltaEtaInEE", .009)),
36  maxHOEEE_(ps.getUntrackedParameter<double>("maxHOEEB_", .10)),
37  maxSigmaiEiEEE_(ps.getUntrackedParameter<double>("maxSigmaiEiEEE", .03)),
38  maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 10)),
39  maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
40  maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
41  minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
42  minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
43  maxIso_(ps.getUntrackedParameter<double>("maxIso", 0.3)),
44  minPtHighest_(ps.getUntrackedParameter<double>("minPtHighest", 24)),
45  minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 60)),
46  maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 120))
47 {
48 }
49 
51  // Read Electron Collection
53  iEvent.getByToken(electronToken_, electronColl);
54 
56  iEvent.getByToken(bsToken_, beamSpot);
57 
58  std::vector<TLorentzVector> list;
59  std::vector<int> chrgeList;
60 
61  if (electronColl.isValid()) {
62  for (auto const& ele: *electronColl) {
63  if (!ele.ecalDriven()) continue;
64  if (ele.pt() < minPt_) continue;
65  // set a max Eta cut
66  if (!(ele.isEB() || ele.isEE())) 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_
76  && fabs(deltaEtaIn) >= maxDeltaEtaInEB_
77  && hOverE >= maxHOEEB_
78  && sigmaee >= maxSigmaiEiEEB_) continue;
79  }
80  else if (ele.isEE()) {
81  if (fabs(deltaPhiIn) >= maxDeltaPhiInEE_
82  && fabs(deltaEtaIn) >= maxDeltaEtaInEE_
83  && hOverE >= maxHOEEE_
84  && sigmaee >= maxSigmaiEiEEE_) continue;
85  }
86 
87  reco::GsfTrackRef trk = ele.gsfTrack();
88  if (!trk.isNonnull() ) continue; // only electrons with tracks
89  double chi2 = trk->chi2();
90  double ndof = trk->ndof();
91  double chbyndof = (ndof > 0) ? chi2/ndof : 0;
92  if (chbyndof >= maxNormChi2_) continue;
93 
94  double trkd0 = trk->d0();
95  if (beamSpot.isValid()) {
96  trkd0 = -(trk->dxy(beamSpot->position()));
97  }
98  else {
99  edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get BeamSpot for label: "
100  << bsTag_;
101  }
102  if (std::fabs(trkd0) >= maxD0_) continue;
103 
104  const reco::HitPattern& hitp = trk->hitPattern();
105  int nPixelHits = hitp.numberOfValidPixelHits();
106  if (nPixelHits < minPixelHits_) continue;
107 
108  int nStripHits = hitp.numberOfValidStripHits();
109  if (nStripHits < minStripHits_) continue;
110 
111  // PF Isolation
112  reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
113  float absiso = pfIso.sumChargedHadronPt + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
114  float eiso = absiso/(ele.pt());
115  if (eiso > maxIso_) continue;
116 
117  TLorentzVector lv;
118  lv.SetPtEtaPhiE(ele.pt(), ele.eta(), ele.phi(), ele.energy());
119  list.push_back(lv);
120  chrgeList.push_back(ele.charge());
121  }
122  }
123  else {
124  edm::LogError("ZtoEEEventSelector") << "Error >> Failed to get ElectronCollection for label: "
125  << electronTag_;
126  }
127  if (list.size() < 2) return false;
128  if (chrgeList[0] + chrgeList[1] != 0) return false;
129 
130  if (list[0].Pt() < minPtHighest_) return false;
131  TLorentzVector zv = list[0] + list[1];
132  if (zv.M() < minInvMass_ || zv.M() > maxInvMass_) return false;
133 
134  return true;
135 }
136 // Define this as a plug-in
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
const edm::InputTag electronTag_
const double maxSigmaiEiEEE_
const double maxDeltaPhiInEB_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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:646
int numberOfValidStripHits() const
Definition: HitPattern.h:931
const double minInvMass_
const double maxDeltaEtaInEE_
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:641
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:640
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
const double maxDeltaEtaInEB_
bool isValid() const
Definition: HandleBase.h:74
const double maxInvMass_
bool filter(edm::Event &, edm::EventSetup const &) override
const edm::InputTag bsTag_
fixed size matrix
HLT enums.
ZtoEEEventSelector(const edm::ParameterSet &)
int numberOfValidPixelHits() const
Definition: HitPattern.h:916
const Point & position() const
position
Definition: BeamSpot.h:62
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:639
const double maxSigmaiEiEEB_
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger list("!*","!HLTx*"if it matches 2 triggers or more) will accept the event if all the matching triggers are FAIL.It will reject the event if any of the triggers are PASS or EXCEPTION(this matches the behavior of"!*"before the partial wildcard feature was incorporated).Triggers which are in the READY state are completely ignored.(READY should never be returned since the trigger paths have been run
const double maxNormChi2_