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 "RecoEgamma/EgammaTools/interface/EleEnergyRetriever.h"
0004 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0005
0006 class GsfEleEmHadD1IsoRhoCut : public CutApplicatorWithEventContentBase {
0007 public:
0008 GsfEleEmHadD1IsoRhoCut(const edm::ParameterSet& c);
0009
0010 result_type operator()(const reco::GsfElectronPtr&) const final;
0011
0012 void setConsumes(edm::ConsumesCollector&) final;
0013 void getEventContent(const edm::EventBase&) final;
0014
0015 double value(const reco::CandidatePtr& cand) const final;
0016
0017 CandidateType candidateType() const final { return ELECTRON; }
0018
0019 private:
0020 float rhoConstant_;
0021 EBEECutValues slopeTerm_;
0022 EBEECutValues slopeStart_;
0023 EBEECutValues constTerm_;
0024 EleEnergyRetriever energyRetriever_;
0025
0026 edm::Handle<double> rhoHandle_;
0027 };
0028
0029 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleEmHadD1IsoRhoCut, "GsfEleEmHadD1IsoRhoCut");
0030
0031 GsfEleEmHadD1IsoRhoCut::GsfEleEmHadD1IsoRhoCut(const edm::ParameterSet& params)
0032 : CutApplicatorWithEventContentBase(params),
0033 rhoConstant_(params.getParameter<double>("rhoConstant")),
0034 slopeTerm_(params, "slopeTerm"),
0035 slopeStart_(params, "slopeStart"),
0036 constTerm_(params, "constTerm"),
0037 energyRetriever_(params.getParameter<std::string>("energyType")) {
0038 edm::InputTag rhoTag = params.getParameter<edm::InputTag>("rho");
0039 contentTags_.emplace("rho", rhoTag);
0040 }
0041
0042 void GsfEleEmHadD1IsoRhoCut::setConsumes(edm::ConsumesCollector& cc) {
0043 auto rho = cc.consumes<double>(contentTags_["rho"]);
0044 contentTokens_.emplace("rho", rho);
0045 }
0046
0047 void GsfEleEmHadD1IsoRhoCut::getEventContent(const edm::EventBase& ev) {
0048 ev.getByLabel(contentTags_["rho"], rhoHandle_);
0049 }
0050
0051 CutApplicatorBase::result_type GsfEleEmHadD1IsoRhoCut::operator()(const reco::GsfElectronPtr& cand) const {
0052 const double rho = (*rhoHandle_);
0053
0054 const float isolEmHadDepth1 = cand->dr03EcalRecHitSumEt() + cand->dr03HcalTowerSumEt(1);
0055
0056 const float sinTheta = cand->p() != 0. ? cand->pt() / cand->p() : 0.;
0057 const float et = energyRetriever_(*cand) * sinTheta;
0058
0059 const float cutValue =
0060 et > slopeStart_(cand) ? slopeTerm_(cand) * (et - slopeStart_(cand)) + constTerm_(cand) : constTerm_(cand);
0061 return isolEmHadDepth1 < cutValue + rhoConstant_ * rho;
0062 }
0063
0064 double GsfEleEmHadD1IsoRhoCut::value(const reco::CandidatePtr& cand) const {
0065 reco::GsfElectronPtr ele(cand);
0066 return ele->dr03EcalRecHitSumEt() + ele->dr03HcalTowerSumEt(1);
0067 }