CMS 3D CMS Logo

CutBasedElectronID Class Reference

#include <RecoEgamma/ElectronIdentification/interface/CutBasedElectronID.h>

Inheritance diagram for CutBasedElectronID:

ElectronIDAlgo

List of all members.

Public Member Functions

int classify (const reco::GsfElectron *)
 CutBasedElectronID ()
double result (const reco::GsfElectron *, const edm::Event &, const edm::EventSetup &)
void setup (const edm::ParameterSet &conf)
virtual ~CutBasedElectronID ()

Private Attributes

edm::ParameterSet cuts_
std::string quality_


Detailed Description

Definition at line 6 of file CutBasedElectronID.h.


Constructor & Destructor Documentation

CutBasedElectronID::CutBasedElectronID (  )  [inline]

Definition at line 10 of file CutBasedElectronID.h.

00010 {};

virtual CutBasedElectronID::~CutBasedElectronID (  )  [inline, virtual]

Definition at line 12 of file CutBasedElectronID.h.

00012 {};


Member Function Documentation

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

Definition at line 23 of file CutBasedElectronID.cc.

References reco::GsfElectron::eSuperClusterOverP(), eta, reco::Particle::p4(), reco::GsfElectron::trackMomentumAtVtx(), and reco::GsfElectron::trackMomentumOut().

Referenced by result().

00023                                                                 {
00024   
00025   double eta = electron->p4().Eta();
00026   double eOverP = electron->eSuperClusterOverP();
00027   double pin  = electron->trackMomentumAtVtx().R(); 
00028   double pout = electron->trackMomentumOut().R(); 
00029   double fBrem = (pin-pout)/pin;
00030   
00031   int cat;
00032   
00033   if((fabs(eta)<1.479 && fBrem<0.06) || (fabs(eta)>1.479 && fBrem<0.1)) 
00034     cat=1;
00035   else if (eOverP < 1.2 && eOverP > 0.8) 
00036     cat=0;
00037   else 
00038     cat=2;
00039   
00040   return cat;
00041 }

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

Reimplemented from ElectronIDAlgo.

Definition at line 43 of file CutBasedElectronID.cc.

