Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "CommonTools/Egamma/interface/EffectiveAreas.h"
0004 
0005 class GsfEleRelPFIsoScaledCut : public CutApplicatorWithEventContentBase {
0006 public:
0007   GsfEleRelPFIsoScaledCut(const edm::ParameterSet& c);
0008 
0009   result_type operator()(const reco::GsfElectronPtr&) const final;
0010 
0011   void setConsumes(edm::ConsumesCollector&) final;
0012   void getEventContent(const edm::EventBase&) final;
0013 
0014   double value(const reco::CandidatePtr& cand) const final;
0015 
0016   CandidateType candidateType() const final { return ELECTRON; }
0017 
0018 private:
0019   const float barrelC0_, endcapC0_, barrelCpt_, endcapCpt_, barrelCutOff_;
0020   EffectiveAreas effectiveAreas_;
0021   edm::Handle<double> rhoHandle_;
0022 };
0023 
0024 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleRelPFIsoScaledCut, "GsfEleRelPFIsoScaledCut");
0025 
0026 GsfEleRelPFIsoScaledCut::GsfEleRelPFIsoScaledCut(const edm::ParameterSet& c)
0027     : CutApplicatorWithEventContentBase(c),
0028       barrelC0_(c.getParameter<double>("barrelC0")),
0029       endcapC0_(c.getParameter<double>("endcapC0")),
0030       barrelCpt_(c.getParameter<double>("barrelCpt")),
0031       endcapCpt_(c.getParameter<double>("endcapCpt")),
0032       barrelCutOff_(c.getParameter<double>("barrelCutOff")),
0033       effectiveAreas_((c.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath()) {
0034   edm::InputTag rhoTag = c.getParameter<edm::InputTag>("rho");
0035   contentTags_.emplace("rho", rhoTag);
0036 }
0037 
0038 void GsfEleRelPFIsoScaledCut::setConsumes(edm::ConsumesCollector& cc) {
0039   auto rho = cc.consumes<double>(contentTags_["rho"]);
0040   contentTokens_.emplace("rho", rho);
0041 }
0042 
0043 void GsfEleRelPFIsoScaledCut::getEventContent(const edm::EventBase& ev) {
0044   ev.getByLabel(contentTags_["rho"], rhoHandle_);
0045 }
0046 
0047 CutApplicatorBase::result_type GsfEleRelPFIsoScaledCut::operator()(const reco::GsfElectronPtr& cand) const {
0048   // Establish the cut value
0049   double absEta = std::abs(cand->superCluster()->eta());
0050 
0051   const float C0 = (absEta < barrelCutOff_ ? barrelC0_ : endcapC0_);
0052   const float Cpt = (absEta < barrelCutOff_ ? barrelCpt_ : endcapCpt_);
0053   const float isoCut = C0 + Cpt / cand->pt();
0054 
0055   return value(cand) < isoCut;
0056 }
0057 
0058 double GsfEleRelPFIsoScaledCut::value(const reco::CandidatePtr& cand) const {
0059   // Establish the cut value
0060   reco::GsfElectronPtr ele(cand);
0061   double absEta = std::abs(ele->superCluster()->eta());
0062 
0063   // Compute the combined isolation with effective area correction
0064   auto pfIso = ele->pfIsolationVariables();
0065   const float chad = pfIso.sumChargedHadronPt;
0066   const float nhad = pfIso.sumNeutralHadronEt;
0067   const float pho = pfIso.sumPhotonEt;
0068   const float eA = effectiveAreas_.getEffectiveArea(absEta);
0069   const float rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0;  // std::max likes float arguments
0070   const float iso = chad + std::max(0.0f, nhad + pho - rho * eA);
0071   return iso / cand->pt();
0072 }