CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
WtoLNuSelector.cc
Go to the documentation of this file.
1 // W->lnu Event Selector
2 
18 
19 #include "TLorentzVector.h"
20 #include "TH1.h"
21 #include "TMath.h"
22 
24 
26  electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
27  bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
28  muonTag_(ps.getUntrackedParameter<edm::InputTag>("muonInputTag", edm::InputTag("muons"))),
29  //pfmetTag_(ps.getUntrackedParameter<edm::InputTag>("pfmetTag", edm::InputTag("pfMet"))),
30  pfmetTag_(ps.getUntrackedParameter<edm::InputTag>("pfmetTag", edm::InputTag("pfMetT1T2Txy"))),
31  electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
32  bsToken_(consumes<reco::BeamSpot>(bsTag_)),
33  muonToken_(consumes<reco::MuonCollection>(muonTag_)),
34  pfmetToken_(consumes<reco::PFMETCollection>(pfmetTag_))
35 {
36 }
37 
39  // Read Electron Collection
41  iEvent.getByToken(electronToken_, electronColl);
42 
44  iEvent.getByToken(bsToken_, beamSpot);
45 
46  std::vector<TLorentzVector> eleList;
47  if (electronColl.isValid()) {
48  for (auto const& ele : *electronColl) {
49  if (!ele.ecalDriven() ) continue;
50 
51  double hOverE = ele.hadronicOverEm();
52  double sigmaee = ele.sigmaIetaIeta();
53  double deltaPhiIn = ele.deltaPhiSuperClusterTrackAtVtx();
54  double deltaEtaIn = ele.deltaEtaSuperClusterTrackAtVtx();
55 
56  // separate cut for barrel and endcap
57  if (ele.isEB()) {
58  if (std::fabs(deltaPhiIn) >= .15
59  && std::fabs(deltaEtaIn) >= .007
60  && hOverE >= .12
61  && sigmaee >= .01) continue;
62  }
63  else if (ele.isEE()) {
64  if (std::fabs(deltaPhiIn) >= .10
65  && std::fabs(deltaEtaIn) >= .009
66  && hOverE >= .10
67  && sigmaee >= .03) continue;
68  }
69 
70  reco::GsfTrackRef trk = ele.gsfTrack();
71  if (!trk.isNonnull() ) continue; // only electrons wd tracks
72  double chi2 = trk->chi2();
73  double ndof = trk->ndof();
74  double chbyndof = (ndof > 0) ? chi2/ndof : 0;
75  double trkd0 = trk->d0();
76  double trkdz = trk->dz();
77  if (beamSpot.isValid()) {
78  trkd0 = -(trk->dxy(beamSpot->position()));
79  trkdz = trk->dz(beamSpot->position());
80  }
81  else {
82  edm::LogError("WtoLNuSelector") << "Error >> Failed to get BeamSpot for label: "
83  << bsTag_;
84  }
85  if (chbyndof >= 10 || std::fabs(trkd0) >= 0.02 || std::fabs(trkdz) >= 20) continue;
86  const reco::HitPattern& hitp = trk->hitPattern();
87  int nPixelHits = hitp.numberOfValidPixelHits();
88  int nStripHits = hitp.numberOfValidStripHits();
89  if (nPixelHits < 1 || nStripHits < 8) continue;
90 
91  // PF Isolation
92  reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
93  float absiso = pfIso.sumChargedHadronPt
94  + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
95  float eiso = absiso/ele.pt();
96  if (eiso > 0.3) continue;
97 
98  TLorentzVector le;
99  le.SetPtEtaPhiE(ele.pt(), ele.eta(), ele.phi(), ele.energy());
100  eleList.push_back(le);
101  }
102  }
103  else {
104  edm::LogError("WtoLNuSelector") << "Error >> Failed to get ElectronCollection for label: "
105  << electronTag_;
106  }
107 
108  // Read Muon Collection
110  iEvent.getByToken(muonToken_, muonColl);
111 
112  std::vector<TLorentzVector> muList;
113  if (muonColl.isValid()) {
114  for (auto const& mu : *muonColl) {
115  if (!mu.isGlobalMuon()
116  || !mu.isPFMuon()
117  || std::fabs(mu.eta()) > 2.1
118  || mu.pt() <= 5) continue;
119 
120  reco::TrackRef gtrkref = mu.globalTrack();
121  if (!gtrkref.isNonnull()) continue;
122  const reco::Track* gtk = &(*gtrkref);
123  double chi2 = gtk->chi2();
124  double ndof = gtk->ndof();
125  double chbyndof = (ndof > 0) ? chi2/ndof : 0;
126 
127  const reco::HitPattern& hitp = gtk->hitPattern();
128  int nPixelHits = hitp.numberOfValidPixelHits();
129  int nStripHits = hitp.numberOfValidStripHits();
130 
131  reco::TrackRef itrkref = mu.innerTrack(); // tracker segment only
132  if (!itrkref.isNonnull()) continue;
133  const reco::Track* tk = &(*itrkref);
134  double trkd0 = tk->d0();
135  double trkdz = tk->dz();
136  if (beamSpot.isValid()) {
137  trkd0 = -(tk->dxy(beamSpot->position()));
138  trkdz = tk->dz(beamSpot->position());
139  }
140  // Hits/section in the muon chamber
141  int nChambers = mu.numberOfChambers();
142  int nMatches = mu.numberOfMatches();
143  int nMatchedStations = mu.numberOfMatchedStations();
144 
145  // PF Isolation
146  const reco::MuonPFIsolation& pfIso04 = mu.pfIsolationR04();
147  double absiso = pfIso04.sumChargedParticlePt + std::max(0.0, pfIso04.sumNeutralHadronEt + pfIso04.sumPhotonEt - 0.5 * pfIso04.sumPUPt);
148 
149  // Final selection
150  if (chbyndof < 10
151  && std::fabs(trkd0) < 0.02
152  && std::fabs(trkdz) < 20.0
153  && nPixelHits > 1
154  && nStripHits > 8
155  && nChambers > 2
156  && nMatches > 2
157  && nMatchedStations > 2
158  && absiso/mu.pt() < 0.3)
159  {
160  TLorentzVector lm;
161  lm.SetPtEtaPhiE(mu.pt(), mu.eta(), mu.phi(), mu.energy());
162  muList.push_back(lm);
163  }
164  }
165  }
166  else {
167  edm::LogError("WtoLNuSelector") << "Error >> Failed to get MuonCollection for label: " << muonTag_;
168  }
169 
170  // Require either a high pt electron or muon
171  if (eleList.size() < 1 && muList.size() < 1) return false;
172 
173  // Both should not be present at the same time
174  if ((eleList.size() > 0 && eleList[0].Pt() > 20) &&
175  (muList.size() > 0 && muList[0].Pt() > 20)) return false;
176 
177  // find the high pt lepton
178  TLorentzVector vlep;
179  if (eleList.size() > 0 && muList.size() > 0) {
180  vlep = (eleList[0].Pt() > muList[0].Pt()) ? eleList[0] : muList[0];
181  }
182  else if (eleList.size() > 0) {
183  vlep = eleList[0];
184  }
185  else {
186  vlep = muList[0];
187  }
188  if (vlep.Pt() < 20) return false;
189 
191  iEvent.getByToken(pfmetToken_, pfColl);
192 
193  if (pfColl.isValid()) {
194  double mt = getMt(vlep, pfColl->front());
195  if (mt < 60 || mt > 80) return false;
196  }
197  else {
198  edm::LogError("WtoLNuSelector") << "Error >> Failed to get PFMETCollection for label: " << pfmetTag_;
199  return false;
200  }
201 
202  return true;
203 }
204 double WtoLNuSelector::getMt(const TLorentzVector& vlep, const reco::PFMET& obj) {
205  double met = obj.et();
206  double phi = obj.phi();
207 
208  TLorentzVector vmet;
209  double metx = met * std::cos(phi);
210  double mety = met * std::sin(phi);
211  vmet.SetPxPyPzE(metx, mety, 0.0, met);
212 
213  // transverse mass
214  TLorentzVector vw = vlep + vmet;
215 
216  return std::sqrt(2 * vlep.Et() * met * (1 - std::cos(deltaPhi(vlep.Phi(), phi))));
217 }
218 // Define this as a plug-in
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
double getMt(const TLorentzVector &vlep, const reco::PFMET &obj)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:592
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
Definition: GsfElectron.h:580
int numberOfValidStripHits() const
Definition: HitPattern.h:853
virtual double phi() const final
momentum azimuthal angle
const edm::InputTag electronTag_
const edm::EDGetTokenT< reco::GsfElectronCollection > electronToken_
float sumPhotonEt
sum pt of PF photons
WtoLNuSelector(const edm::ParameterSet &)
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
float sumNeutralHadronEt
sum pt of neutral hadrons
float sumChargedParticlePt
sum-pt of charged Particles(inludes e/mu)
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
int iEvent
Definition: GenABIO.cc:230
const edm::EDGetTokenT< reco::MuonCollection > muonToken_
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:544
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:550
T sqrt(T t)
Definition: SSEVec.h:18
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:575
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:574
const int mu
Definition: Constants.h:22
bool isValid() const
Definition: HandleBase.h:75
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:604
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
const edm::InputTag bsTag_
bool filter(edm::Event &, edm::EventSetup const &) override
const edm::InputTag muonTag_
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
virtual double et() const final
transverse energy
const edm::InputTag pfmetTag_
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:586
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:573
const edm::EDGetTokenT< reco::PFMETCollection > pfmetToken_
Collection of PF MET.