Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:17:57

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "DataFormats/EgammaCandidates/interface/ConversionFwd.h"
0004 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
0005 #include "CommonTools/Egamma/interface/ConversionTools.h"
0006 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0007 
0008 class GsfEleTrkPtIsoRhoCut : public CutApplicatorWithEventContentBase {
0009 public:
0010   GsfEleTrkPtIsoRhoCut(const edm::ParameterSet& c);
0011 
0012   result_type operator()(const reco::GsfElectronPtr&) const final;
0013 
0014   void setConsumes(edm::ConsumesCollector&) final;
0015   void getEventContent(const edm::EventBase&) final;
0016 
0017   double value(const reco::CandidatePtr& cand) const final;
0018 
0019   CandidateType candidateType() const final { return ELECTRON; }
0020 
0021 private:
0022   EBEECutValues slopeTerm_;
0023   EBEECutValues slopeStart_;
0024   EBEECutValues constTerm_;
0025   EBEECutValues rhoEtStart_;
0026   EBEECutValues rhoEA_;
0027 
0028   edm::Handle<double> rhoHandle_;
0029 };
0030 
0031 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleTrkPtIsoRhoCut, "GsfEleTrkPtIsoRhoCut");
0032 
0033 GsfEleTrkPtIsoRhoCut::GsfEleTrkPtIsoRhoCut(const edm::ParameterSet& params)
0034     : CutApplicatorWithEventContentBase(params),
0035       slopeTerm_(params, "slopeTerm"),
0036       slopeStart_(params, "slopeStart"),
0037       constTerm_(params, "constTerm"),
0038       rhoEtStart_(params, "rhoEtStart"),
0039       rhoEA_(params, "rhoEA") {
0040   edm::InputTag rhoTag = params.getParameter<edm::InputTag>("rho");
0041   contentTags_.emplace("rho", rhoTag);
0042 }
0043 
0044 void GsfEleTrkPtIsoRhoCut::setConsumes(edm::ConsumesCollector& cc) {
0045   auto rho = cc.consumes<double>(contentTags_["rho"]);
0046   contentTokens_.emplace("rho", rho);
0047 }
0048 
0049 void GsfEleTrkPtIsoRhoCut::getEventContent(const edm::EventBase& ev) { ev.getByLabel(contentTags_["rho"], rhoHandle_); }
0050 
0051 CutApplicatorBase::result_type GsfEleTrkPtIsoRhoCut::operator()(const reco::GsfElectronPtr& cand) const {
0052   const double rho = (*rhoHandle_);
0053   const float isolTrkPt = cand->dr03TkSumPt();
0054 
0055   const float et = cand->et();
0056   const float cutValue =
0057       et > slopeStart_(cand) ? slopeTerm_(cand) * (et - slopeStart_(cand)) + constTerm_(cand) : constTerm_(cand);
0058 
0059   const float rhoCutValue = et >= rhoEtStart_(cand) ? rhoEA_(cand) * rho : 0.;
0060 
0061   return isolTrkPt < cutValue + rhoCutValue;
0062 }
0063 
0064 double GsfEleTrkPtIsoRhoCut::value(const reco::CandidatePtr& cand) const {
0065   reco::GsfElectronPtr ele(cand);
0066   return ele->dr03TkSumPt();
0067 }