Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:07

0001 #include "CutBasedElectronID.h"
0002 
0003 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0005 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0006 
0007 #include <algorithm>
0008 
0009 CutBasedElectronID::CutBasedElectronID(const edm::ParameterSet& conf, edm::ConsumesCollector& iC) {
0010   verticesCollection_ = iC.consumes<std::vector<reco::Vertex> >(conf.getParameter<edm::InputTag>("verticesCollection"));
0011 }
0012 
0013 void CutBasedElectronID::setup(const edm::ParameterSet& conf) {
0014   // Get all the parameters
0015   //baseSetup(conf);
0016 
0017   type_ = conf.getParameter<std::string>("electronIDType");
0018   quality_ = conf.getParameter<std::string>("electronQuality");
0019   version_ = conf.getParameter<std::string>("electronVersion");
0020   //verticesCollection_ = conf.getParameter<edm::InputTag>("verticesCollection");
0021 
0022   if (type_ == "classbased" and (version_ == "V06")) {
0023     newCategories_ = conf.getParameter<bool>("additionalCategories");
0024   }
0025 
0026   if (type_ == "classbased" and (version_ == "V03" or version_ == "V04" or version_ == "V05" or version_.empty())) {
0027     wantBinning_ = conf.getParameter<bool>("etBinning");
0028     newCategories_ = conf.getParameter<bool>("additionalCategories");
0029   }
0030 
0031   if (type_ == "robust" || type_ == "classbased") {
0032     std::string stringCut = type_ + quality_ + "EleIDCuts" + version_;
0033     cuts_ = conf.getParameter<edm::ParameterSet>(stringCut);
0034   } else {
0035     throw cms::Exception("Configuration")
0036         << "Invalid electronType parameter in CutBasedElectronID: must be robust or classbased\n";
0037   }
0038 }
0039 
0040 double CutBasedElectronID::result(const reco::GsfElectron* electron, const edm::Event& e, const edm::EventSetup& es) {
0041   if (type_ == "classbased")
0042     return cicSelection(electron, e, es);
0043   else if (type_ == "robust")
0044     return robustSelection(electron, e, es);
0045 
0046   return 0;
0047 }
0048 
0049 int CutBasedElectronID::classify(const reco::GsfElectron* electron) {
0050   double eta = fabs(electron->superCluster()->eta());
0051   double eOverP = electron->eSuperClusterOverP();
0052   double fBrem = electron->fbrem();
0053 
0054   int cat = -1;
0055   if (version_ == "V00" || version_ == "V01") {
0056     if ((electron->isEB() && fBrem < 0.06) || (electron->isEE() && fBrem < 0.1))
0057       cat = 1;
0058     else if (eOverP < 1.2 && eOverP > 0.8)
0059       cat = 0;
0060     else
0061       cat = 2;
0062 
0063     return cat;
0064 
0065   } else if (version_ == "V02") {
0066     if (electron->isEB()) {  // BARREL
0067       if (fBrem < 0.12)
0068         cat = 1;
0069       else if (eOverP < 1.2 && eOverP > 0.9)
0070         cat = 0;
0071       else
0072         cat = 2;
0073     } else {  // ENDCAP
0074       if (fBrem < 0.2)
0075         cat = 1;
0076       else if (eOverP < 1.22 && eOverP > 0.82)
0077         cat = 0;
0078       else
0079         cat = 2;
0080     }
0081 
0082     return cat;
0083 
0084   } else {
0085     if (electron->isEB()) {
0086       if ((fBrem >= 0.12) and (eOverP > 0.9) and (eOverP < 1.2))
0087         cat = 0;
0088       else if (((eta > .445 and eta < .45) or (eta > .79 and eta < .81) or (eta > 1.137 and eta < 1.157) or
0089                 (eta > 1.47285 and eta < 1.4744)) and
0090                newCategories_)
0091         cat = 6;
0092       else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
0093         cat = 8;
0094       else if (fBrem < 0.12)
0095         cat = 1;
0096       else
0097         cat = 2;
0098     } else {
0099       if ((fBrem >= 0.2) and (eOverP > 0.82) and (eOverP < 1.22))
0100         cat = 3;
0101       else if (eta > 1.5 and eta < 1.58 and newCategories_)
0102         cat = 7;
0103       else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
0104         cat = 8;
0105       else if (fBrem < 0.2)
0106         cat = 4;
0107       else
0108         cat = 5;
0109     }
0110 
0111     return cat;
0112   }
0113 
0114   return -1;
0115 }
0116 
0117 double CutBasedElectronID::cicSelection(const reco::GsfElectron* electron,
0118                                         const edm::Event& e,
0119                                         const edm::EventSetup& es) {
0120   double scTheta = (2 * atan(exp(-electron->superCluster()->eta())));
0121   double scEt = electron->superCluster()->energy() * sin(scTheta);
0122 
0123   double eta = fabs(electron->superCluster()->eta());
0124 
0125   double eOverP = electron->eSuperClusterOverP();
0126   double eSeedOverPin = electron->eSeedClusterOverP();
0127   double fBrem = electron->fbrem();
0128   double hOverE = electron->hadronicOverEm();
0129   double sigmaee = electron->sigmaIetaIeta();  //sqrt(vLocCov[0]);
0130   double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
0131   double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
0132 
0133   double ip = 0;
0134   int mishits = electron->gsfTrack()->missingInnerHits();
0135   double tkIso = electron->dr03TkSumPt();
0136   double ecalIso = electron->dr04EcalRecHitSumEt();
0137   double hcalIso = electron->dr04HcalTowerSumEt();
0138 
0139   if (version_ == "V00") {
0140     sigmaee = electron->sigmaEtaEta();  //sqrt(vCov[0]);
0141     if (electron->isEE())
0142       sigmaee = sigmaee - 0.02 * (fabs(eta) - 2.3);  //correct sigmaetaeta dependence on eta in endcap
0143   }
0144 
0145   if (version_ != "V01" or version_ != "V00") {
0146     edm::Handle<reco::VertexCollection> vtxH;
0147     e.getByToken(verticesCollection_, vtxH);
0148     if (!vtxH->empty()) {
0149       reco::VertexRef vtx(vtxH, 0);
0150       ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(), vtx->y(), vtx->z())));
0151     } else
0152       ip = fabs(electron->gsfTrack()->dxy());
0153 
0154     if (electron->isEB()) {
0155       sigmaee = electron->sigmaIetaIeta();  //sqrt(vCov[0]);
0156     }
0157   }
0158 
0159   std::vector<double> cut;
0160 
0161   int cat = classify(electron);
0162   int eb;
0163 
0164   if (electron->isEB())
0165     eb = 0;
0166   else
0167     eb = 1;
0168 
0169   // LOOSE and TIGHT Selections
0170   if (type_ == "classbased" && (version_ == "V01" || version_ == "V00")) {
0171     if ((eOverP < 0.8) && (fBrem < 0.2))
0172       return 0.;
0173 
0174     cut = cuts_.getParameter<std::vector<double> >("hOverE");
0175     if (hOverE > cut[cat + 4 * eb])
0176       return 0.;
0177 
0178     cut = cuts_.getParameter<std::vector<double> >("sigmaEtaEta");
0179     if (sigmaee > cut[cat + 4 * eb])
0180       return 0.;
0181 
0182     cut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
0183     if (eOverP < 1.5) {
0184       if (fabs(deltaPhiIn) > cut[cat + 4 * eb])
0185         return 0.;
0186     } else {
0187       if (fabs(deltaPhiIn) > cut[3 + 4 * eb])
0188         return 0.;
0189     }
0190 
0191     cut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
0192     if (fabs(deltaEtaIn) > cut[cat + 4 * eb])
0193       return 0.;
0194 
0195     cut = cuts_.getParameter<std::vector<double> >("eSeedOverPin");
0196     if (eSeedOverPin < cut[cat + 4 * eb])
0197       return 0.;
0198 
0199     if (quality_ == "tight")
0200       if (eOverP < 0.9 * (1 - fBrem))
0201         return 0.;
0202 
0203     return 1.;
0204   }
0205 
0206   if (type_ == "classbased" and version_ == "V02") {
0207     double result = 0.;
0208 
0209     int bin = 0;
0210 
0211     if (scEt < 20.)
0212       bin = 2;
0213     else if (scEt > 30.)
0214       bin = 0;
0215     else
0216       bin = 1;
0217 
0218     if (fBrem > 0)
0219       eSeedOverPin = eSeedOverPin + fBrem;
0220 
0221     if (bin != 2) {
0222       tkIso = tkIso * pow(40. / scEt, 2);
0223       ecalIso = ecalIso * pow(40. / scEt, 2);
0224       hcalIso = hcalIso * pow(40. / scEt, 2);
0225     }
0226 
0227     std::vector<double> cutTk = cuts_.getParameter<std::vector<double> >("cutisotk");
0228     std::vector<double> cutEcal = cuts_.getParameter<std::vector<double> >("cutisoecal");
0229     std::vector<double> cutHcal = cuts_.getParameter<std::vector<double> >("cutisohcal");
0230     if ((tkIso > cutTk[cat + 3 * eb + bin * 6]) || (ecalIso > cutEcal[cat + 3 * eb + bin * 6]) ||
0231         (hcalIso > cutHcal[cat + 3 * eb + bin * 6]))
0232       result = 0.;
0233     else
0234       result = 2.;
0235 
0236     if (fBrem > -2) {
0237       std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
0238       std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
0239       std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
0240       std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
0241       std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
0242       std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
0243       std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip");
0244       std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutmishits");
0245       if ((hOverE < cuthoe[cat + 3 * eb + bin * 6]) and (sigmaee < cutsee[cat + 3 * eb + bin * 6]) and
0246           (fabs(deltaPhiIn) < cutdphi[cat + 3 * eb + bin * 6]) and
0247           (fabs(deltaEtaIn) < cutdeta[cat + 3 * eb + bin * 6]) and (eSeedOverPin > cuteopin[cat + 3 * eb + bin * 6]) and
0248           (ip < cutip[cat + 3 * eb + bin * 6]) and (mishits < cutmishits[cat + 3 * eb + bin * 6]))
0249         result = result + 1.;
0250     }
0251     return result;
0252   }
0253 
0254   if (version_ == "V03" or version_ == "V04" or version_ == "V05") {
0255     double result = 0.;
0256 
0257     int bin = 0;
0258 
0259     if (wantBinning_) {
0260       if (scEt < 20.)
0261         bin = 2;
0262       else if (scEt > 30.)
0263         bin = 0;
0264       else
0265         bin = 1;
0266     }
0267 
0268     if (fBrem > 0)
0269       eSeedOverPin = eSeedOverPin + fBrem;
0270 
0271     float iso_sum = tkIso + ecalIso + hcalIso;
0272     float iso_sum_corrected = iso_sum * pow(40. / scEt, 2);
0273 
0274     std::vector<double> cutIsoSum = cuts_.getParameter<std::vector<double> >("cutiso_sum");
0275     std::vector<double> cutIsoSumCorr = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
0276     if ((iso_sum < cutIsoSum[cat + bin * 9]) and (iso_sum_corrected < cutIsoSumCorr[cat + bin * 9]))
0277       result += 2.;
0278 
0279     if (fBrem > -2) {
0280       std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
0281       std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
0282       std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
0283       std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
0284       std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
0285       std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
0286 
0287       if ((hOverE < cuthoe[cat + bin * 9]) and (sigmaee < cutsee[cat + bin * 9]) and
0288           (fabs(deltaPhiIn) < cutdphi[cat + bin * 9]) and (fabs(deltaEtaIn) < cutdeta[cat + bin * 9]) and
0289           (eSeedOverPin > cuteopin[cat + bin * 9]) and (scEt > cutet[cat + bin * 9]))
0290         result += 1.;
0291     }
0292 
0293     std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip_gsf");
0294     if (ip < cutip[cat + bin * 9])
0295       result += 8;
0296 
0297     std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutfmishits");
0298     std::vector<double> cutdcotdist = cuts_.getParameter<std::vector<double> >("cutdcotdist");
0299 
0300     float dist = (electron->convDist() == -9999. ? 9999 : electron->convDist());
0301     float dcot = (electron->convDcot() == -9999. ? 9999 : electron->convDcot());
0302 
0303     float dcotdistcomb =
0304         ((0.04 - std::max(fabs(dist), fabs(dcot))) > 0 ? (0.04 - std::max(fabs(dist), fabs(dcot))) : 0);
0305 
0306     if ((mishits < cutmishits[cat + bin * 9]) and (dcotdistcomb < cutdcotdist[cat + bin * 9]))
0307       result += 4;
0308 
0309     return result;
0310   }
0311 
0312   if (type_ == "classbased" && (version_ == "V06" || version_.empty())) {
0313     std::vector<double> cutIsoSum = cuts_.getParameter<std::vector<double> >("cutiso_sum");
0314     std::vector<double> cutIsoSumCorr = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
0315     std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
0316     std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
0317     std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
0318     std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
0319     std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
0320     std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutfmishits");
0321     std::vector<double> cutdcotdist = cuts_.getParameter<std::vector<double> >("cutdcotdist");
0322     std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip_gsf");
0323     std::vector<double> cutIsoSumCorrl = cuts_.getParameter<std::vector<double> >("cutiso_sumoetl");
0324     std::vector<double> cuthoel = cuts_.getParameter<std::vector<double> >("cuthoel");
0325     std::vector<double> cutseel = cuts_.getParameter<std::vector<double> >("cutseel");
0326     std::vector<double> cutdphil = cuts_.getParameter<std::vector<double> >("cutdphiinl");
0327     std::vector<double> cutdetal = cuts_.getParameter<std::vector<double> >("cutdetainl");
0328     std::vector<double> cutipl = cuts_.getParameter<std::vector<double> >("cutip_gsfl");
0329 
0330     int result = 0;
0331 
0332     const int ncuts = 10;
0333     std::vector<bool> cut_results(ncuts, false);
0334 
0335     float iso_sum = tkIso + ecalIso + hcalIso;
0336     float scEta = electron->superCluster()->eta();
0337     if (fabs(scEta) > 1.5)
0338       iso_sum += (fabs(scEta) - 1.5) * 1.09;
0339 
0340     float iso_sumoet = iso_sum * (40. / scEt);
0341 
0342     float eseedopincor = eSeedOverPin + fBrem;
0343     if (fBrem < 0)
0344       eseedopincor = eSeedOverPin;
0345 
0346     float dist = (electron->convDist() == -9999. ? 9999 : electron->convDist());
0347     float dcot = (electron->convDcot() == -9999. ? 9999 : electron->convDcot());
0348 
0349     float dcotdistcomb =
0350         ((0.04 - std::max(fabs(dist), fabs(dcot))) > 0 ? (0.04 - std::max(fabs(dist), fabs(dcot))) : 0);
0351 
0352     for (int cut = 0; cut < ncuts; cut++) {
0353       switch (cut) {
0354         case 0:
0355           cut_results[cut] = compute_cut(fabs(deltaEtaIn), scEt, cutdetal[cat], cutdeta[cat]);
0356           break;
0357         case 1:
0358           cut_results[cut] = compute_cut(fabs(deltaPhiIn), scEt, cutdphil[cat], cutdphi[cat]);
0359           break;
0360         case 2:
0361           cut_results[cut] = (eseedopincor > cuteopin[cat]);
0362           break;
0363         case 3:
0364           cut_results[cut] = compute_cut(hOverE, scEt, cuthoel[cat], cuthoe[cat]);
0365           break;
0366         case 4:
0367           cut_results[cut] = compute_cut(sigmaee, scEt, cutseel[cat], cutsee[cat]);
0368           break;
0369         case 5:
0370           cut_results[cut] = compute_cut(iso_sumoet, scEt, cutIsoSumCorrl[cat], cutIsoSumCorr[cat]);
0371           break;
0372         case 6:
0373           cut_results[cut] = (iso_sum < cutIsoSum[cat]);
0374           break;
0375         case 7:
0376           cut_results[cut] = compute_cut(fabs(ip), scEt, cutipl[cat], cutip[cat]);
0377           break;
0378         case 8:
0379           cut_results[cut] = (mishits < cutmishits[cat]);
0380           break;
0381         case 9:
0382           cut_results[cut] = (dcotdistcomb < cutdcotdist[cat]);
0383           break;
0384       }
0385     }
0386 
0387     // ID part
0388     if (cut_results[0] && cut_results[1] && cut_results[2] && cut_results[3] && cut_results[4])
0389       result = result + 1;
0390 
0391     // ISO part
0392     if (cut_results[5] && cut_results[6])
0393       result = result + 2;
0394 
0395     // IP part
0396     if (cut_results[7])
0397       result = result + 8;
0398 
0399     // Conversion part
0400     if (cut_results[8] && cut_results[9])
0401       result = result + 4;
0402 
0403     return result;
0404   }
0405 
0406   return -1.;
0407 }
0408 
0409 bool CutBasedElectronID::compute_cut(double x, double et, double cut_min, double cut_max, bool gtn) {
0410   float et_min = 10;
0411   float et_max = 40;
0412 
0413   bool accept = false;
0414   float cut = cut_max;  //  the cut at et=40 GeV
0415 
0416   if (et < et_max) {
0417     cut = cut_min + (1 / et_min - 1 / et) * (cut_max - cut_min) / (1 / et_min - 1 / et_max);
0418   }
0419 
0420   if (et < et_min) {
0421     cut = cut_min;
0422   }
0423 
0424   if (gtn) {  // useful for e/p cut which is gt
0425     accept = (x >= cut);
0426   } else {
0427     accept = (x <= cut);
0428   }
0429 
0430   //std::cout << x << " " << cut_min << " " << cut << " " << cut_max << " " << et << " " << accept << std::endl;
0431   return accept;
0432 }
0433 
0434 double CutBasedElectronID::robustSelection(const reco::GsfElectron* electron,
0435                                            const edm::Event& e,
0436                                            const edm::EventSetup& es) {
0437   double scTheta = (2 * atan(exp(-electron->superCluster()->eta())));
0438   double scEt = electron->superCluster()->energy() * sin(scTheta);
0439   double eta = electron->p4().Eta();
0440   double eOverP = electron->eSuperClusterOverP();
0441   double hOverE = electron->hadronicOverEm();
0442   double sigmaee = electron->sigmaIetaIeta();
0443   double e25Max = electron->e2x5Max();
0444   double e15 = electron->e1x5();
0445   double e55 = electron->e5x5();
0446   double e25Maxoe55 = e25Max / e55;
0447   double e15oe55 = e15 / e55;
0448   double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
0449   double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
0450 
0451   double ip = 0;
0452   int mishits = electron->gsfTrack()->missingInnerHits();
0453   double tkIso = electron->dr03TkSumPt();
0454   double ecalIso = electron->dr04EcalRecHitSumEt();
0455   double ecalIsoPed = (electron->isEB()) ? std::max(0., ecalIso - 1.) : ecalIso;
0456   double hcalIso = electron->dr04HcalTowerSumEt();
0457   double hcalIso1 = electron->dr04HcalTowerSumEt(1);
0458   double hcalIso2 = electron->dr04HcalTowerSumEt(2);
0459 
0460   if (version_ == "V00") {
0461     sigmaee = electron->sigmaEtaEta();
0462     if (electron->isEE())
0463       sigmaee = sigmaee - 0.02 * (fabs(eta) - 2.3);  //correct sigmaetaeta dependence on eta in endcap
0464   }
0465 
0466   if (version_ == "V03" or version_ == "V04") {
0467     edm::Handle<reco::BeamSpot> pBeamSpot;
0468     // uses the same name for the vertex collection to avoid adding more new names
0469     e.getByToken(verticesCollection_, pBeamSpot);
0470     if (pBeamSpot.isValid()) {
0471       const reco::BeamSpot* bspot = pBeamSpot.product();
0472       const math::XYZPoint& bspotPosition = bspot->position();
0473       ip = fabs(electron->gsfTrack()->dxy(bspotPosition));
0474     } else
0475       ip = fabs(electron->gsfTrack()->dxy());
0476   }
0477 
0478   if (version_ == "V04" or version_ == "V05") {
0479     ecalIso = electron->dr03EcalRecHitSumEt();
0480     ecalIsoPed = (electron->isEB()) ? std::max(0., ecalIso - 1.) : ecalIso;
0481     hcalIso = electron->dr03HcalTowerSumEt();
0482     hcalIso1 = electron->dr03HcalTowerSumEt(1);
0483     hcalIso2 = electron->dr03HcalTowerSumEt(2);
0484   }
0485 
0486   if (version_ == "V05") {
0487     edm::Handle<reco::VertexCollection> vtxH;
0488     e.getByToken(verticesCollection_, vtxH);
0489     if (!vtxH->empty()) {
0490       reco::VertexRef vtx(vtxH, 0);
0491       ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(), vtx->y(), vtx->z())));
0492     } else
0493       ip = fabs(electron->gsfTrack()->dxy());
0494   }
0495 
0496   // .....................................................................................
0497   std::vector<double> cut;
0498   // ROBUST Selection
0499   if (type_ == "robust") {
0500     double result = 0;
0501 
0502     // hoe, sigmaEtaEta, dPhiIn, dEtaIn
0503     if (electron->isEB())
0504       cut = cuts_.getParameter<std::vector<double> >("barrel");
0505     else
0506       cut = cuts_.getParameter<std::vector<double> >("endcap");
0507     // check isolations: if only isolation passes result = 2
0508     if (quality_ == "highenergy") {
0509       if ((tkIso > cut[6] || hcalIso2 > cut[11]) ||
0510           (electron->isEB() && ((ecalIso + hcalIso1) > cut[7] + cut[8] * scEt)) ||
0511           (electron->isEE() && (scEt >= 50.) && ((ecalIso + hcalIso1) > cut[7] + cut[8] * (scEt - 50))) ||
0512           (electron->isEE() && (scEt < 50.) && ((ecalIso + hcalIso1) > cut[9] + cut[10] * (scEt - 50))))
0513         result = 0;
0514       else
0515         result = 2;
0516     } else {
0517       if ((tkIso > cut[6]) || (ecalIso > cut[7]) || (hcalIso > cut[8]) || (hcalIso1 > cut[9]) || (hcalIso2 > cut[10]) ||
0518           (tkIso / electron->p4().Pt() > cut[11]) || (ecalIso / electron->p4().Pt() > cut[12]) ||
0519           (hcalIso / electron->p4().Pt() > cut[13]) || ((tkIso + ecalIso + hcalIso) > cut[14]) ||
0520           (((tkIso + ecalIso + hcalIso) / electron->p4().Pt()) > cut[15]) ||
0521           ((tkIso + ecalIsoPed + hcalIso) > cut[16]) ||
0522           (((tkIso + ecalIsoPed + hcalIso) / electron->p4().Pt()) > cut[17]))
0523         result = 0.;
0524       else
0525         result = 2.;
0526     }
0527 
0528     if ((hOverE < cut[0]) && (sigmaee < cut[1]) && (fabs(deltaPhiIn) < cut[2]) && (fabs(deltaEtaIn) < cut[3]) &&
0529         (e25Maxoe55 > cut[4] && e15oe55 > cut[5]) && (sigmaee >= cut[18]) && (eOverP > cut[19] && eOverP < cut[20])) {
0530       result = result + 1;
0531     }
0532 
0533     if (ip > cut[21])
0534       return result;
0535     if (mishits > cut[22])  // expected missing hits
0536       return result;
0537     // positive cut[23] means to demand a valid hit in 1st layer PXB
0538     if (cut[23] > 0 &&
0539         !electron->gsfTrack()->hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1))
0540       return result;
0541 
0542     // cut[24]: Dist cut[25]: dcot
0543     float dist = fabs(electron->convDist());
0544     float dcot = fabs(electron->convDcot());
0545     bool isConversion = (cut[24] > 99. || cut[25] > 99.) ? false : (dist < cut[24] && dcot < cut[25]);
0546     if (isConversion)
0547       return result;
0548 
0549     result += 4;
0550 
0551     return result;
0552   }
0553 
0554   return -1.;
0555 }