CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes
ZElectronsSelector Class Reference

Public Member Functions

const float getEffectiveArea (float eta) const
 
void newEvent (const edm::Event &, const edm::EventSetup &)
 
bool operator() (const reco::GsfElectron &) const
 
void printEffectiveAreas () const
 
 ZElectronsSelector (const edm::ParameterSet &, edm::ConsumesCollector &iC)
 

Public Attributes

edm::Handle< double > _rhoHandle
 
std::vector< double > absEtaMax_
 
std::vector< double > absEtaMin_
 
vector< double > dEtaInSeedCut
 
vector< double > dPhiInCut
 
std::vector< double > effectiveAreaValues_
 
vector< double > EInvMinusPInv
 
edm::ParameterSet eleIDWP
 
vector< double > hOverECut
 
vector< int > missHits
 
vector< double > relCombIso
 
vector< double > sigmaIEtaIEtaCut
 
edm::EDGetTokenT< reco::GsfElectronCollectiontheGsfEToken
 
edm::EDGetTokenT< double > theRhoToken
 

Detailed Description

Definition at line 44 of file ZElectronsSelector.cc.

Constructor & Destructor Documentation

ZElectronsSelector::ZElectronsSelector ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 103 of file ZElectronsSelector.cc.

References absEtaMax_, absEtaMin_, dEtaInSeedCut, dPhiInCut, effectiveAreaValues_, EInvMinusPInv, eleIDWP, edm::ParameterSet::getParameter(), hOverECut, missHits, relCombIso, and sigmaIEtaIEtaCut.

103  :
104  theRhoToken(iC.consumes <double> (cfg.getParameter<edm::InputTag>("rho")))
105 {
106  absEtaMin_ = cfg.getParameter<std::vector<double> >("absEtaMin");
107  absEtaMax_ = cfg.getParameter<std::vector<double> >("absEtaMax");
108  effectiveAreaValues_= cfg.getParameter<std::vector<double> >("effectiveAreaValues");
109  //printEffectiveAreas();
110  eleIDWP = cfg.getParameter<edm::ParameterSet>("eleID");
111 
112  missHits = eleIDWP.getParameter<std::vector<int> >("missingHitsCut");
113  sigmaIEtaIEtaCut= eleIDWP.getParameter<std::vector<double> >("full5x5_sigmaIEtaIEtaCut");
114  dEtaInSeedCut = eleIDWP.getParameter<std::vector<double> >("dEtaInSeedCut");
115  dPhiInCut = eleIDWP.getParameter<std::vector<double> >("dPhiInCut");
116  hOverECut = eleIDWP.getParameter<std::vector<double> >("hOverECut");
117  relCombIso = eleIDWP.getParameter<std::vector<double> >("relCombIsolationWithEACut");
118  EInvMinusPInv = eleIDWP.getParameter<std::vector<double> >("EInverseMinusPInverseCut");
119 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
vector< double > EInvMinusPInv
std::vector< double > absEtaMax_
vector< double > dEtaInSeedCut
std::vector< double > effectiveAreaValues_
edm::EDGetTokenT< double > theRhoToken
vector< double > relCombIso
edm::ParameterSet eleIDWP
vector< double > dPhiInCut
std::vector< double > absEtaMin_
vector< double > hOverECut
vector< double > sigmaIEtaIEtaCut

Member Function Documentation

const float ZElectronsSelector::getEffectiveArea ( float  eta) const

Definition at line 87 of file ZElectronsSelector.cc.

References funct::abs(), fftjetcommon_cfi::nEtaBins, and parallelization::uint().

Referenced by operator()().

87  {
88 
89  float effArea = 0;
90  uint nEtaBins = absEtaMin_.size();
91  for(uint iEta = 0; iEta<nEtaBins; iEta++){
92  if( std::abs(eta) >= absEtaMin_[iEta]
93  && std::abs(eta) < absEtaMax_[iEta] ){
94  effArea = effectiveAreaValues_[iEta];
95  break;
96  }
97  }
98 
99  return effArea;
100 }
std::vector< double > absEtaMax_
std::vector< double > effectiveAreaValues_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
def uint(string)
std::vector< double > absEtaMin_
void ZElectronsSelector::newEvent ( const edm::Event ev,
const edm::EventSetup  
)

Definition at line 123 of file ZElectronsSelector.cc.

References _rhoHandle, edm::Event::getByToken(), and theRhoToken.

123  {
124 
126 
127 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
edm::EDGetTokenT< double > theRhoToken
edm::Handle< double > _rhoHandle
bool ZElectronsSelector::operator() ( const reco::GsfElectron el) const

Definition at line 129 of file ZElectronsSelector.cc.

References _rhoHandle, funct::abs(), reco::GsfElectron::deltaEtaSeedClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), dEtaInSeedCut, dPhiInCut, reco::GsfElectron::ecalEnergy(), EInvMinusPInv, reco::GsfElectron::eSuperClusterOverP(), objects.autophobj::float, reco::GsfElectron::full5x5_sigmaIetaIeta(), edm::Ref< C, T, F >::get(), getEffectiveArea(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hadronicOverEm(), hOverECut, reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), edm::HandleBase::isValid(), SiStripPI::max, missHits, reco::HitPattern::MISSING_INNER_HITS, reco::GsfElectron::pfIsolationVariables(), position, edm::Handle< T >::product(), reco::LeafCandidate::pt(), relCombIso, rho, sigmaIEtaIEtaCut, reco::GsfElectron::PflowIsolationVariables::sumChargedHadronPt, reco::GsfElectron::PflowIsolationVariables::sumNeutralHadronEt, reco::GsfElectron::PflowIsolationVariables::sumPhotonEt, and reco::GsfElectron::superCluster().

129  {
130 
131  float pt_e = el.pt();
132  unsigned int ind=0;
133  auto etrack = el.gsfTrack();
134  float abseta = fabs((el.superCluster().get())->position().eta());
135 
136  if (el.isEB()){
137  if (abseta> 1.479) return false; // check if it is really needed
138  }
139  if (el.isEE()){
140  ind=1;
141  if (abseta< 1.479) return false; // check if it is really needed
142  if (abseta>=2.5 ) return false; // check if it is really needed
143  }
144 
145  if (etrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS)>missHits[ind]) return false;
146  if (el.full5x5_sigmaIetaIeta() > sigmaIEtaIEtaCut[ind]) return false;
147  if (fabs(el.deltaPhiSuperClusterTrackAtVtx())> dPhiInCut[ind]) return false;
148  if (fabs(el.deltaEtaSeedClusterTrackAtVtx())> dEtaInSeedCut[ind]) return false;
149  if (el.hadronicOverEm()>hOverECut[ind]) return false;
150  const float eA = getEffectiveArea( abseta );
151  const float rho = _rhoHandle.isValid() ? (float)(*_rhoHandle.product()) : 0;
155  ) > relCombIso[ind]*pt_e) return false;
156  const float ecal_energy_inverse = 1.0/el.ecalEnergy();
157  const float eSCoverP = el.eSuperClusterOverP();
158  if ( std::abs(1.0 - eSCoverP)*ecal_energy_inverse > EInvMinusPInv[ind]) return false;
159 
160  return true;
161 }
const PflowIsolationVariables & pfIsolationVariables() const
Definition: GsfElectron.h:677
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:185
float eSuperClusterOverP() const
Definition: GsfElectron.h:245
vector< double > EInvMinusPInv
vector< double > dEtaInSeedCut
double pt() const final
transverse momentum
bool isEE() const
Definition: GsfElectron.h:353
bool isEB() const
Definition: GsfElectron.h:352
float full5x5_sigmaIetaIeta() const
Definition: GsfElectron.h:458
float hadronicOverEm() const
Definition: GsfElectron.h:491
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:631
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:252
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:630
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:243
bool isValid() const
Definition: HandleBase.h:74
float deltaEtaSeedClusterTrackAtVtx() const
Definition: GsfElectron.h:255
vector< double > relCombIso
T const * product() const
Definition: Handle.h:81
edm::Handle< double > _rhoHandle
float ecalEnergy() const
Definition: GsfElectron.h:848
static int position[264][3]
Definition: ReadPGInfo.cc:509
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:184
vector< double > dPhiInCut
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:629
const float getEffectiveArea(float eta) const
vector< double > hOverECut
vector< double > sigmaIEtaIEtaCut
void ZElectronsSelector::printEffectiveAreas ( ) const

