CMS 3D CMS Logo

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") {
45  edm::InputTag rhoTag = params.getParameter<edm::InputTag>("rho");
46  contentTags_.emplace("rho", rhoTag);
47 }
48 
50  auto rho = cc.consumes<double>(contentTags_["rho"]);
51  contentTokens_.emplace("rho", rho);
52 }
53 
55  ev.getByLabel(contentTags_["rho"], rhoHandle_);
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 }
electrons_cff.bool
bool
Definition: electrons_cff.py:393
CutApplicatorBase::PHOTON
Definition: CutApplicatorBase.h:47
PhoGenericRhoPtScaledCut::constTerm_
EBEECutValues constTerm_
Definition: PhoGenericRhoPtScaledCut.cc:27
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
EBEECutValues.h
contentValuesFiles.fullPath
fullPath
Definition: contentValuesFiles.py:64
edm
HLT enums.
Definition: AlignableModifier.h:19
PhoGenericRhoPtScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition: PhoGenericRhoPtScaledCut.cc:73
EffectiveAreas.h
PhoGenericRhoPtScaledCut::linearRhoTerm_
EffectiveAreas linearRhoTerm_
Definition: PhoGenericRhoPtScaledCut.cc:28
watchdog.const
const
Definition: watchdog.py:83
edm::Handle< double >
PhoGenericRhoPtScaledCut::candidateType
CandidateType candidateType() const final
Definition: PhoGenericRhoPtScaledCut.cc:19
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition: CutApplicatorWithEventContentBase.h:35
HLT_FULL_cff.rhoTag
rhoTag
Definition: HLT_FULL_cff.py:14986
candidate_functions::CandidateCut::result_type
bool result_type
Definition: CandidateCut.h:11
trigObjTnPSource_cfi.var
var
Definition: trigObjTnPSource_cfi.py:21
Photon.h
EffectiveAreas
Definition: EffectiveAreas.h:7
CutApplicatorWithEventContentBase.h
PhoGenericRhoPtScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition: PhoGenericRhoPtScaledCut.cc:49
PhoGenericRhoPtScaledCut
Definition: PhoGenericRhoPtScaledCut.cc:8
PhoGenericRhoPtScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition: PhoGenericRhoPtScaledCut.cc:32
PhoGenericRhoPtScaledCut::varFunc_
ThreadSafeFunctor< StringObjectFunction< reco::Photon > > varFunc_
Definition: PhoGenericRhoPtScaledCut.cc:22
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
ThreadSafeFunctor
Definition: ThreadSafeFunctor.h:12
DDAxes::rho
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PhoGenericRhoPtScaledCut::quadraticPtTerm_
EBEECutValues quadraticPtTerm_
Definition: PhoGenericRhoPtScaledCut.cc:30
edm::ParameterSet
Definition: ParameterSet.h:47
CutApplicatorBase::CandidateType
CandidateType
Definition: CutApplicatorBase.h:47
edmplugin::PluginFactory
Definition: PluginFactory.h:34
cand
Definition: decayParser.h:32
EgHLTOffHistBins_cfi.et
et
Definition: EgHLTOffHistBins_cfi.py:8
ThreadSafeFunctor.h
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
cc
edm::Ptr
Definition: AssociationVector.h:31
std
Definition: JetResolutionObject.h:76
PhoGenericRhoPtScaledCut::linearPtTerm_
EBEECutValues linearPtTerm_
Definition: PhoGenericRhoPtScaledCut.cc:29
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition: CutApplicatorWithEventContentBase.h:40
edm::EventBase
Definition: EventBase.h:46
PhoGenericRhoPtScaledCut::operator()
result_type operator()(const reco::PhotonPtr &) const final
Definition: PhoGenericRhoPtScaledCut.cc:58
EffectiveAreas::getEffectiveArea
const float getEffectiveArea(float eta) const
Definition: EffectiveAreas.cc:44
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
PhoGenericRhoPtScaledCut(const edm::ParameterSet &c)
Definition: PhoGenericRhoPtScaledCut.cc:37
EBEECutValuesT< double >
PhoGenericRhoPtScaledCut::lessThan_
bool lessThan_
Definition: PhoGenericRhoPtScaledCut.cc:23
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PhoGenericRhoPtScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition: PhoGenericRhoPtScaledCut.cc:54
StringObjectFunction.h
edm::InputTag
Definition: InputTag.h:15
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
CutApplicatorWithEventContentBase
Definition: CutApplicatorWithEventContentBase.h:19