CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PhoGenericRhoPtScaledCut.cc
Go to the documentation of this file.
7 
9 public:
11 
12  result_type operator()(const reco::PhotonPtr&) const final;
13 
15  void getEventContent(const edm::EventBase&) final;
16 
17  double value(const reco::CandidatePtr& cand) const final;
18 
20 
21 private:
23  bool lessThan_;
24  //cut value is constTerm + linearRhoTerm_*rho + linearPtTerm*pt + quadraticPtTerm*pt*pt
25  //note EBEECutValues & Effective areas are conceptually the same thing, both are eta
26  //binned constants, just Effective areas have arbitary rather than barrel/endcap binnng
31 
33 };
34 
36 
39  varFunc_(params.getParameter<std::string>("cutVariable")),
40  lessThan_(params.getParameter<bool>("lessThan")),
41  constTerm_(params, "constTerm"),
42  linearRhoTerm_(params.getParameter<edm::FileInPath>("effAreasConfigFile").fullPath()),
43  linearPtTerm_(params, "linearPtTerm"),
44  quadraticPtTerm_(params, "quadPtTerm") {
46  contentTags_.emplace("rho", rhoTag);
47 }
48 
50  auto rho = cc.consumes<double>(contentTags_["rho"]);
51  contentTokens_.emplace("rho", rho);
52 }
53 
56 }
57 
59  const double rho = (*rhoHandle_);
60 
61  const float var = varFunc_(*pho);
62 
63  const float et = pho->et();
64  const float absEta = std::abs(pho->superCluster()->eta());
65  const float cutValue = constTerm_(pho) + linearRhoTerm_.getEffectiveArea(absEta) * rho + linearPtTerm_(pho) * et +
66  quadraticPtTerm_(pho) * et * et;
67  if (lessThan_)
68  return var < cutValue;
69  else
70  return var >= cutValue;
71 }
72 
74  reco::PhotonPtr pho(cand);
75  return varFunc_(*pho);
76 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EventSetup & c
std::unordered_map< std::string, edm::InputTag > contentTags_
void getEventContent(const edm::EventBase &) final
bool ev
PhoGenericRhoPtScaledCut(const edm::ParameterSet &c)
const float getEffectiveArea(float eta) const
edm::Handle< double > rhoHandle_
double value(const reco::CandidatePtr &cand) const final
list var
if using global norm cols_to_minmax = [&#39;t_delta&#39;, &#39;t_hmaxNearP&#39;,&#39;t_emaxNearP&#39;, &#39;t_hAnnular&#39;, &#39;t_eAnnular&#39;,&#39;t_pt&#39;,&#39;t_nVtx&#39;,&#39;t_ieta&#39;,&#39;t_eHcal10&#39;, &#39;t_eHcal30&#39;,&#39;t_rhoh&#39;,&#39;t_eHcal&#39;] df[cols_to_minmax] = df[cols_to_minmax].apply(lambda x: (x - x.min()) / (x.max() - x.min()) if (x.max() - x.min() &gt; 0) else 1.0/200.0)
CandidateType candidateType() const final
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setConsumes(edm::ConsumesCollector &) final
result_type operator()(const reco::PhotonPtr &) const final
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
ThreadSafeFunctor< StringObjectFunction< reco::Photon > > varFunc_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool getByLabel(InputTag const &, Handle< T > &) const
Definition: EventBase.h:92
#define DEFINE_EDM_PLUGIN(factory, type, name)