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 "DataFormats/Common/interface/ValueMap.h"
0004 #include "RecoEgamma/EgammaTools/interface/EBEECutValues.h"
0005 
0006 class GsfEleValueMapIsoRhoCut : public CutApplicatorWithEventContentBase {
0007 public:
0008   GsfEleValueMapIsoRhoCut(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 getRhoCorr(const reco::GsfElectronPtr& cand) const;
0021 
0022   EBEECutValues slopeTerm_;
0023   EBEECutValues slopeStart_;
0024   EBEECutValues constTerm_;
0025   EBEECutValues rhoEtStart_;
0026   EBEECutValues rhoEA_;
0027 
0028   bool useRho_;
0029 
0030   edm::Handle<double> rhoHandle_;
0031   edm::Handle<edm::ValueMap<float> > valueHandle_;
0032 };
0033 
0034 DEFINE_EDM_PLUGIN(CutApplicatorFactory, GsfEleValueMapIsoRhoCut, "GsfEleValueMapIsoRhoCut");
0035 
0036 GsfEleValueMapIsoRhoCut::GsfEleValueMapIsoRhoCut(const edm::ParameterSet& params)
0037     : CutApplicatorWithEventContentBase(params),
0038       slopeTerm_(params, "slopeTerm"),
0039       slopeStart_(params, "slopeStart"),
0040       constTerm_(params, "constTerm"),
0041       rhoEtStart_(params, "rhoEtStart"),
0042       rhoEA_(params, "rhoEA") {
0043   auto rho = params.getParameter<edm::InputTag>("rho");
0044   if (!rho.label().empty()) {
0045     useRho_ = true;
0046     contentTags_.emplace("rho", rho);
0047   } else
0048     useRho_ = false;
0049 
0050   contentTags_.emplace("value", params.getParameter<edm::InputTag>("value"));
0051 }
0052 
0053 void GsfEleValueMapIsoRhoCut::setConsumes(edm::ConsumesCollector& cc) {
0054   if (useRho_)
0055     contentTokens_.emplace("rho", cc.consumes<double>(contentTags_["rho"]));
0056   contentTokens_.emplace("value", cc.consumes<edm::ValueMap<float> >(contentTags_["value"]));
0057 }
0058 
0059 void GsfEleValueMapIsoRhoCut::getEventContent(const edm::EventBase& ev) {
0060   if (useRho_)
0061     ev.getByLabel(contentTags_["rho"], rhoHandle_);
0062   ev.getByLabel(contentTags_["value"], valueHandle_);
0063 }
0064 
0065 CutApplicatorBase::result_type GsfEleValueMapIsoRhoCut::operator()(const reco::GsfElectronPtr& cand) const {
0066   const float val = (*valueHandle_)[cand];
0067 
0068   const float et = cand->et();
0069   const float cutValue =
0070       et > slopeStart_(cand) ? slopeTerm_(cand) * (et - slopeStart_(cand)) + constTerm_(cand) : constTerm_(cand);
0071   const float rhoCutValue = getRhoCorr(cand);
0072 
0073   return val < cutValue + rhoCutValue;
0074 }
0075 
0076 double GsfEleValueMapIsoRhoCut::value(const reco::CandidatePtr& cand) const {
0077   reco::GsfElectronPtr ele(cand);
0078   return (*valueHandle_)[cand];
0079 }
0080 
0081 float GsfEleValueMapIsoRhoCut::getRhoCorr(const reco::GsfElectronPtr& cand) const {
0082   if (!useRho_)
0083     return 0.;
0084   else {
0085     const double rho = (*rhoHandle_);
0086     return cand->et() >= rhoEtStart_(cand) ? rhoEA_(cand) * rho : 0.;
0087   }
0088 }