References classify(), EcalClusterLazyTools::covariances(), EgammaValidation_Zee_cff::cut, cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::eSuperClusterOverP(), eta, ElectronIDAlgo::getClusterShape(), edm::ParameterSet::getParameter(), reco::GsfElectron::hadronicOverEm(), reco::Particle::p4(), quality_, funct::sqrt(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackMomentumAtVtx(), and reco::GsfElectron::trackMomentumOut().

00045                                                            { 
00046   
00047   double eta = fabs(electron->p4().Eta());
00048   //EcalClusterTools shapeRef = getClusterShape(electron, e);
00049   
00050   double eOverP = electron->eSuperClusterOverP();
00051   double eSeed = electron->superCluster()->seed()->energy();
00052   double pin  = electron->trackMomentumAtVtx().R();   
00053   double eSeedOverPin = eSeed/pin; 
00054   double pout = electron->trackMomentumOut().R(); 
00055   double fBrem = (pin-pout)/pin;
00056   
00057   double hOverE = electron->hadronicOverEm();
00058   EcalClusterLazyTools lazyTools = getClusterShape(e,es);
00059   std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
00060   double sigmaee = sqrt(vCov[0]);
00061   double deltaPhiIn = electron->deltaPhiSuperClusterTrackAtVtx();
00062   double deltaEtaIn = electron->deltaEtaSuperClusterTrackAtVtx();
00063   
00064   int eb;
00065   if (eta < 1.479) 
00066     eb = 0;
00067   else {
00068     eb = 1; 
00069     sigmaee = sigmaee - 0.02*(fabs(eta) - 2.3);   //correct sigmaetaeta dependence on eta in endcap
00070   }
00071 
00072   std::vector<double> cut;
00073     
00074   // ROBUST Selection
00075   if (quality_ == "robust") {
00076 
00077     // hoe, sigmaEtaEta, dPhiIn, dEtaIn
00078     if (eta < 1.479)
00079       cut = cuts_.getParameter<std::vector<double> >("barrel");
00080     else
00081       cut = cuts_.getParameter<std::vector<double> >("endcap");
00082 
00083     if (hOverE > cut[0]) 
00084       return 0.;    
00085 
00086     if (sigmaee > cut[1]) 
00087       return 0.;    
00088 
00089     if (fabs(deltaPhiIn) > cut[2]) 
00090       return 0.;    
00091 
00092     if (fabs(deltaEtaIn) > cut[3]) 
00093       return 0.;    
00094     
00095     return 1.;
00096   }
00097   
00098   int cat = classify(electron);
00099 
00100   // TIGHT Selection
00101   if (quality_ == "tight") {
00102     
00103     if ((eOverP < 0.8) && (fBrem < 0.2)) 
00104       return 0.;
00105     
00106     if (eOverP < 0.9*(1-fBrem))
00107       return 0.;
00108     
00109     cut = cuts_.getParameter<std::vector<double> >("hOverE");
00110     if (hOverE > cut[cat+4*eb]) 
00111       return 0.;    
00112     
00113     cut = cuts_.getParameter<std::vector<double> >("sigmaEtaEta");
00114     if (sigmaee > cut[cat+4*eb]) 
00115       return 0.;    
00116     
00117     cut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
00118     if (eOverP < 1.5) {
00119       if (fabs(deltaPhiIn) > cut[cat+4*eb]) 
00120         return 0.;    
00121     } else {
00122       if (fabs(deltaPhiIn) > cut[3+4*eb])
00123         return 0.;
00124     }
00125     
00126     cut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
00127     if (fabs(deltaEtaIn) > cut[cat+4*eb]) 
00128       return 0.;    
00129     
00130     cut = cuts_.getParameter<std::vector<double> >("eSeedOverPinMin");
00131     if (eSeedOverPin < cut[cat+4*eb]) 
00132       return 0.;  
00133 
00134     cut = cuts_.getParameter<std::vector<double> >("eSeedOverPinMax");
00135     if (eSeedOverPin > cut[cat+4*eb]) 
00136       return 0.;  
00137 
00138     return 1.;
00139   }
00140   
00141     // LOOSE Selection
00142   if (quality_ == "loose") {
00143     
00144     if ((eOverP < 0.8) && (fBrem < 0.2)) 
00145       return 0.;
00146     
00147     cut = cuts_.getParameter<std::vector<double> >("hOverE");
00148     if (hOverE > cut[cat+4*eb]) 
00149       return 0.;    
00150     
00151     cut = cuts_.getParameter<std::vector<double> >("sigmaEtaEta");
00152     if (sigmaee > cut[cat+4*eb]) 
00153       return 0.;    
00154     
00155     cut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
00156     if (eOverP < 1.5) {
00157       if (fabs(deltaPhiIn) > cut[cat+4*eb]) 
00158         return 0.;    
00159     } else {
00160       if (fabs(deltaPhiIn) > cut[3+4*eb])
00161         return 0.;
00162     }
00163     
00164     cut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
00165     if (fabs(deltaEtaIn) > cut[cat+4*eb]) 
00166       return 0.;    
00167     
00168     cut = cuts_.getParameter<std::vector<double> >("eSeedOverPin");
00169     if (eSeedOverPin < cut[cat+4*eb]) 
00170       return 0.;  
00171     
00172     return 1.; 
00173   }
00174   
00175   return 0.;
00176 }

void CutBasedElectronID::setup ( const edm::ParameterSet conf  )  [virtual]

Reimplemented from ElectronIDAlgo.

Definition at line 4 of file CutBasedElectronID.cc.

References ElectronIDAlgo::baseSetup(), cuts_, cmsRelvalreport::exit, edm::ParameterSet::getParameter(), and quality_.

00004                                                           {
00005   
00006   // Get all the parameters
00007   baseSetup(conf);
00008   
00009   quality_ = conf.getParameter<std::string>("electronQuality");
00010   
00011   if (quality_ == "tight") {
00012     cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
00013   } else if (quality_=="robust") {
00014     cuts_ = conf.getParameter<edm::ParameterSet>("robustEleIDCuts");
00015   } else if (quality_=="loose") {
00016     cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
00017   } else {
00018     edm::LogError("CutBasedElectronID") << "Invalid electronQuality parameter: must be tight, loose or robust." ;
00019     exit (1);
00020   }
00021 }


Member Data Documentation

edm::ParameterSet CutBasedElectronID::cuts_ [private]

Definition at line 21 of file CutBasedElectronID.h.

Referenced by result(), and setup().

std::string CutBasedElectronID::quality_ [private]

Definition at line 20 of file CutBasedElectronID.h.

Referenced by result(), and setup().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:17:36 2009 for CMSSW by  doxygen 1.5.4