CMS 3D CMS Logo

PTDRElectronID Class Reference

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

Inheritance diagram for PTDRElectronID:

ElectronIDAlgo

List of all members.

Public Member Functions

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

Private Attributes

std::vector< intacceptCracks_
edm::ParameterSet cuts_
std::string quality_
std::vector< intuseBremFraction_
std::vector< intuseDeltaEtaIn_
std::vector< intuseDeltaPhiIn_
std::vector< intuseDeltaPhiOut_
std::vector< intuseE9overE25_
std::vector< intuseEoverPIn_
std::vector< intuseEoverPOut_
std::vector< intuseHoverE_
std::vector< intuseInvEMinusInvP_
std::vector< intuseSigmaEtaEta_
std::vector< intuseSigmaPhiPhi_
int variables_


Detailed Description

Definition at line 6 of file PTDRElectronID.h.


Constructor & Destructor Documentation

PTDRElectronID::PTDRElectronID (  )  [inline]

Definition at line 10 of file PTDRElectronID.h.

00010 {};

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

Definition at line 12 of file PTDRElectronID.h.

00012 {};


Member Function Documentation

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

Reimplemented from ElectronIDAlgo.

Definition at line 39 of file PTDRElectronID.cc.

References acceptCracks_, reco::GsfElectron::caloEnergy(), reco::GsfElectron::classification(), EcalClusterLazyTools::covariances(), cuts_, reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSeedClusterTrackAtCalo(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), EcalClusterLazyTools::e3x3(), EcalClusterLazyTools::e5x5(), reco::GsfElectron::eSeedClusterOverPout(), reco::GsfElectron::eSuperClusterOverP(), ElectronIDAlgo::getClusterShape(), edm::ParameterSet::getParameter(), reco::GsfElectron::hadronicOverEm(), funct::sqrt(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::trackMomentumOut(), useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, value, and variables_.

00041                                                        {
00042 
00043   //determine which element of the cut arrays in electronId.cfi to read
00044   //depending on the electron classification
00045   int icut=0;
00046   switch (electron->classification()) {
00047   case 0:   icut=0; break;
00048   case 10:  icut=1; break;
00049   case 20:  icut=2; break;
00050   case 30:  icut=3; break;
00051   case 31:  icut=3; break;
00052   case 32:  icut=3; break;
00053   case 33:  icut=3; break;
00054   case 34:  icut=3; break;
00055   case 40:  icut=8; break;
00056   case 100: icut=4; break;
00057   case 110: icut=5; break;
00058   case 120: icut=6; break;
00059   case 130: icut=7; break;
00060   case 131: icut=7; break;
00061   case 132: icut=7; break;
00062   case 133: icut=7; break;
00063   case 134: icut=7; break;
00064 
00065   default:
00066     edm::LogError("CutBasedElectronID") << "Error: unrecognized electron classification ";
00067     break;
00068   }
00069 
00070   if (acceptCracks_[variables_])
00071     if (electron->classification()==40) return 1.;
00072   
00073   if (useEoverPIn_[variables_]) {
00074     double value = electron->eSuperClusterOverP();
00075     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPInMax");
00076     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPInMin");
00077     if (value<mincut[icut] || value>maxcut[icut]) return 0.;
00078   }
00079 
00080   if (useDeltaEtaIn_[variables_]) {
00081     double value = electron->deltaEtaSuperClusterTrackAtVtx();
00082     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaEtaIn");
00083     if (fabs(value)>maxcut[icut]) return 0.;
00084   }
00085 
00086   if (useDeltaPhiIn_[variables_]) {
00087     double value = electron->deltaPhiSuperClusterTrackAtVtx();
00088     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiIn");
00089     if (fabs(value)>maxcut[icut]) return 0.;
00090   }
00091 
00092   if (useHoverE_[variables_]) {
00093     double value = electron->hadronicOverEm();
00094     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("HoverE");
00095     if (fabs(value)>maxcut[icut]) return 0.;
00096   }
00097 
00098   if (useEoverPOut_[variables_]) {
00099     double value = electron->eSeedClusterOverPout();
00100     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("EoverPOutMax");
00101     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("EoverPOutMin");
00102     if (value<mincut[icut] || value>maxcut[icut]) return 0.;
00103   }
00104 
00105   if (useDeltaPhiOut_[variables_]) {
00106     double value = electron->deltaPhiSeedClusterTrackAtCalo();
00107     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("deltaPhiOut");
00108     if (fabs(value)>maxcut[icut]) return 0.;
00109   }
00110 
00111   if (useInvEMinusInvP_[variables_]) {
00112     double value = (1./electron->caloEnergy())-(1./electron->trackMomentumAtVtx().R());
00113     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("invEMinusInvP");
00114     if (value>maxcut[icut]) return 0.;
00115   }
00116 
00117   if (useBremFraction_[variables_]) {
00118     double value = electron->trackMomentumAtVtx().R()-electron->trackMomentumOut().R();
00119     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("bremFraction");
00120     if (value<mincut[icut]) return 0.;
00121   }
00122 
00123   EcalClusterLazyTools lazyTools = getClusterShape(e,es);
00124   std::vector<float> vCov = lazyTools.covariances(*(electron->superCluster()->seed())) ;
00125   
00126   if (useE9overE25_[variables_]) {
00127     double value = lazyTools.e3x3(*(electron->superCluster()->seed()))/lazyTools.e5x5(*(electron->superCluster()->seed()));
00128     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("E9overE25");
00129     if (fabs(value)<mincut[icut]) return 0.;
00130   }
00131 
00132   if (useSigmaEtaEta_[variables_]) {
00133     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMax");
00134     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaEtaEtaMin");
00135     if (sqrt(vCov[0])<mincut[icut] || sqrt(vCov[0])>maxcut[icut]) return 0.;
00136   }
00137 
00138   if (useSigmaPhiPhi_[variables_]) {
00139     std::vector<double> mincut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMin");
00140     std::vector<double> maxcut = cuts_.getParameter<std::vector<double> >("sigmaPhiPhiMax");
00141     if (sqrt(vCov[1])<mincut[icut] || sqrt(vCov[1])>maxcut[icut]) return 0.;
00142   }
00143 
00144   return 1.;
00145 
00146 }

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

Reimplemented from ElectronIDAlgo.

Definition at line 3 of file PTDRElectronID.cc.

References acceptCracks_, ElectronIDAlgo::baseSetup(), cuts_, cmsRelvalreport::exit, edm::ParameterSet::getParameter(), quality_, useBremFraction_, useDeltaEtaIn_, useDeltaPhiIn_, useDeltaPhiOut_, useE9overE25_, useEoverPIn_, useEoverPOut_, useHoverE_, useInvEMinusInvP_, useSigmaEtaEta_, useSigmaPhiPhi_, and variables_.

00003                                                       {
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     edm::LogError("PTDRElectronID") << "Invalid electronQuality parameter: must be tight, medium or loose." ;
00034     exit (1);
00035   }
00036   
00037 }


