CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfEleEmHadD1IsoRhoCut.cc
Go to the documentation of this file.
6 
8 
10 public:
12 
13  result_type operator()(const reco::GsfElectronPtr&) const override final;
14 
15  void setConsumes(edm::ConsumesCollector&) override final;
16  void getEventContent(const edm::EventBase&) override final;
17 
18  CandidateType candidateType() const override final {
19  return ELECTRON;
20  }
21 
22 private:
23  float rhoConstant_;
27 
28 
30 
31 };
32 
35  "GsfEleEmHadD1IsoRhoCut");
36 
39  rhoConstant_(params.getParameter<double>("rhoConstant")),
40  slopeTerm_(params,"slopeTerm"),
41  slopeStart_(params,"slopeStart"),
42  constTerm_(params,"constTerm"){
43  edm::InputTag rhoTag = params.getParameter<edm::InputTag>("rho");
44  contentTags_.emplace("rho",rhoTag);
45 
46 }
47 
49  auto rho = cc.consumes<double>(contentTags_["rho"]);
50  contentTokens_.emplace("rho",rho);
51 }
52 
55 }
56 
57 CutApplicatorBase::result_type
59 operator()(const reco::GsfElectronPtr& cand) const{
60  const double rho = (*rhoHandle_);
61 
62  const float isolEmHadDepth1 = cand->dr03EcalRecHitSumEt() + cand->dr03HcalDepth1TowerSumEt();
63 
64  const float et = cand->et();
65  const float cutValue = et > slopeStart_(cand) ? slopeTerm_(cand)*(et-slopeStart_(cand)) + constTerm_(cand) : constTerm_(cand);
66  return isolEmHadDepth1 < cutValue + rhoConstant_*rho;
67 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
CandidateType candidateType() const overridefinal
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition: DDAxes.h:10
bool ev
void getEventContent(const edm::EventBase &) overridefinal
std::unordered_map< std::string, edm::InputTag > contentTags_
edm::Handle< double > rhoHandle_
GsfEleEmHadD1IsoRhoCut(const edm::ParameterSet &c)
string const
Definition: compareJSON.py:14
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:87
#define DEFINE_EDM_PLUGIN(factory, type, name)
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
void setConsumes(edm::ConsumesCollector &) overridefinal