Definition at line 75 of file ZElectronsSelector.cc.

References fftjetcommon_cfi::nEtaBins, and parallelization::uint().

75  {
76 
77 
78  printf(" eta_min eta_max effective area\n");
79  uint nEtaBins = absEtaMin_.size();
80  for(uint iEta = 0; iEta<nEtaBins; iEta++){
81  printf(" %8.4f %8.4f %8.5f\n",
82  absEtaMin_[iEta], absEtaMax_[iEta],
83  effectiveAreaValues_[iEta]);
84  }
85 
86 }
std::vector< double > absEtaMax_
std::vector< double > effectiveAreaValues_
def uint(string)
std::vector< double > absEtaMin_

Member Data Documentation

edm::Handle< double > ZElectronsSelector::_rhoHandle

Definition at line 55 of file ZElectronsSelector.cc.

Referenced by newEvent(), and operator()().

std::vector<double> ZElectronsSelector::absEtaMax_

Definition at line 58 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

std::vector<double> ZElectronsSelector::absEtaMin_

Definition at line 57 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::dEtaInSeedCut

Definition at line 65 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::dPhiInCut

Definition at line 66 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

std::vector<double> ZElectronsSelector::effectiveAreaValues_

Definition at line 59 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::EInvMinusPInv

Definition at line 69 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

edm::ParameterSet ZElectronsSelector::eleIDWP

Definition at line 61 of file ZElectronsSelector.cc.

Referenced by ZElectronsSelector().

vector<double> ZElectronsSelector::hOverECut

Definition at line 67 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<int> ZElectronsSelector::missHits

Definition at line 63 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::relCombIso

Definition at line 68 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

vector<double> ZElectronsSelector::sigmaIEtaIEtaCut

Definition at line 64 of file ZElectronsSelector.cc.

Referenced by operator()(), and ZElectronsSelector().

edm::EDGetTokenT<reco::GsfElectronCollection> ZElectronsSelector::theGsfEToken

Definition at line 54 of file ZElectronsSelector.cc.

edm::EDGetTokenT<double> ZElectronsSelector::theRhoToken

Definition at line 53 of file ZElectronsSelector.cc.

Referenced by newEvent().