CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoEgamma/ElectronIdentification/src/CutBasedElectronID.cc

Go to the documentation of this file.
00001 #include "RecoEgamma/ElectronIdentification/interface/CutBasedElectronID.h"
00002 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00003 #include "DataFormats/VertexReco/interface/Vertex.h"
00004 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00005 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00006 //#include "RecoEgamma/EgammaTools/interface/ConversionFinder.h"
00007 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00008 
00009 #include <algorithm>
00010 
00011 void CutBasedElectronID::setup(const edm::ParameterSet& conf) {
00012 
00013   // Get all the parameters
00014   baseSetup(conf);
00015 
00016   type_ = conf.getParameter<std::string>("electronIDType");
00017   quality_ = conf.getParameter<std::string>("electronQuality");
00018   version_ = conf.getParameter<std::string>("electronVersion");
00019   verticesCollection_ = conf.getParameter<edm::InputTag>("verticesCollection");
00020 
00021   if (type_ == "classbased" and (version_ == "V03" or version_ == "V04" or version_ == "V05" or version_ == "")) {
00022     wantBinning_ = conf.getParameter<bool>("etBinning");
00023     newCategories_ = conf.getParameter<bool>("additionalCategories");
00024   }
00025 
00026   if (type_ == "robust" || type_ == "classbased") {
00027     std::string stringCut = type_+quality_+"EleIDCuts"+version_;
00028     cuts_ = conf.getParameter<edm::ParameterSet>(stringCut);
00029   }
00030   else {
00031     throw cms::Exception("Configuration")
00032       << "Invalid electronType parameter in CutBasedElectronID: must be robust or classbased\n";
00033   }
00034 }
00035 
00036 double CutBasedElectronID::result(const reco::GsfElectron* electron ,
00037                                   const edm::Event& e ,
00038                                   const edm::EventSetup& es) {
00039 
00040   if (type_ == "classbased")
00041     return cicSelection(electron, e, es);
00042   else if (type_ == "robust")
00043     return robustSelection(electron, e, es);
00044 
00045   return 0;
00046 
00047 }
00048 
00049 int CutBasedElectronID::classify(const reco::GsfElectron* electron) {
00050 
00051   double eta = fabs(electron->superCluster()->eta());
00052   double eOverP = electron->eSuperClusterOverP();
00053   double fBrem = electron->fbrem();
00054 
00055   int cat = -1;
00056   if (version_ == "V00" || version_ == "V01") {
00057     if((electron->isEB() && fBrem<0.06) || (electron->isEE() && fBrem<0.1))
00058       cat=1;
00059     else if (eOverP < 1.2 && eOverP > 0.8)
00060       cat=0;
00061     else
00062       cat=2;
00063 
00064     return cat;
00065 
00066   } else if (version_ == "V02") {
00067     if (electron->isEB()) {       // BARREL
00068       if(fBrem < 0.12)
00069         cat=1;
00070       else if (eOverP < 1.2 && eOverP > 0.9)
00071         cat=0;
00072       else
00073         cat=2;
00074     } else {                     // ENDCAP
00075       if(fBrem < 0.2)
00076         cat=1;
00077       else if (eOverP < 1.22 && eOverP > 0.82)
00078         cat=0;
00079       else
00080         cat=2;
00081     }
00082 
00083     return cat;
00084 
00085   } else {
00086     if (electron->isEB()) {
00087       if ((fBrem >= 0.12) and (eOverP > 0.9) and (eOverP < 1.2))
00088         cat = 0;
00089       else if (((eta >  .445   and eta <  .45  ) or
00090                 (eta >  .79    and eta <  .81  ) or
00091                 (eta > 1.137   and eta < 1.157 ) or
00092                 (eta > 1.47285 and eta < 1.4744)) and newCategories_)
00093         cat = 6;
00094       else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
00095         cat = 8;
00096       else if (fBrem < 0.12)
00097         cat = 1;
00098       else
00099         cat = 2;
00100     } else {
00101       if ((fBrem >= 0.2) and (eOverP > 0.82) and (eOverP < 1.22))
00102         cat = 3;
00103       else if (eta > 1.5 and eta <  1.58 and newCategories_)
00104         cat = 7;
00105       else if (electron->trackerDrivenSeed() and !electron->ecalDrivenSeed() and newCategories_)
00106         cat = 8;
00107       else if (fBrem < 0.2)
00108         cat = 4;
00109       else
00110         cat = 5;
00111     }
00112 
00113     return cat;
00114   }
00115 
00116   return -1;
00117 }
00118 
00119 double CutBasedElectronID::cicSelection(const reco::GsfElectron* electron,
00120                                         const edm::Event& e,
00121                                         const edm::EventSetup& es) {
00122 
00123   double scTheta = (2*atan(exp(-electron->superCluster()->eta())));
00124   double scEt = electron->superCluster()->energy()*sin(scTheta);
00125 
00126   double eta = electron->p4().Eta();
00127   double eOverP = electron->eSuperClusterOverP();
00128   double eSeedOverPin = electron->eSeedClusterOverP();
00129   double fBrem = electron->fbrem();
00130   double hOverE = electron->hadronicOverEm();
00131   double sigmaee = electron->sigmaIetaIeta(); //sqrt(vLocCov[0]);
00132   double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
00133   double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
00134 
00135   double ip = 0;
00136   int mishits = electron->gsfTrack()->trackerExpectedHitsInner().numberOfHits();
00137   double tkIso = electron->dr03TkSumPt();
00138   double ecalIso = electron->dr04EcalRecHitSumEt();
00139   double hcalIso = electron->dr04HcalTowerSumEt();
00140 
00141   if (version_ == "V00") {
00142     sigmaee = electron->sigmaEtaEta();//sqrt(vCov[0]);
00143     if (electron->isEE())
00144       sigmaee = sigmaee - 0.02*(fabs(eta) - 2.3);   //correct sigmaetaeta dependence on eta in endcap
00145   }
00146 
00147   if (version_ != "V01" or version_ != "V00") {
00148     edm::Handle<reco::VertexCollection> vtxH;
00149     e.getByLabel(verticesCollection_, vtxH);
00150     if (vtxH->size() != 0) {
00151       reco::VertexRef vtx(vtxH, 0);
00152       ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(),vtx->y(),vtx->z())));
00153     } else
00154       ip = fabs(electron->gsfTrack()->dxy());
00155 
00156     if (electron->isEB()) {
00157       sigmaee = electron->sigmaIetaIeta(); //sqrt(vCov[0]);
00158     }
00159   }
00160 
00161   std::vector<double> cut;
00162 
00163   int cat = classify(electron);
00164   int eb;
00165 
00166   if (electron->isEB())
00167     eb = 0;
00168   else
00169     eb = 1;
00170 
00171   // LOOSE and TIGHT Selections
00172   if (type_ == "classbased" && (version_ == "V01" || version_ == "V00")) {
00173 
00174     if ((eOverP < 0.8) && (fBrem < 0.2))
00175       return 0.;
00176 
00177     cut = cuts_.getParameter<std::vector<double> >("hOverE");
00178     if (hOverE > cut[cat+4*eb])
00179       return 0.;
00180 
00181     cut = cuts_.getParameter<std::vector<double> >("sigmaEtaEta");
00182     if (sigmaee > cut[cat+4*eb])
00183       return 0.;
00184 
00185     cut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
00186     if (eOverP < 1.5) {
00187       if (fabs(deltaPhiIn) > cut[cat+4*eb])
00188         return 0.;
00189     } else {
00190       if (fabs(deltaPhiIn) > cut[3+4*eb])
00191         return 0.;
00192     }
00193 
00194     cut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
00195     if (fabs(deltaEtaIn) > cut[cat+4*eb])
00196       return 0.;
00197 
00198     cut = cuts_.getParameter<std::vector<double> >("eSeedOverPin");
00199     if (eSeedOverPin < cut[cat+4*eb])
00200       return 0.;
00201 
00202     if (quality_ == "tight")
00203       if (eOverP < 0.9*(1-fBrem))
00204         return 0.;
00205 
00206     return 1.;
00207   }
00208 
00209   if (type_ == "classbased" and version_ == "V02") {
00210     double result = 0.;
00211 
00212     int bin = 0;
00213 
00214     if (scEt < 20.)
00215       bin = 2;
00216     else if (scEt > 30.)
00217       bin = 0;
00218     else
00219       bin = 1;
00220 
00221     if (fBrem > 0)
00222       eSeedOverPin = eSeedOverPin + fBrem;
00223 
00224     if (bin != 2) {
00225       tkIso = tkIso*pow(40./scEt, 2);
00226       ecalIso = ecalIso*pow(40./scEt, 2);
00227       hcalIso = hcalIso*pow(40./scEt, 2);
00228     }
00229 
00230     std::vector<double> cutTk = cuts_.getParameter<std::vector<double> >("cutisotk");
00231     std::vector<double> cutEcal = cuts_.getParameter<std::vector<double> >("cutisoecal");
00232     std::vector<double> cutHcal = cuts_.getParameter<std::vector<double> >("cutisohcal");
00233     if ((tkIso > cutTk[cat+3*eb+bin*6]) ||
00234         (ecalIso > cutEcal[cat+3*eb+bin*6]) ||
00235         (hcalIso > cutHcal[cat+3*eb+bin*6]))
00236       result = 0.;
00237     else
00238       result = 2.;
00239 
00240     if (fBrem > -2) {
00241       std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
00242       std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
00243       std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
00244       std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
00245       std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
00246       std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
00247       std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip");
00248       std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutmishits");
00249       if ((hOverE < cuthoe[cat+3*eb+bin*6]) and
00250           (sigmaee < cutsee[cat+3*eb+bin*6]) and
00251           (fabs(deltaPhiIn) < cutdphi[cat+3*eb+bin*6]) and
00252           (fabs(deltaEtaIn) < cutdeta[cat+3*eb+bin*6]) and
00253           (eSeedOverPin > cuteopin[cat+3*eb+bin*6]) and
00254           (ip < cutip[cat+3*eb+bin*6]) and
00255           (mishits < cutmishits[cat+3*eb+bin*6]))
00256         result = result + 1.;
00257     }
00258     return result;
00259   }
00260 
00261   if (version_ == "V03" or version_ == "V04" or version_ == "V05") {
00262     double result = 0.;
00263 
00264     int bin = 0;
00265 
00266     if (wantBinning_) {
00267       if (scEt < 20.)
00268         bin = 2;
00269       else if (scEt > 30.)
00270         bin = 0;
00271       else
00272         bin = 1;
00273     }
00274 
00275     if (fBrem > 0)
00276       eSeedOverPin = eSeedOverPin + fBrem;
00277 
00278     float iso_sum = tkIso + ecalIso + hcalIso;
00279     float iso_sum_corrected = iso_sum*pow(40./scEt, 2);
00280 
00281     std::vector<double> cutIsoSum = cuts_.getParameter<std::vector<double> >("cutiso_sum");
00282     std::vector<double> cutIsoSumCorr = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
00283     if ((iso_sum < cutIsoSum[cat+bin*9]) and
00284         (iso_sum_corrected < cutIsoSumCorr[cat+bin*9]))
00285       result += 2.;
00286 
00287     if (fBrem > -2) {
00288       std::vector<double> cuthoe = cuts_.getParameter<std::vector<double> >("cuthoe");
00289       std::vector<double> cutsee = cuts_.getParameter<std::vector<double> >("cutsee");
00290       std::vector<double> cutdphi = cuts_.getParameter<std::vector<double> >("cutdphiin");
00291       std::vector<double> cutdeta = cuts_.getParameter<std::vector<double> >("cutdetain");
00292       std::vector<double> cuteopin = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
00293       std::vector<double> cutet = cuts_.getParameter<std::vector<double> >("cutet");
00294 
00295       if ((hOverE < cuthoe[cat+bin*9]) and
00296           (sigmaee < cutsee[cat+bin*9]) and
00297           (fabs(deltaPhiIn) < cutdphi[cat+bin*9]) and
00298           (fabs(deltaEtaIn) < cutdeta[cat+bin*9]) and
00299           (eSeedOverPin > cuteopin[cat+bin*9]) and
00300           (scEt > cutet[cat+bin*9]))
00301         result += 1.;
00302     }
00303 
00304     std::vector<double> cutip = cuts_.getParameter<std::vector<double> >("cutip_gsf");
00305     if (ip < cutip[cat+bin*9])
00306       result += 8;
00307     
00308     std::vector<double> cutmishits = cuts_.getParameter<std::vector<double> >("cutfmishits");
00309     std::vector<double> cutdcotdist = cuts_.getParameter<std::vector<double> >("cutdcotdist");
00310     
00311     float dist = (electron->convDist() == -9999.? 9999:electron->convDist());
00312     float dcot = (electron->convDcot() == -9999.? 9999:electron->convDcot());
00313 
00314     float dcotdistcomb = ((0.04 - std::max(dist, dcot)) > 0?(0.04 - std::max(dist, dcot)):0);
00315 
00316     if ((mishits < cutmishits[cat+bin*9]) and
00317         (dcotdistcomb < cutdcotdist[cat+bin*9]))
00318       result += 4;
00319 
00320     return result;
00321   }
00322 
00323   if (type_ == "classbased" && (version_ == "V06" || version_ == "")) { 
00324     std::vector<double> cutIsoSum      = cuts_.getParameter<std::vector<double> >("cutiso_sum");
00325     std::vector<double> cutIsoSumCorr  = cuts_.getParameter<std::vector<double> >("cutiso_sumoet");
00326     std::vector<double> cuthoe         = cuts_.getParameter<std::vector<double> >("cuthoe");
00327     std::vector<double> cutsee         = cuts_.getParameter<std::vector<double> >("cutsee");
00328     std::vector<double> cutdphi        = cuts_.getParameter<std::vector<double> >("cutdphiin");
00329     std::vector<double> cutdeta        = cuts_.getParameter<std::vector<double> >("cutdetain");
00330     std::vector<double> cuteopin       = cuts_.getParameter<std::vector<double> >("cuteseedopcor");
00331     std::vector<double> cutmishits     = cuts_.getParameter<std::vector<double> >("cutfmishits");
00332     std::vector<double> cutdcotdist    = cuts_.getParameter<std::vector<double> >("cutdcotdist");
00333     std::vector<double> cutip          = cuts_.getParameter<std::vector<double> >("cutip_gsf");
00334     std::vector<double> cutIsoSumCorrl = cuts_.getParameter<std::vector<double> >("cutiso_sumoetl");
00335     std::vector<double> cuthoel        = cuts_.getParameter<std::vector<double> >("cuthoel");
00336     std::vector<double> cutseel        = cuts_.getParameter<std::vector<double> >("cutseel");
00337     std::vector<double> cutdphil       = cuts_.getParameter<std::vector<double> >("cutdphiinl");
00338     std::vector<double> cutdetal       = cuts_.getParameter<std::vector<double> >("cutdetainl");
00339     std::vector<double> cutipl         = cuts_.getParameter<std::vector<double> >("cutip_gsfl");
00340     
00341     int result = 0;
00342     
00343     const int ncuts = 10;
00344     std::vector<bool> cut_results(ncuts, false);
00345     
00346     float iso_sum = tkIso + ecalIso + hcalIso;
00347     float scEta = electron->superCluster()->eta();
00348     if(fabs(scEta)>1.5) 
00349       iso_sum += (fabs(scEta)-1.5)*1.09;
00350     
00351     float iso_sumoet = iso_sum*(40./scEt);
00352     
00353     float eseedopincor = eSeedOverPin + fBrem;
00354     if(fBrem < 0)
00355       eseedopincor = eSeedOverPin;
00356 
00357     float dist = (electron->convDist() == -9999.? 9999:electron->convDist());
00358     float dcot = (electron->convDcot() == -9999.? 9999:electron->convDcot());
00359 
00360     float dcotdistcomb = ((0.04 - std::max(dist, dcot)) > 0?(0.04 - std::max(dist, dcot)):0);
00361 
00362     for (int cut=0; cut<ncuts; cut++) {
00363       switch (cut) {
00364       case 0:
00365         cut_results[cut] = compute_cut(fabs(deltaEtaIn), scEt, cutdetal[cat], cutdeta[cat]);
00366         break;
00367       case 1:
00368         cut_results[cut] = compute_cut(fabs(deltaPhiIn), scEt, cutdphil[cat], cutdphi[cat]);
00369         break;
00370       case 2:
00371         cut_results[cut] = (eseedopincor > cuteopin[cat]);
00372         break;
00373       case 3:
00374         cut_results[cut] = compute_cut(hOverE, scEt, cuthoel[cat], cuthoe[cat]);
00375         break;
00376       case 4:
00377         cut_results[cut] = compute_cut(sigmaee, scEt, cutseel[cat], cutsee[cat]);
00378         break;
00379       case 5:
00380         cut_results[cut] = compute_cut(iso_sumoet, scEt, cutIsoSumCorrl[cat], cutIsoSumCorr[cat]);
00381         break;
00382       case 6:
00383         cut_results[cut] = (iso_sum < cutIsoSum[cat]);
00384         break;
00385       case 7:
00386         cut_results[cut] = compute_cut(fabs(ip), scEt, cutipl[cat], cutip[cat]);
00387         break;
00388       case 8:
00389         cut_results[cut] = (mishits < cutmishits[cat]);
00390         break;
00391       case 9:
00392         cut_results[cut] = (dcotdistcomb < cutdcotdist[cat]);
00393         break;
00394       }
00395     }
00396     
00397     // ID part
00398     if (cut_results[0] & cut_results[1] & cut_results[2] & cut_results[3] & cut_results[4])
00399       result = result + 1;
00400     
00401     // ISO part
00402     if (cut_results[5] & cut_results[6])
00403       result = result + 2;
00404     
00405     // IP part
00406     if (cut_results[7])
00407       result = result + 8;
00408     
00409     // Conversion part
00410     if (cut_results[8] and cut_results[9])
00411       result = result + 4;
00412     
00413     return result;
00414   }
00415 
00416   return -1.;
00417 }
00418 
00419 
00420 bool CutBasedElectronID::compute_cut(double x, double et, double cut_min, double cut_max, bool gtn) {
00421 
00422   float et_min = 10;
00423   float et_max = 40;
00424 
00425   bool accept = false;
00426   float cut = cut_max; //  the cut at et=40 GeV
00427 
00428   if(et < et_max) {
00429     cut = cut_min + (1/et_min - 1/et)*(cut_max - cut_min)/(1/et_min - 1/et_max);
00430   } 
00431   
00432   if(et < et_min) {
00433     cut = cut_min;
00434   } 
00435 
00436   if(gtn) {   // useful for e/p cut which is gt
00437     accept = (x >= cut);
00438   } 
00439   else {
00440     accept = (x <= cut);
00441   }
00442 
00443   //std::cout << x << " " << cut_min << " " << cut << " " << cut_max << " " << et << " " << accept << std::endl;
00444   return accept;
00445 }
00446 
00447 double CutBasedElectronID::robustSelection(const reco::GsfElectron* electron ,
00448                                            const edm::Event& e ,
00449                                            const edm::EventSetup& es) {
00450 
00451   double scTheta = (2*atan(exp(-electron->superCluster()->eta())));
00452   double scEt = electron->superCluster()->energy()*sin(scTheta);
00453   double eta = electron->p4().Eta();
00454   double eOverP = electron->eSuperClusterOverP();
00455   double hOverE = electron->hadronicOverEm();
00456   double sigmaee = electron->sigmaIetaIeta();
00457   double e25Max = electron->e2x5Max();
00458   double e15 = electron->e1x5();
00459   double e55 = electron->e5x5();
00460   double e25Maxoe55 = e25Max/e55;
00461   double e15oe55 = e15/e55 ;
00462   double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
00463   double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
00464 
00465   double ip = 0;
00466   int mishits = electron->gsfTrack()->trackerExpectedHitsInner().numberOfHits();
00467   double tkIso = electron->dr03TkSumPt();
00468   double ecalIso = electron->dr04EcalRecHitSumEt();
00469   double ecalIsoPed = (electron->isEB())?std::max(0.,ecalIso-1.):ecalIso;
00470   double hcalIso = electron->dr04HcalTowerSumEt();
00471   double hcalIso1 = electron->dr04HcalDepth1TowerSumEt();
00472   double hcalIso2 = electron->dr04HcalDepth2TowerSumEt();
00473 
00474   if (version_ == "V00") {
00475     sigmaee = electron->sigmaEtaEta();
00476      if (electron->isEE())
00477        sigmaee = sigmaee - 0.02*(fabs(eta) - 2.3);   //correct sigmaetaeta dependence on eta in endcap
00478   }
00479 
00480   if (version_ == "V03" or version_ == "V04") {
00481     edm::Handle<reco::BeamSpot> pBeamSpot;
00482     // uses the same name for the vertex collection to avoid adding more new names
00483     e.getByLabel(verticesCollection_, pBeamSpot);
00484     if (pBeamSpot.isValid()) {
00485       const reco::BeamSpot *bspot = pBeamSpot.product();
00486       const math::XYZPoint bspotPosition = bspot->position();
00487       ip = fabs(electron->gsfTrack()->dxy(bspotPosition));
00488     } else
00489       ip = fabs(electron->gsfTrack()->dxy());
00490   }
00491 
00492   if (version_ == "V04" or version_ == "V05") {
00493     ecalIso = electron->dr03EcalRecHitSumEt();
00494     ecalIsoPed = (electron->isEB())?std::max(0.,ecalIso-1.):ecalIso;
00495     hcalIso = electron->dr03HcalTowerSumEt();
00496     hcalIso1 = electron->dr03HcalDepth1TowerSumEt();
00497     hcalIso2 = electron->dr03HcalDepth2TowerSumEt();
00498   }
00499 
00500   if (version_ == "V05") {
00501     edm::Handle<reco::VertexCollection> vtxH;
00502     e.getByLabel(verticesCollection_, vtxH);
00503     if (vtxH->size() != 0) {
00504       reco::VertexRef vtx(vtxH, 0);
00505       ip = fabs(electron->gsfTrack()->dxy(math::XYZPoint(vtx->x(),vtx->y(),vtx->z())));
00506     } else
00507       ip = fabs(electron->gsfTrack()->dxy());
00508   }
00509 
00510   // .....................................................................................
00511   std::vector<double> cut;
00512   // ROBUST Selection
00513   if (type_ == "robust") {
00514 
00515     double result = 0;
00516 
00517     // hoe, sigmaEtaEta, dPhiIn, dEtaIn
00518     if (electron->isEB())
00519       cut = cuts_.getParameter<std::vector<double> >("barrel");
00520     else
00521       cut = cuts_.getParameter<std::vector<double> >("endcap");
00522     // check isolations: if only isolation passes result = 2
00523     if (quality_ == "highenergy") {
00524       if ((tkIso > cut[6] || hcalIso2 > cut[11]) ||
00525           (electron->isEB() && ((ecalIso + hcalIso1) > cut[7]+cut[8]*scEt)) ||
00526           (electron->isEE() && (scEt >= 50.) && ((ecalIso + hcalIso1) > cut[7]+cut[8]*(scEt-50))) ||
00527           (electron->isEE() && (scEt < 50.) && ((ecalIso + hcalIso1) > cut[9]+cut[10]*(scEt-50))))
00528         result = 0;
00529       else
00530         result = 2;
00531     } else {
00532       if ((tkIso > cut[6]) || (ecalIso > cut[7]) || (hcalIso > cut[8]) || (hcalIso1 > cut[9]) || (hcalIso2 > cut[10]) ||
00533           (tkIso/electron->p4().Pt() > cut[11]) || (ecalIso/electron->p4().Pt() > cut[12]) || (hcalIso/electron->p4().Pt() > cut[13]) ||
00534           ((tkIso+ecalIso+hcalIso)>cut[14]) || (((tkIso+ecalIso+hcalIso)/ electron->p4().Pt()) > cut[15]) ||
00535           ((tkIso+ecalIsoPed+hcalIso)>cut[16]) || (((tkIso+ecalIsoPed+hcalIso)/ electron->p4().Pt()) > cut[17])  )
00536         result = 0.;
00537       else
00538         result = 2.;
00539     }
00540 
00541     if ((hOverE < cut[0]) && (sigmaee < cut[1]) && (fabs(deltaPhiIn) < cut[2]) &&
00542         (fabs(deltaEtaIn) < cut[3]) && (e25Maxoe55 > cut[4] && e15oe55 > cut[5]) &&
00543         (sigmaee >= cut[18]) && (eOverP > cut[19] &&  eOverP < cut[20]) )
00544      { result = result + 1 ; }
00545 
00546     if (ip > cut[21])
00547       return result;
00548     if (mishits > cut[22]) // expected missing hits
00549       return result;
00550     // positive cut[23] means to demand a valid hit in 1st layer PXB
00551     if (cut[23] >0 && not (electron->gsfTrack()->hitPattern().hasValidHitInFirstPixelBarrel()))
00552       return result;
00553 
00554     // cut[24]: Dist cut[25]: dcot
00555     float dist = fabs(electron->convDist());
00556     float dcot = fabs(electron->convDcot());
00557     bool isConversion = (cut[24]>99. || cut[25]>99.)?false:(dist < cut[24] && dcot < cut[25]);
00558     if (isConversion)
00559       return result ;
00560     
00561     result += 4 ;
00562 
00563     return result ;
00564    }
00565 
00566   return -1. ;
00567  }