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 "DataFormats/Common/interface/ValueMap.h"
0004 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0005 
0006 class GsfEleFull5x5E2x5OverE5x5WithSatCut : public CutApplicatorBase {
0007 public:
0008   GsfEleFull5x5E2x5OverE5x5WithSatCut(const edm::ParameterSet& c);
0009 
0010   result_type operator()(const reco::GsfElectronPtr&) const final;
0011 
0012   double value(const reco::CandidatePtr& cand) const final;
0013 
0014   CandidateType candidateType() const final { return ELECTRON; }
0015 
0016 private:
0017   EBEECutValues minE1x5OverE5x5Cut_;
0018   EBEECutValues minE2x5OverE5x5Cut_;
0019   EBEECutValuesInt maxNrSatCrysIn5x5Cut_;
0020 };
0021 
0022 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleFull5x5E2x5OverE5x5WithSatCut, "GsfEleFull5x5E2x5OverE5x5WithSatCut");
0023 
0024 GsfEleFull5x5E2x5OverE5x5WithSatCut::GsfEleFull5x5E2x5OverE5x5WithSatCut(const edm::ParameterSet& params)
0025     : CutApplicatorBase(params),
0026       minE1x5OverE5x5Cut_(params, "minE1x5OverE5x5"),
0027       minE2x5OverE5x5Cut_(params, "minE2x5OverE5x5"),
0028       maxNrSatCrysIn5x5Cut_(params, "maxNrSatCrysIn5x5") {}
0029 
0030 CutApplicatorBase::result_type GsfEleFull5x5E2x5OverE5x5WithSatCut::operator()(const reco::GsfElectronPtr& cand) const {
0031   if (cand->nSaturatedXtals() > maxNrSatCrysIn5x5Cut_(cand))
0032     return true;
0033 
0034   const double e5x5 = cand->full5x5_e5x5();
0035   const double e2x5OverE5x5 = e5x5 != 0 ? cand->full5x5_e2x5Max() / e5x5 : 0;
0036   const double e1x5OverE5x5 = e5x5 != 0 ? cand->full5x5_e1x5() / e5x5 : 0;
0037 
0038   return e1x5OverE5x5 > minE1x5OverE5x5Cut_(cand) || e2x5OverE5x5 > minE2x5OverE5x5Cut_(cand);
0039 }
0040 
0041 double GsfEleFull5x5E2x5OverE5x5WithSatCut::value(const reco::CandidatePtr& cand) const {
0042   reco::GsfElectronPtr ele(cand);
0043   //btw we broke somebodies nice model of assuming every cut is 1D....
0044   //what this is returning is fairly meaningless...
0045   return ele->full5x5_e1x5() ? ele->full5x5_e2x5Max() / ele->full5x5_e1x5() : 0.;
0046 }