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 
0004 class GsfEleHadronicOverEMEnergyScaledCut : public CutApplicatorWithEventContentBase {
0005 public:
0006   GsfEleHadronicOverEMEnergyScaledCut(const edm::ParameterSet& c)
0007       : CutApplicatorWithEventContentBase(c),
0008         barrelC0_(c.getParameter<double>("barrelC0")),
0009         barrelCE_(c.getParameter<double>("barrelCE")),
0010         barrelCr_(c.getParameter<double>("barrelCr")),
0011         endcapC0_(c.getParameter<double>("endcapC0")),
0012         endcapCE_(c.getParameter<double>("endcapCE")),
0013         endcapCr_(c.getParameter<double>("endcapCr")),
0014         barrelCutOff_(c.getParameter<double>("barrelCutOff")) {
0015     edm::InputTag rhoTag = c.getParameter<edm::InputTag>("rho");
0016     contentTags_.emplace("rho", rhoTag);
0017   }
0018 
0019   result_type operator()(const reco::GsfElectronPtr&) const final;
0020 
0021   void setConsumes(edm::ConsumesCollector&) final;
0022   void getEventContent(const edm::EventBase&) final;
0023 
0024   double value(const reco::CandidatePtr& cand) const final;
0025 
0026   CandidateType candidateType() const final { return ELECTRON; }
0027 
0028 private:
0029   const float barrelC0_, barrelCE_, barrelCr_, endcapC0_, endcapCE_, endcapCr_, barrelCutOff_;
0030   edm::Handle<double> rhoHandle_;
0031 };
0032 
0033 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleHadronicOverEMEnergyScaledCut, "GsfEleHadronicOverEMEnergyScaledCut");
0034 
0035 void GsfEleHadronicOverEMEnergyScaledCut::setConsumes(edm::ConsumesCollector& cc) {
0036   auto rho = cc.consumes<double>(contentTags_["rho"]);
0037   contentTokens_.emplace("rho", rho);
0038 }
0039 
0040 void GsfEleHadronicOverEMEnergyScaledCut::getEventContent(const edm::EventBase& ev) {
0041   ev.getByLabel(contentTags_["rho"], rhoHandle_);
0042 }
0043 
0044 CutApplicatorBase::result_type GsfEleHadronicOverEMEnergyScaledCut::operator()(const reco::GsfElectronPtr& cand) const {
0045   const double rho = rhoHandle_.isValid() ? (float)(*rhoHandle_) : 0;
0046   const float energy = cand->superCluster()->energy();
0047   const float c0 = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelC0_ : endcapC0_);
0048   const float cE = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCE_ : endcapCE_);
0049   const float cR = (std::abs(cand->superCluster()->position().eta()) < barrelCutOff_ ? barrelCr_ : endcapCr_);
0050   return cand->hadronicOverEm() < c0 + cE / energy + cR * rho / energy;
0051 }
0052 
0053 double GsfEleHadronicOverEMEnergyScaledCut::value(const reco::CandidatePtr& cand) const {
0054   reco::GsfElectronPtr ele(cand);
0055   return ele->hadronicOverEm();
0056 }