Member Data Documentation

std::vector<int> PTDRElectronID::acceptCracks_ [private]

Definition at line 32 of file PTDRElectronID.h.

Referenced by result(), and setup().

edm::ParameterSet PTDRElectronID::cuts_ [private]

Definition at line 34 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::string PTDRElectronID::quality_ [private]

Definition at line 19 of file PTDRElectronID.h.

Referenced by setup().

std::vector<int> PTDRElectronID::useBremFraction_ [private]

Definition at line 29 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaEtaIn_ [private]

Definition at line 22 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaPhiIn_ [private]

Definition at line 23 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useDeltaPhiOut_ [private]

Definition at line 27 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useE9overE25_ [private]

Definition at line 25 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useEoverPIn_ [private]

Definition at line 21 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useEoverPOut_ [private]

Definition at line 26 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useHoverE_ [private]

Definition at line 24 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useInvEMinusInvP_ [private]

Definition at line 28 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useSigmaEtaEta_ [private]

Definition at line 30 of file PTDRElectronID.h.

Referenced by result(), and setup().

std::vector<int> PTDRElectronID::useSigmaPhiPhi_ [private]

Definition at line 31 of file PTDRElectronID.h.

Referenced by result(), and setup().

int PTDRElectronID::variables_ [private]

Definition at line 36 of file PTDRElectronID.h.

Referenced by result(), and setup().


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