CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
CutBasedElectronID Class Reference

#include <CutBasedElectronID.h>

Inheritance diagram for CutBasedElectronID:
ElectronIDAlgo

Public Member Functions

double cicSelection (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
 
int classify (const reco::GsfElectron *)
 
bool compute_cut (double x, double et, double cut_min, double cut_max, bool gtn=false)
 
 CutBasedElectronID ()
 
double result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
 
double robustSelection (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
 
void setup (const edm::ParameterSet &conf)
 
virtual ~CutBasedElectronID ()
 
- Public Member Functions inherited from ElectronIDAlgo
void baseSetup (const edm::ParameterSet &conf)
 
 ElectronIDAlgo ()
 
virtual ~ElectronIDAlgo ()
 

Private Attributes

edm::ParameterSet cuts_
 
bool newCategories_
 
std::string quality_
 
std::string type_
 
std::string version_
 
edm::InputTag verticesCollection_
 
bool wantBinning_
 

Additional Inherited Members

- Protected Member Functions inherited from ElectronIDAlgo
EcalClusterLazyTools getClusterShape (const edm::Event &, const edm::EventSetup &)
 
- Protected Attributes inherited from ElectronIDAlgo
edm::InputTag reducedBarrelRecHitCollection_
 
edm::InputTag reducedEndcapRecHitCollection_
 

Detailed Description

Definition at line 6 of file CutBasedElectronID.h.

Constructor & Destructor Documentation

CutBasedElectronID::CutBasedElectronID ( )
inline

Definition at line 10 of file CutBasedElectronID.h.

10 {};
virtual CutBasedElectronID::~CutBasedElectronID ( )
inlinevirtual

Definition at line 12 of file CutBasedElectronID.h.

12 {};

Member Function Documentation

double CutBasedElectronID::cicSelection ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
)

Definition at line 119 of file CutBasedElectronID.cc.

References newFWLiteAna::bin, classify(), compute_cut(), reco::GsfElectron::convDcot(), reco::GsfElectron::convDist(), align_tpl::cut, cutsInCategoriesElectronIdentificationV06_DataTuning_cfi::cutdcotdist, cutsInCategoriesElectronIdentificationV06_cfi::cuthoe, cutsInCategoriesElectronIdentificationV06_cfi::cuthoel, cuts_, cutsInCategoriesElectronIdentificationV06_cfi::cutsee, cutsInCategoriesElectronIdentificationV06_cfi::cutseel, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::dr03TkSumPt(), reco::GsfElectron::dr04EcalRecHitSumEt(), reco::GsfElectron::dr04HcalTowerSumEt(), reco::GsfElectron::eSeedClusterOverP(), reco::GsfElectron::eSuperClusterOverP(), eta(), funct::exp(), reco::GsfElectron::fbrem(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), max(), or, reco::LeafCandidate::p4(), funct::pow(), quality_, result(), reco::GsfElectron::sigmaEtaEta(), reco::GsfElectron::sigmaIetaIeta(), funct::sin(), reco::GsfElectron::superCluster(), type_, version_, verticesCollection_, and wantBinning_.

Referenced by result().

121  {
122 
123  double scTheta = (2*atan(exp(-electron->superCluster()->eta())));
124  double scEt = electron->superCluster()->energy()*sin(scTheta);
125 
126  double eta = electron->p4().Eta();
127  double eOverP = electron->eSuperClusterOverP();
128  double eSeedOverPin = electron->eSeedClusterOverP();
129  double fBrem = electron->fbrem();
130  double hOverE = electron->hadronicOverEm();
131  double sigmaee = electron->sigmaIetaIeta(); //sqrt(vLocCov[0]);
132  double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
133  double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
134 
135  double ip = 0;
136  int mishits = electron->gsfTrack()->trackerExpectedHitsInner().numberOfHits();
137  double tkIso = electron->dr03TkSumPt();
138  double ecalIso = electron->dr04EcalRecHitSumEt();
139  double hcalIso = electron->dr04HcalTowerSumEt();
140 
141  if (version_ == "V00") {
142  sigmaee = electron->sigmaEtaEta();//sqrt(vCov[0]);
143  if (electron->isEE())
144  sigmaee = sigmaee - 0.02*(fabs(eta) - 2.3); //correct sigmaetaeta dependence on eta in endcap
145  }
146 
147  if (version_ != "V01" or version_ != "V00") {
150  if (vtxH->size() != 0) {
151  reco::VertexRef vtx(vtxH, 0);
152  ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(),vtx->y(),vtx->z())));
153  } else
154  ip = fabs(electron->gsfTrack()->dxy());
155 
156  if (electron->isEB()) {
157  sigmaee = electron->sigmaIetaIeta(); //sqrt(vCov[0]);
158  }
159  }
160 
161  std::vector<double> cut;
162 
163  int cat = classify(electron);
164  int eb;
165 
166  if (electron->isEB())
167  eb = 0;
168  else
169  eb = 1;
170 
171  // LOOSE and TIGHT Selections
172  if (type_ == "classbased" && (version_ == "V01" || version_ == "V00")) {
173 
174  if ((eOverP < 0.8) && (fBrem < 0.2))
175  return 0.;
176 
177  cut = cuts_.getParameter<std::vector<double> >("hOverE");
178  if (hOverE > cut[cat+4*eb])
179  return 0.;
180 
181  cut = cuts_.getParameter<std::vector<double> >("sigmaEtaEta");
182  if (sigmaee > cut[cat+4*eb])
183  return 0.;
184 
185  cut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
186  if (eOverP < 1.5) {
187  if (fabs(deltaPhiIn) > cut[cat+4*eb])
188  return 0.;
189  } else {
190  if (fabs(deltaPhiIn) > cut[3+4*eb])
191  return 0.;
192  }
193 
194  cut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
195  if (fabs(deltaEtaIn) > cut[cat+4*eb])
196  return 0.;
197 
198  cut = cuts_.getParameter<std::vector<double> >("eSeedOverPin");
199  if (eSeedOverPin < cut[cat+4*eb])
200  return 0.;
201 
202  if (quality_ == "tight")
203  if (eOverP < 0.9*(1-fBrem))
204  return 0.;
205 
206  return 1.;
207  }
208 
209  if (type_ == "classbased" and version_ == "V02") {
210  double result = 0.;
211 
212  int bin = 0;
213 
214  if (scEt < 20.)
215  bin = 2;
216  else if (scEt > 30.)
217  bin = 0;
218  else
219  bin = 1;
220 
221  if (fBrem > 0)
222  eSeedOverPin = eSeedOverPin + fBrem;
223 
224  if (bin != 2) {
225  tkIso = tkIso*pow(40./scEt, 2);
226  ecalIso = ecalIso*pow(40./scEt, 2);
227  hcalIso = hcalIso*pow(40./scEt, 2);
228  }
229 
230  std::vector<double> cutTk = cuts_.getParameter<std::vector<double> >("cutisotk");
231  std::vector<double> cutEcal = cuts_.getParameter<std::vector<double> >("cutisoecal");
232  std::vector<double> cutHcal = cuts_.getParameter<std::vector<double> >("cutisohcal");
233  if ((tkIso > cutTk[cat+3*eb+bin*6]) ||
234  (ecalIso > cutEcal[cat+3*eb+bin*6]) ||
235  (hcalIso > cutHcal[cat+3*eb+bin*6]))
236  result = 0.;
237  else
238  result = 2.;
239 
240  if (fBrem > -2) {
241  std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
242  std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
243  std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
244  std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
245  std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
246  std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
247  std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip");
248  std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutmishits");
249  if ((hOverE < cuthoe[cat+3*eb+bin*6]) and
250  (sigmaee < cutsee[cat+3*eb+bin*6]) and
251  (fabs(deltaPhiIn) < cutdphi[cat+3*eb+bin*6]) and
252  (fabs(deltaEtaIn) < cutdeta[cat+3*eb+bin*6]) and
253  (eSeedOverPin > cuteopin[cat+3*eb+bin*6]) and
254  (ip < cutip[cat+3*eb+bin*6]) and
255  (mishits < cutmishits[cat+3*eb+bin*6]))
256  result = result + 1.;
257  }
258  return result;
259  }
260 
261  if (version_ == "V03" or version_ == "V04" or version_ == "V05") {
262  double result = 0.;
263 
264  int bin = 0;
265 
266  if (wantBinning_) {
267  if (scEt < 20.)
268  bin = 2;
269  else if (scEt > 30.)
270  bin = 0;
271  else
272  bin = 1;
273  }
274 
275  if (fBrem > 0)
276  eSeedOverPin = eSeedOverPin + fBrem;
277 
278  float iso_sum = tkIso + ecalIso + hcalIso;
279  float iso_sum_corrected = iso_sum*pow(40./scEt, 2);
280 
281  std::vector<double> cutIsoSum = cuts_.getParameter<std::vector<double> >("cutiso_sum");
282  std::vector<double> cutIsoSumCorr = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
283  if ((iso_sum < cutIsoSum[cat+bin*9]) and
284  (iso_sum_corrected < cutIsoSumCorr[cat+bin*9]))
285  result += 2.;
286 
287  if (fBrem > -2) {
288  std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
289  std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
290  std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
291  std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
292  std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
293  std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
294 
295  if ((hOverE < cuthoe[cat+bin*9]) and
296  (sigmaee < cutsee[cat+bin*9]) and
297  (fabs(deltaPhiIn) < cutdphi[cat+bin*9]) and
298  (fabs(deltaEtaIn) < cutdeta[cat+bin*9]) and
299  (eSeedOverPin > cuteopin[cat+bin*9]) and
300  (scEt > cutet[cat+bin*9]))
301  result += 1.;
302  }
303 
304  std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip_gsf");
305  if (ip < cutip[cat+bin*9])
306  result += 8;
307 
308  std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutfmishits");
309  std::vector<double> cutdcotdist = cuts_.getParameter<std::vector<double> >("cutdcotdist");
310 
311  float dist = (electron->convDist() == -9999.? 9999:electron->convDist());
312  float dcot = (electron->convDcot() == -9999.? 9999:electron->convDcot());
313 
314  float dcotdistcomb = ((0.04 - std::max(dist, dcot)) > 0?(0.04 - std::max(dist, dcot)):0);
315 
316  if ((mishits < cutmishits[cat+bin*9]) and
317  (dcotdistcomb < cutdcotdist[cat+bin*9]))
318  result += 4;
319 
320  return result;
321  }
322 
323  if (type_ == "classbased" && (version_ == "V06" || version_ == "")) {
324  std::vector<double> cutIsoSum = cuts_.getParameter<std::vector<double> >("cutiso_sum");
325  std::vector<double> cutIsoSumCorr = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
326  std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
327  std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
328  std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
329  std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
330  std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
331  std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutfmishits");
332  std::vector<double> cutdcotdist = cuts_.getParameter<std::vector<double> >("cutdcotdist");
333  std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip_gsf");
334  std::vector<double> cutIsoSumCorrl = cuts_.getParameter<std::vector<double> >("cutiso_sumoetl");
335  std::vector<double> cuthoel = cuts_.getParameter<std::vector<double> >("cuthoel");
336  std::vector<double> cutseel = cuts_.getParameter<std::vector<double> >("cutseel");
337  std::vector<double> cutdphil = cuts_.getParameter<std::vector<double> >("cutdphiinl");
338  std::vector<double> cutdetal = cuts_.getParameter<std::vector<double> >("cutdetainl");
339  std::vector<double> cutipl = cuts_.getParameter<std::vector<double> >("cutip_gsfl");
340 
341  int result = 0;
342 
343  const int ncuts = 10;
344  std::vector<bool> cut_results(ncuts, false);
345 
346  float iso_sum = tkIso + ecalIso + hcalIso;
347  float scEta = electron->superCluster()->eta();
348  if(fabs(scEta)>1.5)
349  iso_sum += (fabs(scEta)-1.5)*1.09;
350 
351  float iso_sumoet = iso_sum*(40./scEt);
352 
353  float eseedopincor = eSeedOverPin + fBrem;
354  if(fBrem < 0)
355  eseedopincor = eSeedOverPin;
356 
357  float dist = (electron->convDist() == -9999.? 9999:electron->convDist());
358  float dcot = (electron->convDcot() == -9999.? 9999:electron->convDcot());
359 
360  float dcotdistcomb = ((0.04 - std::max(dist, dcot)) > 0?(0.04 - std::max(dist, dcot)):0);
361 
362  for (int cut=0; cut<ncuts; cut++) {
363  switch (cut) {
364  case 0:
365  cut_results[cut] = compute_cut(fabs(deltaEtaIn), scEt, cutdetal[cat], cutdeta[cat]);
366  break;
367  case 1:
368  cut_results[cut] = compute_cut(fabs(deltaPhiIn), scEt, cutdphil[cat], cutdphi[cat]);
369  break;
370  case 2:
371  cut_results[cut] = (eseedopincor > cuteopin[cat]);
372  break;
373  case 3:
374  cut_results[cut] = compute_cut(hOverE, scEt, cuthoel[cat], cuthoe[cat]);
375  break;
376  case 4:
377  cut_results[cut] = compute_cut(sigmaee, scEt, cutseel[cat], cutsee[cat]);
378  break;
379  case 5:
380  cut_results[cut] = compute_cut(iso_sumoet, scEt, cutIsoSumCorrl[cat], cutIsoSumCorr[cat]);
381  break;
382  case 6:
383  cut_results[cut] = (iso_sum < cutIsoSum[cat]);
384  break;
385  case 7:
386  cut_results[cut] = compute_cut(fabs(ip), scEt, cutipl[cat], cutip[cat]);
387  break;
388  case 8:
389  cut_results[cut] = (mishits < cutmishits[cat]);
390  break;
391  case 9:
392  cut_results[cut] = (dcotdistcomb < cutdcotdist[cat]);
393  break;
394  }
395  }
396 
397  // ID part
398  if (cut_results[0] & cut_results[1] & cut_results[2] & cut_results[3] & cut_results[4])
399  result = result + 1;
400 
401  // ISO part
402  if (cut_results[5] & cut_results[6])
403  result = result + 2;
404 
405  // IP part
406  if (cut_results[7])
407  result = result + 8;
408 
409  // Conversion part
410  if (cut_results[8] and cut_results[9])
411  result = result + 4;
412 
413  return result;
414  }
415 
416  return -1.;
417 }
T getParameter(std::string const &) const
float dr04HcalTowerSumEt() const
Definition: GsfElectron.h:434
float eSuperClusterOverP() const
Definition: GsfElectron.h:201
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
edm::ParameterSet cuts_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
float fbrem() const
Definition: GsfElectron.h:528
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:166
float convDist() const
Definition: GsfElectron.h:478
double result(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
T eta() const
bool isEE() const
Definition: GsfElectron.h:335
bool isEB() const
Definition: GsfElectron.h:334
return((rh^lh)&mask)
float convDcot() const
Definition: GsfElectron.h:479
int classify(const reco::GsfElectron *)
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:205
float sigmaIetaIeta() const
Definition: GsfElectron.h:378
float hadronicOverEm() const
Definition: GsfElectron.h:393
float eSeedClusterOverP() const
Definition: GsfElectron.h:202
float dr04EcalRecHitSumEt() const
Definition: GsfElectron.h:431
const T & max(const T &a, const T &b)
edm::InputTag verticesCollection_
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:208
float dr03TkSumPt() const
Definition: GsfElectron.h:422
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
tuple cut
Definition: align_tpl.py:88
GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:167
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
bool compute_cut(double x, double et, double cut_min, double cut_max, bool gtn=false)
perl if(1 lt scalar(@::datatypes))
Definition: edlooper.cc:31
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float sigmaEtaEta() const
Definition: GsfElectron.h:377
int CutBasedElectronID::classify ( const reco::GsfElectron electron)

Definition at line 49 of file CutBasedElectronID.cc.

References reco::GsfElectron::ecalDrivenSeed(), reco::GsfElectron::eSuperClusterOverP(), eta(), reco::GsfElectron::fbrem(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), newCategories_, or, reco::GsfElectron::superCluster(), reco::GsfElectron::trackerDrivenSeed(), and version_.

Referenced by cicSelection().

49  {
50 
51  double eta = fabs(electron->superCluster()->eta());
52  double eOverP = electron->eSuperClusterOverP();
53  double fBrem = electron->fbrem();
54 
55  int cat = -1;
56  if (version_ == "V00" || version_ == "V01") {
57  if((electron->isEB() && fBrem<0.06) || (electron->isEE() && fBrem<0.1))
58  cat=1;
59  else if (eOverP < 1.2 && eOverP > 0.8)
60  cat=0;
61  else
62  cat=2;
63 
64  return cat;
65 
66  } else if (version_ == "V02") {
67  if (electron->isEB()) { // BARREL
68  if(fBrem < 0.12)
69  cat=1;
70  else if (eOverP < 1.2 && eOverP > 0.9)
71  cat=0;
72  else
73  cat=2;
74  } else { // ENDCAP
75  if(fBrem < 0.2)
76  cat=1;
77  else if (eOverP < 1.22 && eOverP > 0.82)
78  cat=0;
79  else
80  cat=2;
81  }
82 
83  return cat;
84 
85  } else {
86  if (electron->isEB()) {
87  if ((fBrem >= 0.12) and (eOverP > 0.9) and (eOverP < 1.2))
88  cat = 0;
89  else if (((eta > .445 and eta < .45 ) or
90  (eta > .79 and eta < .81 ) or
91  (eta > 1.137 and eta < 1.157 ) or
92  (eta > 1.47285 and eta < 1.4744)) and newCategories_)
93  cat = 6;
94  else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
95  cat = 8;
96  else if (fBrem < 0.12)
97  cat = 1;
98  else
99  cat = 2;
100  } else {
101  if ((fBrem >= 0.2) and (eOverP > 0.82) and (eOverP < 1.22))
102  cat = 3;
103  else if (eta > 1.5 and eta < 1.58 and newCategories_)
104  cat = 7;
105  else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
106  cat = 8;
107  else if (fBrem < 0.2)
108  cat = 4;
109  else
110  cat = 5;
111  }
112 
113  return cat;
114  }
115 
116  return -1;
117 }
float eSuperClusterOverP() const
Definition: GsfElectron.h:201
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
float fbrem() const
Definition: GsfElectron.h:528
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:166
T eta() const
bool isEE() const
Definition: GsfElectron.h:335
bool isEB() const
Definition: GsfElectron.h:334
bool trackerDrivenSeed() const
Definition: GsfElectron.h:169
perl if(1 lt scalar(@::datatypes))
Definition: edlooper.cc:31
bool ecalDrivenSeed() const
Definition: GsfElectron.h:168
bool CutBasedElectronID::compute_cut ( double  x,
double  et,
double  cut_min,
double  cut_max,
bool  gtn = false 
)

Definition at line 420 of file CutBasedElectronID.cc.

References accept(), align_tpl::cut, and ExpressReco_HICollisions_FallBack::et.

Referenced by cicSelection().

420  {
421 
422  float et_min = 10;
423  float et_max = 40;
424 
425  bool accept = false;
426  float cut = cut_max; // the cut at et=40 GeV
427 
428  if(et < et_max) {
429  cut = cut_min + (1/et_min - 1/et)*(cut_max - cut_min)/(1/et_min - 1/et_max);
430  }
431 
432  if(et < et_min) {
433  cut = cut_min;
434  }
435 
436  if(gtn) { // useful for e/p cut which is gt
437  accept = (x >= cut);
438  }
439  else {
440  accept = (x <= cut);
441  }
442 
443  //std::cout << x << " " << cut_min << " " << cut << " " << cut_max << " " << et << " " << accept << std::endl;
444  return accept;
445 }
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:21
tuple cut
Definition: align_tpl.py:88
double CutBasedElectronID::result ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
)
virtual

