CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoEgamma/ElectronIdentification/src/PTDRElectronID.cc

Go to the documentation of this file.
00001 #include "RecoEgamma/ElectronIdentification/interface/PTDRElectronID.h"
00002 
00003 void PTDRElectronID::setup(const edm::ParameterSet& conf) {
00004 
00005   // Get all the parameters
00006   baseSetup(conf);
00007   
00008   quality_ =  conf.getParameter<std::string>("electronQuality");
00009   
00010   useEoverPIn_ = conf.getParameter<std::vector<int> >("useEoverPIn");
00011   useDeltaEtaIn_ = conf.getParameter<std::vector<int> >("useDeltaEtaIn");
00012   useDeltaPhiIn_ = conf.getParameter<std::vector<int> >("useDeltaPhiIn");
00013   useHoverE_ = conf.getParameter<std::vector<int> >("useHoverE");
00014   useE9overE25_ = conf.getParameter<std::vector<int> >("useE9overE25");
00015   useEoverPOut_ = conf.getParameter<std::vector<int> >("useEoverPOut");
00016   useDeltaPhiOut_ = conf.getParameter<std::vector<int> >("useDeltaPhiOut");
00017   useInvEMinusInvP_ = conf.getParameter<std::vector<int> >("useInvEMinusInvP");
00018   useBremFraction_ = conf.getParameter<std::vector<int> >("useBremFraction");
00019   useSigmaEtaEta_ = conf.getParameter<std::vector<int> >("useSigmaEtaEta");
00020   useSigmaPhiPhi_ = conf.getParameter<std::vector<int> >("useSigmaPhiPhi");
00021   acceptCracks_ = conf.getParameter<std::vector<int> >("acceptCracks");
00022   
00023   if (quality_=="tight") {
00024     cuts_ = conf.getParameter<edm::ParameterSet>("tightEleIDCuts");
00025     variables_ = 2 ;
00026   } else if (quality_=="medium") {
00027     cuts_ = conf.getParameter<edm::ParameterSet>("mediumEleIDCuts");
00028     variables_ = 1 ;
00029   } else if (quality_=="loose") {
00030     cuts_ = conf.getParameter<edm::ParameterSet>("looseEleIDCuts");
00031     variables_ = 0 ;
00032   } else {
00033      throw cms::Exception("Configuration") << "Invalid electronQuality parameter in PTDElectronID: must be tight, medium or loose." ;
00034   }  
00035 }
00036 
00037 double PTDRElectronID::result(const reco::GsfElectron* electron,
00038                               const edm::Event& e ,
00039                               const edm::EventSetup& es) {
00040 
00041   //determine which element of the cut arrays in cfi file to read
00042   //depending on the electron classification
00043   int icut=0;
00044   int elClass = electron->classification() ;
00045   if (electron->isEB()) //barrel
00046      {
00047        if (elClass == reco::GsfElectron::GOLDEN)    icut=0;
00048        if (elClass == reco::GsfElectron::BIGBREM)   icut=1;
00049        //if (elClass == reco::GsfElectron::NARROW)    icut=2;
00050        if (elClass == reco::GsfElectron::SHOWERING) icut=3;
00051        if (elClass == reco::GsfElectron::GAP)       icut=8;
00052      }
00053   if (electron->isEE()) //endcap
00054      {
00055        if (elClass == reco::GsfElectron::GOLDEN)    icut=4;
00056        if (elClass == reco::GsfElectron::BIGBREM)   icut=5;
00057        //if (elClass == reco::GsfElectron::NARROW)    icut=6;
00058        if (elClass == reco::GsfElectron::SHOWERING) icut=7;
00059        if (elClass == reco::GsfElectron::GAP)       icut=8;
00060      }
00061   if (elClass == reco::GsfElectron::UNKNOWN) 
00062      {
00063        edm::LogError("PTDRElectronID") << "Error: unrecognized electron classification ";
00064        return 1.;
00065      }
00066 
00067   if (acceptCracks_[variables_])
00068     if (elClass == reco::GsfElectron::GAP) return 1.;
00069   
00070   if (useEoverPIn_[variables_]) {
00071     double value = electron->eSuperClusterOverP();
00072     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
00073     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
00074     if (value<mincut[icut] || value>maxcut[icut]) return 0.;
00075   }
00076 
00077   if (useDeltaEtaIn_[variables_]) {
00078     double value = electron->deltaEtaSuperClusterTrackAtVtx();
00079     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
00080     if (fabs(value)>maxcut[icut]) return 0.;
00081   }
00082 
00083   if (useDeltaPhiIn_[variables_]) {
00084     double value = electron->deltaPhiSuperClusterTrackAtVtx();
00085     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
00086     if (fabs(value)>maxcut[icut]) return 0.;
00087   }
00088 
00089   if (useHoverE_[variables_]) {
00090     double value = electron->hadronicOverEm();
00091     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
00092     if (value>maxcut[icut]) return 0.;
00093   }
00094 
00095   if (useEoverPOut_[variables_]) {
00096     double value = electron->eSeedClusterOverPout();
00097     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
00098     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
00099     if (value<mincut[icut] || value>maxcut[icut]) return 0.;
00100   }
00101 
00102   if (useDeltaPhiOut_[variables_]) {
00103     double value = electron->deltaPhiSeedClusterTrackAtCalo();
00104     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
00105     if (fabs(value)>maxcut[icut]) return 0.;
00106   }
00107 
00108   if (useInvEMinusInvP_[variables_]) {
00109     double value = (1./electron->caloEnergy())-(1./electron->trackMomentumAtVtx().R());
00110     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
00111     if (value>maxcut[icut]) return 0.;
00112   }
00113 
00114   if (useBremFraction_[variables_]) {
00115     double value = electron->trackMomentumAtVtx().R()-electron->trackMomentumOut().R();
00116     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
00117     if (value<mincut[icut]) return 0.;
00118   }
00119 
00120   EcalClusterLazyTools lazyTools = getClusterShape(e,es);
00121   std::vector<float> vCov = lazyTools.localCovariances(*(electron->superCluster()->seed())) ;
00122   //std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
00123     
00124   if (useE9overE25_[variables_]) {
00125     double value = lazyTools.e3x3(*(electron->superCluster()->seed()))/lazyTools.e5x5(*(electron->superCluster()->seed()));
00126     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
00127     if (fabs(value)<mincut[icut]) return 0.;
00128   }
00129 
00130   if (useSigmaEtaEta_[variables_]) {
00131     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
00132     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
00133     if (sqrt(vCov[0])<mincut[icut] || sqrt(vCov[0])>maxcut[icut]) return 0.;
00134   }
00135 
00136   if (useSigmaPhiPhi_[variables_]) {
00137     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
00138     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
00139     if (sqrt(vCov[1])<mincut[icut] || sqrt(vCov[1])>maxcut[icut]) return 0.;
00140   }
00141 
00142   return 1.;
00143 
00144 }