CMS 3D CMS Logo

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 (const edm::ParameterSet &conf, edm::ConsumesCollector &iC)
 
double result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &) override
 
double robustSelection (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
 
void setup (const edm::ParameterSet &conf) override
 
 ~CutBasedElectronID () override
 
- Public Member Functions inherited from ElectronIDAlgo
 ElectronIDAlgo ()
 
virtual ~ElectronIDAlgo ()
 

Private Attributes

edm::ParameterSet cuts_
 
bool newCategories_
 
std::string quality_
 
std::string type_
 
std::string version_
 
edm::EDGetTokenT< std::vector< reco::Vertex > > verticesCollection_
 
bool wantBinning_
 

Additional Inherited Members

- Protected Attributes inherited from ElectronIDAlgo
edm::InputTag reducedBarrelRecHitCollection_
 
edm::InputTag reducedEndcapRecHitCollection_
 

Detailed Description

Definition at line 11 of file CutBasedElectronID.h.

Constructor & Destructor Documentation

◆ CutBasedElectronID()

CutBasedElectronID::CutBasedElectronID ( const edm::ParameterSet conf,
edm::ConsumesCollector iC 
)

Definition at line 9 of file CutBasedElectronID.cc.

References edm::ConsumesCollector::consumes(), edm::ParameterSet::getParameter(), and verticesCollection_.

9  {
10  verticesCollection_ = iC.consumes<std::vector<reco::Vertex> >(conf.getParameter<edm::InputTag>("verticesCollection"));
11 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< std::vector< reco::Vertex > > verticesCollection_

◆ ~CutBasedElectronID()

CutBasedElectronID::~CutBasedElectronID ( )
inlineoverride

Definition at line 15 of file CutBasedElectronID.h.

15 {};

Member Function Documentation

◆ cicSelection()

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

Definition at line 117 of file CutBasedElectronID.cc.

References newFWLiteAna::bin, eostools::cat(), classify(), compute_cut(), DMR_cfg::cut, electronIdCutBasedExt_cfi::cutdcotdist, electronIdCutBasedExt_cfi::cutdeta, electronIdCutBasedExt_cfi::cutdphi, electronIdCutBasedExt_cfi::cuteopin, electronIdCutBasedExt_cfi::cutet, electronIdCutBasedExt_cfi::cuthoe, electronIdCutBasedExt_cfi::cuthoel, electronIdCutBasedExt_cfi::cutip, electronIdCutBasedExt_cfi::cutmishits, cuts_, electronIdCutBasedExt_cfi::cutsee, electronIdCutBasedExt_cfi::cutseel, ticlTrackstersMerge_cfi::cutTk, electronIdCutBasedExt_cfi::deltaEtaIn, electronIdCutBasedExt_cfi::deltaPhiIn, MillePedeFileConverter_cfg::e, run3scouting_cff::ecalIso, HPSPFTauProducerPuppi_cfi::electron, EgHLTOffHistBins_cfi::eOverP, electronIdCutBasedExt_cfi::eSeedOverPin, PVValHelper::eta, JetChargeProducer_cfi::exp, edm::ParameterSet::getParameter(), run3scouting_cff::hcalIso, EgHLTOffHistBins_cfi::hOverE, SiStripPI::max, or, funct::pow(), quality_, result(), funct::sin(), type_, version_, verticesCollection_, L1BJetProducer_cff::vtx, and wantBinning_.

Referenced by result().

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

◆ classify()

int CutBasedElectronID::classify ( const reco::GsfElectron electron)

Definition at line 49 of file CutBasedElectronID.cc.

References eostools::cat(), HPSPFTauProducerPuppi_cfi::electron, EgHLTOffHistBins_cfi::eOverP, PVValHelper::eta, newCategories_, or, and version_.

Referenced by cicSelection().

49  {
50  double eta = fabs(electron->superCluster()->eta());
51  double eOverP = electron->eSuperClusterOverP();
52  double fBrem = electron->fbrem();
53 
54  int cat = -1;
55  if (version_ == "V00" || version_ == "V01") {
56  if ((electron->isEB() && fBrem < 0.06) || (electron->isEE() && fBrem < 0.1))
57  cat = 1;
58  else if (eOverP < 1.2 && eOverP > 0.8)
59  cat = 0;
60  else
61  cat = 2;
62 
63  return cat;
64 
65  } else if (version_ == "V02") {
66  if (electron->isEB()) { // BARREL
67  if (fBrem < 0.12)
68  cat = 1;
69  else if (eOverP < 1.2 && eOverP > 0.9)
70  cat = 0;
71  else
72  cat = 2;
73  } else { // ENDCAP
74  if (fBrem < 0.2)
75  cat = 1;
76  else if (eOverP < 1.22 && eOverP > 0.82)
77  cat = 0;
78  else
79  cat = 2;
80  }
81 
82  return cat;
83 
84  } else {
85  if (electron->isEB()) {
86  if ((fBrem >= 0.12) and (eOverP > 0.9) and (eOverP < 1.2))
87  cat = 0;
88  else if (((eta > .445 and eta < .45) or (eta > .79 and eta < .81) or (eta > 1.137 and eta < 1.157) or
89  (eta > 1.47285 and eta < 1.4744)) and
91  cat = 6;
92  else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
93  cat = 8;
94  else if (fBrem < 0.12)
95  cat = 1;
96  else
97  cat = 2;
98  } else {
99  if ((fBrem >= 0.2) and (eOverP > 0.82) and (eOverP < 1.22))
100  cat = 3;
101  else if (eta > 1.5 and eta < 1.58 and newCategories_)
102  cat = 7;
103  else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
104  cat = 8;
105  else if (fBrem < 0.2)
106  cat = 4;
107  else
108  cat = 5;
109  }
110 
111  return cat;
112  }
113 
114  return -1;
115 }
def cat(path)
Definition: eostools.py:401
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::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
if(threadIdxLocalY==0 &&threadIdxLocalX==0)

◆ compute_cut()

bool CutBasedElectronID::compute_cut ( double  x,
double  et,
double  cut_min,
double  cut_max,
bool  gtn = false 
)

Definition at line 409 of file CutBasedElectronID.cc.

References accept(), DMR_cfg::cut, EgHLTOffHistBins_cfi::et, and x.

Referenced by cicSelection().

409  {
410  float et_min = 10;
411  float et_max = 40;
412 
413  bool accept = false;
414  float cut = cut_max; // the cut at et=40 GeV
415 
416  if (et < et_max) {
417  cut = cut_min + (1 / et_min - 1 / et) * (cut_max - cut_min) / (1 / et_min - 1 / et_max);
418  }
419 
420  if (et < et_min) {
421  cut = cut_min;
422  }
423 
424  if (gtn) { // useful for e/p cut which is gt
425  accept = (x >= cut);
426  } else {
427  accept = (x <= cut);
428  }
429 
430  //std::cout << x << " " << cut_min << " " << cut << " " << cut_max << " " << et << " " << accept << std::endl;
431  return accept;
432 }
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31

◆ result()

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

Reimplemented from ElectronIDAlgo.

Definition at line 40 of file CutBasedElectronID.cc.

References cicSelection(), MillePedeFileConverter_cfg::e, HPSPFTauProducerPuppi_cfi::electron, robustSelection(), and type_.

Referenced by cicSelection(), and robustSelection().

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

◆ robustSelection()

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

Definition at line 434 of file CutBasedElectronID.cc.

References DMR_cfg::cut, cuts_, electronIdCutBasedExt_cfi::deltaEtaIn, electronIdCutBasedExt_cfi::deltaPhiIn, MillePedeFileConverter_cfg::e, run3scouting_cff::ecalIso, HPSPFTauProducerPuppi_cfi::electron, EgHLTOffHistBins_cfi::eOverP, PVValHelper::eta, JetChargeProducer_cfi::exp, funct::false, edm::ParameterSet::getParameter(), run3scouting_cff::hcalIso, EgHLTOffHistBins_cfi::hOverE, edm::HandleBase::isValid(), SiStripPI::max, or, GeomDetEnumerators::PixelBarrel, reco::BeamSpot::position(), edm::Handle< T >::product(), quality_, result(), funct::sin(), type_, version_, verticesCollection_, and L1BJetProducer_cff::vtx.

Referenced by result().

436  {
437  double scTheta = (2 * atan(exp(-electron->superCluster()->eta())));
438  double scEt = electron->superCluster()->energy() * sin(scTheta);
439  double eta = electron->p4().Eta();
440  double eOverP = electron->eSuperClusterOverP();
441  double hOverE = electron->hadronicOverEm();
442  double sigmaee = electron->sigmaIetaIeta();
443  double e25Max = electron->e2x5Max();
444  double e15 = electron->e1x5();
445  double e55 = electron->e5x5();
446  double e25Maxoe55 = e25Max / e55;
447  double e15oe55 = e15 / e55;
448  double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
449  double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
450 
451  double ip = 0;
452  int mishits = electron->gsfTrack()->missingInnerHits();
453  double tkIso = electron->dr03TkSumPt();
454  double ecalIso = electron->dr04EcalRecHitSumEt();
455  double ecalIsoPed = (electron->isEB()) ? std::max(0., ecalIso - 1.) : ecalIso;
456  double hcalIso = electron->dr04HcalTowerSumEt();
457  double hcalIso1 = electron->dr04HcalTowerSumEt(1);
458  double hcalIso2 = electron->dr04HcalTowerSumEt(2);
459 
460  if (version_ == "V00") {
461  sigmaee = electron->sigmaEtaEta();
462  if (electron->isEE())
463  sigmaee = sigmaee - 0.02 * (fabs(eta) - 2.3); //correct sigmaetaeta dependence on eta in endcap
464  }
465 
466  if (version_ == "V03" or version_ == "V04") {
467  edm::Handle<reco::BeamSpot> pBeamSpot;
468  // uses the same name for the vertex collection to avoid adding more new names
469  e.getByToken(verticesCollection_, pBeamSpot);
470  if (pBeamSpot.isValid()) {
471  const reco::BeamSpot* bspot = pBeamSpot.product();
472  const math::XYZPoint& bspotPosition = bspot->position();
473  ip = fabs(electron->gsfTrack()->dxy(bspotPosition));
474  } else
475  ip = fabs(electron->gsfTrack()->dxy());
476  }
477 
478  if (version_ == "V04" or version_ == "V05") {
479  ecalIso = electron->dr03EcalRecHitSumEt();
480  ecalIsoPed = (electron->isEB()) ? std::max(0., ecalIso - 1.) : ecalIso;
481  hcalIso = electron->dr03HcalTowerSumEt();
482  hcalIso1 = electron->dr03HcalTowerSumEt(1);
483  hcalIso2 = electron->dr03HcalTowerSumEt(2);
484  }
485 
486  if (version_ == "V05") {
488  e.getByToken(verticesCollection_, vtxH);
489  if (!vtxH->empty()) {
490  reco::VertexRef vtx(vtxH, 0);
491  ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(), vtx->y(), vtx->z())));
492  } else
493  ip = fabs(electron->gsfTrack()->dxy());
494  }
495 
496  // .....................................................................................
497  std::vector<double> cut;
498  // ROBUST Selection
499  if (type_ == "robust") {
500  double result = 0;
501 
502  // hoe, sigmaEtaEta, dPhiIn, dEtaIn
503  if (electron->isEB())
504  cut = cuts_.getParameter<std::vector<double> >("barrel");
505  else
506  cut = cuts_.getParameter<std::vector<double> >("endcap");
507  // check isolations: if only isolation passes result = 2
508  if (quality_ == "highenergy") {
509  if ((tkIso > cut[6] || hcalIso2 > cut[11]) ||
510  (electron->isEB() && ((ecalIso + hcalIso1) > cut[7] + cut[8] * scEt)) ||
511  (electron->isEE() && (scEt >= 50.) && ((ecalIso + hcalIso1) > cut[7] + cut[8] * (scEt - 50))) ||
512  (electron->isEE() && (scEt < 50.) && ((ecalIso + hcalIso1) > cut[9] + cut[10] * (scEt - 50))))
513  result = 0;
514  else
515  result = 2;
516  } else {
517  if ((tkIso > cut[6]) || (ecalIso > cut[7]) || (hcalIso > cut[8]) || (hcalIso1 > cut[9]) || (hcalIso2 > cut[10]) ||
518  (tkIso / electron->p4().Pt() > cut[11]) || (ecalIso / electron->p4().Pt() > cut[12]) ||
519  (hcalIso / electron->p4().Pt() > cut[13]) || ((tkIso + ecalIso + hcalIso) > cut[14]) ||
520  (((tkIso + ecalIso + hcalIso) / electron->p4().Pt()) > cut[15]) ||
521  ((tkIso + ecalIsoPed + hcalIso) > cut[16]) ||
522  (((tkIso + ecalIsoPed + hcalIso) / electron->p4().Pt()) > cut[17]))
523  result = 0.;
524  else
525  result = 2.;
526  }
527 
528  if ((hOverE < cut[0]) && (sigmaee < cut[1]) && (fabs(deltaPhiIn) < cut[2]) && (fabs(deltaEtaIn) < cut[3]) &&
529  (e25Maxoe55 > cut[4] && e15oe55 > cut[5]) && (sigmaee >= cut[18]) && (eOverP > cut[19] && eOverP < cut[20])) {
530  result = result + 1;
531  }
532 
533  if (ip > cut[21])
534  return result;
535  if (mishits > cut[22]) // expected missing hits
536  return result;
537  // positive cut[23] means to demand a valid hit in 1st layer PXB
538  if (cut[23] > 0 &&
539  !electron->gsfTrack()->hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1))
540  return result;
541 
542  // cut[24]: Dist cut[25]: dcot
543  float dist = fabs(electron->convDist());
544  float dcot = fabs(electron->convDcot());
545  bool isConversion = (cut[24] > 99. || cut[25] > 99.) ? false : (dist < cut[24] && dcot < cut[25]);
546  if (isConversion)
547  return result;
548 
549  result += 4;
550 
551  return result;
552  }
553 
554  return -1.;
555 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const Point & position() const
position
Definition: BeamSpot.h:59
double result(const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &) override
edm::ParameterSet cuts_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< std::vector< reco::Vertex > > verticesCollection_
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::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isValid() const
Definition: HandleBase.h:70

