Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:06

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0004 
0005 class GsfEleHadronicOverEMLinearCut : public CutApplicatorBase {
0006 public:
0007   GsfEleHadronicOverEMLinearCut(const edm::ParameterSet& params)
0008       : CutApplicatorBase(params),
0009         slopeTerm_(params, "slopeTerm"),
0010         slopeStart_(params, "slopeStart"),
0011         constTerm_(params, "constTerm") {}
0012 
0013   result_type operator()(const reco::GsfElectronPtr&) const final;
0014 
0015   double value(const reco::CandidatePtr& cand) const final;
0016 
0017   CandidateType candidateType() const final { return ELECTRON; }
0018 
0019 private:
0020   EBEECutValues slopeTerm_;
0021   EBEECutValues slopeStart_;
0022   EBEECutValues constTerm_;
0023 };
0024 
0025 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleHadronicOverEMLinearCut, "GsfEleHadronicOverEMLinearCut");
0026 
0027 CutApplicatorBase::result_type GsfEleHadronicOverEMLinearCut::operator()(const reco::GsfElectronPtr& cand) const {
0028   const float energy = cand->superCluster()->energy();
0029   const float cutValue = energy > slopeStart_(cand) ? slopeTerm_(cand) * (energy - slopeStart_(cand)) + constTerm_(cand)
0030                                                     : constTerm_(cand);
0031 
0032   return cand->hadronicOverEm() * energy < cutValue;
0033 }
0034 
0035 double GsfEleHadronicOverEMLinearCut::value(const reco::CandidatePtr& cand) const {
0036   reco::GsfElectronPtr ele(cand);
0037   const float energy = ele->superCluster()->energy();
0038   return ele->hadronicOverEm() * energy;
0039 }