Reimplemented from ElectronIDAlgo.

Definition at line 36 of file CutBasedElectronID.cc.

References cicSelection(), robustSelection(), and type_.

Referenced by cicSelection(), and robustSelection().

38  {
39 
40  if (type_ == "classbased")
41  return cicSelection(electron, e, es);
42  else if (type_ == "robust")
43  return robustSelection(electron, e, es);
44 
45  return 0;
46 
47 }
double robustSelection(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
double cicSelection(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
double CutBasedElectronID::robustSelection ( const reco::GsfElectron electron,
const edm::Event e,
const edm::EventSetup es 
)

Definition at line 447 of file CutBasedElectronID.cc.

References reco::GsfElectron::convDcot(), reco::GsfElectron::convDist(), align_tpl::cut, cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::dr03EcalRecHitSumEt(), reco::GsfElectron::dr03HcalDepth1TowerSumEt(), reco::GsfElectron::dr03HcalDepth2TowerSumEt(), reco::GsfElectron::dr03HcalTowerSumEt(), reco::GsfElectron::dr03TkSumPt(), reco::GsfElectron::dr04EcalRecHitSumEt(), reco::GsfElectron::dr04HcalDepth1TowerSumEt(), reco::GsfElectron::dr04HcalDepth2TowerSumEt(), reco::GsfElectron::dr04HcalTowerSumEt(), reco::GsfElectron::e1x5(), reco::GsfElectron::e2x5Max(), reco::GsfElectron::e5x5(), reco::GsfElectron::eSuperClusterOverP(), eta(), funct::exp(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), edm::HandleBase::isValid(), max(), or, reco::LeafCandidate::p4(), reco::BeamSpot::position(), edm::Handle< T >::product(), quality_, result(), reco::GsfElectron::sigmaEtaEta(), reco::GsfElectron::sigmaIetaIeta(), funct::sin(), reco::GsfElectron::superCluster(), type_, version_, and verticesCollection_.

Referenced by result().

449  {
450 
451  double scTheta = (2*atan(exp(-electron->superCluster()->eta())));
452  double scEt = electron->superCluster()->energy()*sin(scTheta);
453  double eta = electron->p4().Eta();
454  double eOverP = electron->eSuperClusterOverP();
455  double hOverE = electron->hadronicOverEm();
456  double sigmaee = electron->sigmaIetaIeta();
457  double e25Max = electron->e2x5Max();
458  double e15 = electron->e1x5();
459  double e55 = electron->e5x5();
460  double e25Maxoe55 = e25Max/e55;
461  double e15oe55 = e15/e55 ;
462  double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
463  double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
464 
465  double ip = 0;
466  int mishits = electron->gsfTrack()->trackerExpectedHitsInner().numberOfHits();
467  double tkIso = electron->dr03TkSumPt();
468  double ecalIso = electron->dr04EcalRecHitSumEt();
469  double ecalIsoPed = (electron->isEB())?std::max(0.,ecalIso-1.):ecalIso;
470  double hcalIso = electron->dr04HcalTowerSumEt();
471  double hcalIso1 = electron->dr04HcalDepth1TowerSumEt();
472  double hcalIso2 = electron->dr04HcalDepth2TowerSumEt();
473 
474  if (version_ == "V00") {
475  sigmaee = electron->sigmaEtaEta();
476  if (electron->isEE())
477  sigmaee = sigmaee - 0.02*(fabs(eta) - 2.3); //correct sigmaetaeta dependence on eta in endcap
478  }
479 
480  if (version_ == "V03" or version_ == "V04") {
481  edm::Handle<reco::BeamSpot> pBeamSpot;
482  // uses the same name for the vertex collection to avoid adding more new names
483  e.getByLabel(verticesCollection_, pBeamSpot);
484  if (pBeamSpot.isValid()) {
485  const reco::BeamSpot *bspot = pBeamSpot.product();
486  const math::XYZPoint bspotPosition = bspot->position();
487  ip = fabs(electron->gsfTrack()->dxy(bspotPosition));
488  } else
489  ip = fabs(electron->gsfTrack()->dxy());
490  }
491 
492  if (version_ == "V04" or version_ == "V05") {
493  ecalIso = electron->dr03EcalRecHitSumEt();
494  ecalIsoPed = (electron->isEB())?std::max(0.,ecalIso-1.):ecalIso;
495  hcalIso = electron->dr03HcalTowerSumEt();
496  hcalIso1 = electron->dr03HcalDepth1TowerSumEt();
497  hcalIso2 = electron->dr03HcalDepth2TowerSumEt();
498  }
499 
500  if (version_ == "V05") {
503  if (vtxH->size() != 0) {
504  reco::VertexRef vtx(vtxH, 0);
505  ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(),vtx->y(),vtx->z())));
506  } else
507  ip = fabs(electron->gsfTrack()->dxy());
508  }
509 
510  // .....................................................................................
511  std::vector<double> cut;
512  // ROBUST Selection
513  if (type_ == "robust") {
514 
515  double result = 0;
516 
517  // hoe, sigmaEtaEta, dPhiIn, dEtaIn
518  if (electron->isEB())
519  cut = cuts_.getParameter<std::vector<double> >("barrel");
520  else
521  cut = cuts_.getParameter<std::vector<double> >("endcap");
522  // check isolations: if only isolation passes result = 2
523  if (quality_ == "highenergy") {
524  if ((tkIso > cut[6] || hcalIso2 > cut[11]) ||
525  (electron->isEB() && ((ecalIso + hcalIso1) > cut[7]+cut[8]*scEt)) ||
526  (electron->isEE() && (scEt >= 50.) && ((ecalIso + hcalIso1) > cut[7]+cut[8]*(scEt-50))) ||
527  (electron->isEE() && (scEt < 50.) && ((ecalIso + hcalIso1) > cut[9]+cut[10]*(scEt-50))))
528  result = 0;
529  else
530  result = 2;
531  } else {
532  if ((tkIso > cut[6]) || (ecalIso > cut[7]) || (hcalIso > cut[8]) || (hcalIso1 > cut[9]) || (hcalIso2 > cut[10]) ||
533  (tkIso/electron->p4().Pt() > cut[11]) || (ecalIso/electron->p4().Pt() > cut[12]) || (hcalIso/electron->p4().Pt() > cut[13]) ||
534  ((tkIso+ecalIso+hcalIso)>cut[14]) || (((tkIso+ecalIso+hcalIso)/ electron->p4().Pt()) > cut[15]) ||
535  ((tkIso+ecalIsoPed+hcalIso)>cut[16]) || (((tkIso+ecalIsoPed+hcalIso)/ electron->p4().Pt()) > cut[17]) )
536  result = 0.;
537  else
538  result = 2.;
539  }
540 
541  if ((hOverE < cut[0]) && (sigmaee < cut[1]) && (fabs(deltaPhiIn) < cut[2]) &&
542  (fabs(deltaEtaIn) < cut[3]) && (e25Maxoe55 > cut[4] && e15oe55 > cut[5]) &&
543  (sigmaee >= cut[18]) && (eOverP > cut[19] && eOverP < cut[20]) )
544  { result = result + 1 ; }
545 
546  if (ip > cut[21])
547  return result;
548  if (mishits > cut[22]) // expected missing hits
549  return result;
550  // positive cut[23] means to demand a valid hit in 1st layer PXB
551  if (cut[23] >0 && not (electron->gsfTrack()->hitPattern().hasValidHitInFirstPixelBarrel()))
552  return result;
553 
554  // cut[24]: Dist cut[25]: dcot
555  float dist = fabs(electron->convDist());
556  float dcot = fabs(electron->convDcot());
557  bool isConversion = (cut[24]>99. || cut[25]>99.)?false:(dist < cut[24] && dcot < cut[25]);
558  if (isConversion)
559  return result ;
560 
561  result += 4 ;
562 
563  return result ;
564  }
565 
566  return -1. ;
567  }
T getParameter(std::string const &) const
float dr04HcalTowerSumEt() const
Definition: GsfElectron.h:434
float eSuperClusterOverP() const
Definition: GsfElectron.h:201
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
float dr04HcalDepth1TowerSumEt() const
Definition: GsfElectron.h:432
edm::ParameterSet cuts_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
float dr03HcalDepth2TowerSumEt() const
Definition: GsfElectron.h:425
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:166
float e2x5Max() const
Definition: GsfElectron.h:380
float convDist() const
Definition: GsfElectron.h:478
double result(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
T eta() const
bool isEE() const
Definition: GsfElectron.h:335
bool isEB() const
Definition: GsfElectron.h:334
float convDcot() const
Definition: GsfElectron.h:479
float dr04HcalDepth2TowerSumEt() const
Definition: GsfElectron.h:433
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:205
float sigmaIetaIeta() const
Definition: GsfElectron.h:378
float hadronicOverEm() const
Definition: GsfElectron.h:393
float dr04EcalRecHitSumEt() const
Definition: GsfElectron.h:431
const T & max(const T &a, const T &b)
edm::InputTag verticesCollection_
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:208
float dr03TkSumPt() const
Definition: GsfElectron.h:422
float dr03HcalDepth1TowerSumEt() const
Definition: GsfElectron.h:424
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
tuple cut
Definition: align_tpl.py:88
GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:167
float e1x5() const
Definition: GsfElectron.h:379
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
T const * product() const
Definition: Handle.h:74
float dr03EcalRecHitSumEt() const
Definition: GsfElectron.h:423
float e5x5() const
Definition: GsfElectron.h:381
float dr03HcalTowerSumEt() const
Definition: GsfElectron.h:426
const Point & position() const
position
Definition: BeamSpot.h:63
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
float sigmaEtaEta() const
Definition: GsfElectron.h:377
void CutBasedElectronID::setup ( const edm::ParameterSet conf)
virtual

Reimplemented from ElectronIDAlgo.

Definition at line 11 of file CutBasedElectronID.cc.

References ElectronIDAlgo::baseSetup(), cuts_, edm::hlt::Exception, edm::ParameterSet::getParameter(), newCategories_, or, quality_, type_, version_, verticesCollection_, and wantBinning_.

11  {
12 
13  // Get all the parameters
14  baseSetup(conf);
15 
16  type_ = conf.getParameter<std::string>("electronIDType");
17  quality_ = conf.getParameter<std::string>("electronQuality");
18  version_ = conf.getParameter<std::string>("electronVersion");
19  verticesCollection_ = conf.getParameter<edm::InputTag>("verticesCollection");
20 
21  if (type_ == "classbased" and (version_ == "V03" or version_ == "V04" or version_ == "V05" or version_ == "")) {
22  wantBinning_ = conf.getParameter<bool>("etBinning");
23  newCategories_ = conf.getParameter<bool>("additionalCategories");
24  }
25 
26  if (type_ == "robust" || type_ == "classbased") {
27  std::string stringCut = type_+quality_+"EleIDCuts"+version_;
28  cuts_ = conf.getParameter<edm::ParameterSet>(stringCut);
29  }
30  else {
31  throw cms::Exception("Configuration")
32  << "Invalid electronType parameter in CutBasedElectronID: must be robust or classbased\n";
33  }
34 }
T getParameter(std::string const &) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
edm::ParameterSet cuts_
void baseSetup(const edm::ParameterSet &conf)
edm::InputTag verticesCollection_

Member Data Documentation

edm::ParameterSet CutBasedElectronID::cuts_
private

Definition at line 28 of file CutBasedElectronID.h.

Referenced by cicSelection(), robustSelection(), and setup().

bool CutBasedElectronID::newCategories_
private

Definition at line 23 of file CutBasedElectronID.h.

Referenced by classify(), and setup().

std::string CutBasedElectronID::quality_
private

Definition at line 25 of file CutBasedElectronID.h.

Referenced by cicSelection(), robustSelection(), and setup().

std::string CutBasedElectronID::type_
private
std::string CutBasedElectronID::version_
private

Definition at line 26 of file CutBasedElectronID.h.

Referenced by cicSelection(), classify(), robustSelection(), and setup().

edm::InputTag CutBasedElectronID::verticesCollection_
private

Definition at line 27 of file CutBasedElectronID.h.

Referenced by cicSelection(), robustSelection(), and setup().

bool CutBasedElectronID::wantBinning_
private

Definition at line 22 of file CutBasedElectronID.h.

Referenced by cicSelection(), and setup().