◆ setup()

void CutBasedElectronID::setup ( const edm::ParameterSet conf)
overridevirtual

Reimplemented from ElectronIDAlgo.

Definition at line 13 of file CutBasedElectronID.cc.

References cuts_, Exception, edm::ParameterSet::getParameter(), newCategories_, or, quality_, AlCaHLTBitMon_QueryRunRegistry::string, type_, version_, and wantBinning_.

13  {
14  // Get all the parameters
15  //baseSetup(conf);
16 
17  type_ = conf.getParameter<std::string>("electronIDType");
18  quality_ = conf.getParameter<std::string>("electronQuality");
19  version_ = conf.getParameter<std::string>("electronVersion");
20  //verticesCollection_ = conf.getParameter<edm::InputTag>("verticesCollection");
21 
22  if (type_ == "classbased" and (version_ == "V06")) {
23  newCategories_ = conf.getParameter<bool>("additionalCategories");
24  }
25 
26  if (type_ == "classbased" and (version_ == "V03" or version_ == "V04" or version_ == "V05" or version_.empty())) {
27  wantBinning_ = conf.getParameter<bool>("etBinning");
28  newCategories_ = conf.getParameter<bool>("additionalCategories");
29  }
30 
31  if (type_ == "robust" || type_ == "classbased") {
32  std::string stringCut = type_ + quality_ + "EleIDCuts" + version_;
33  cuts_ = conf.getParameter<edm::ParameterSet>(stringCut);
34  } else {
35  throw cms::Exception("Configuration")
36  << "Invalid electronType parameter in CutBasedElectronID: must be robust or classbased\n";
37  }
38 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ParameterSet cuts_
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::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12

Member Data Documentation

◆ cuts_

edm::ParameterSet CutBasedElectronID::cuts_
private

Definition at line 31 of file CutBasedElectronID.h.

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

◆ newCategories_

bool CutBasedElectronID::newCategories_
private

Definition at line 26 of file CutBasedElectronID.h.

Referenced by classify(), and setup().

◆ quality_

std::string CutBasedElectronID::quality_
private

Definition at line 28 of file CutBasedElectronID.h.

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

◆ type_

std::string CutBasedElectronID::type_
private

◆ version_

std::string CutBasedElectronID::version_
private

Definition at line 29 of file CutBasedElectronID.h.

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

◆ verticesCollection_

edm::EDGetTokenT<std::vector<reco::Vertex> > CutBasedElectronID::verticesCollection_
private

Definition at line 30 of file CutBasedElectronID.h.

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

◆ wantBinning_

bool CutBasedElectronID::wantBinning_
private

Definition at line 25 of file CutBasedElectronID.h.

Referenced by cicSelection(), and setup().