RecoEgamma
PhotonIdentification
plugins
cuts
PhoGenericRhoPtScaledCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorWithEventContentBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/Photon.h
"
3
#include "
RecoEgamma/EgammaTools/interface/EffectiveAreas.h
"
4
#include "
RecoEgamma/EgammaTools/interface/ThreadSafeStringCut.h
"
5
#include "
RecoEgamma/EgammaTools/interface/EBEECutValues.h
"
6
7
class
PhoGenericRhoPtScaledCut
:
public
CutApplicatorWithEventContentBase
{
8
public
:
9
PhoGenericRhoPtScaledCut
(
const
edm::ParameterSet
&
c
);
10
11
result_type
operator()
(
const
reco::PhotonPtr
&)
const
final
;
12
13
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
14
void
getEventContent
(
const
edm::EventBase
&)
final
;
15
16
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
17
18
CandidateType
candidateType
()
const
final {
return
PHOTON
; }
19
20
private
:
21
ThreadSafeStringCut<StringObjectFunction<reco::Photon>
,
reco::Photon
>
varFunc_
;
22
bool
lessThan_
;
23
//cut value is constTerm + linearRhoTerm_*rho + linearPtTerm*pt + quadraticPtTerm*pt*pt
24
//note EBEECutValues & Effective areas are conceptually the same thing, both are eta
25
//binned constants, just Effective areas have arbitary rather than barrel/endcap binnng
26
EBEECutValues
constTerm_
;
27
EffectiveAreas
linearRhoTerm_
;
28
EBEECutValues
linearPtTerm_
;
29
EBEECutValues
quadraticPtTerm_
;
30
31
edm::Handle<double>
rhoHandle_
;
32
};
33
34
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
PhoGenericRhoPtScaledCut
,
"PhoGenericRhoPtScaledCut"
);
35
36
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
(
const
edm::ParameterSet
&
params
)
37
:
CutApplicatorWithEventContentBase
(
params
),
38
varFunc_(
params
.getParameter<
std
::
string
>(
"cutVariable"
)),
39
lessThan_(
params
.getParameter<
bool
>(
"lessThan"
)),
40
constTerm_(
params
,
"constTerm"
),
41
linearRhoTerm_(
params
.getParameter<
edm
::FileInPath>(
"effAreasConfigFile"
).
fullPath
()),
42
linearPtTerm_(
params
,
"linearPtTerm"
),
43
quadraticPtTerm_(
params
,
"quadPtTerm"
) {
44
edm::InputTag
rhoTag
=
params
.getParameter<
edm::InputTag
>(
"rho"
);
45
contentTags_
.emplace(
"rho"
,
rhoTag
);
46
}
47
48
void
PhoGenericRhoPtScaledCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
49
auto
rho
=
cc
.consumes<
double
>(
contentTags_
[
"rho"
]);
50
contentTokens_
.emplace(
"rho"
,
rho
);
51
}
52
53
void
PhoGenericRhoPtScaledCut::getEventContent
(
const
edm::EventBase
&
ev
) {
54
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
);
55
}
56
57
CutApplicatorBase::result_type
PhoGenericRhoPtScaledCut::operator()
(
const
reco::PhotonPtr
& pho)
const
{
58
const
double
rho
= (*rhoHandle_);
59
60
const
float
var
=
varFunc_
(*pho);
61
62
const
float
et
= pho->et();
63
const
float
absEta =
std::abs
(pho->superCluster()->eta());
64
const
float
cutValue =
constTerm_
(pho) +
linearRhoTerm_
.
getEffectiveArea
(absEta) *
rho
+
linearPtTerm_
(pho) *
et
+
65
quadraticPtTerm_
(pho) *
et
*
et
;
66
if
(
lessThan_
)
67
return
var
< cutValue;
68
else
69
return
var
>= cutValue;
70
}
71
72
double
PhoGenericRhoPtScaledCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
73
reco::PhotonPtr
pho(
cand
);
74
return
varFunc_
(*pho);
75
}
electrons_cff.bool
bool
Definition:
electrons_cff.py:372
CutApplicatorBase::PHOTON
Definition:
CutApplicatorBase.h:47
PhoGenericRhoPtScaledCut::constTerm_
EBEECutValues constTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:26
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
HLT_2018_cff.rhoTag
rhoTag
Definition:
HLT_2018_cff.py:13606
PhoGenericRhoPtScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
PhoGenericRhoPtScaledCut.cc:72
PhoGenericRhoPtScaledCut::linearRhoTerm_
EffectiveAreas linearRhoTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:27
watchdog.const
const
Definition:
watchdog.py:83
edm::Handle< double >
ThreadSafeStringCut.h
PhoGenericRhoPtScaledCut::candidateType
CandidateType candidateType() const final
Definition:
PhoGenericRhoPtScaledCut.cc:18
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
ThreadSafeStringCut
Definition:
ThreadSafeStringCut.h:17
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:8
CutApplicatorWithEventContentBase.h
PhoGenericRhoPtScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
PhoGenericRhoPtScaledCut.cc:48
PhoGenericRhoPtScaledCut
Definition:
PhoGenericRhoPtScaledCut.cc:7
PhoGenericRhoPtScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
PhoGenericRhoPtScaledCut.cc:31
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
DDAxes::rho
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
PhoGenericRhoPtScaledCut::quadraticPtTerm_
EBEECutValues quadraticPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:29
edm::ParameterSet
Definition:
ParameterSet.h:36
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
cand
Definition:
decayParser.h:34
EgHLTOffHistBins_cfi.et
et
Definition:
EgHLTOffHistBins_cfi.py:8
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
cc
edm::Ptr
Definition:
AssociationVector.h:31
reco::Photon
Definition:
Photon.h:21
std
Definition:
JetResolutionObject.h:76
PhoGenericRhoPtScaledCut::linearPtTerm_
EBEECutValues linearPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:28
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::varFunc_
ThreadSafeStringCut< StringObjectFunction< reco::Photon >, reco::Photon > varFunc_
Definition:
PhoGenericRhoPtScaledCut.cc:21
PhoGenericRhoPtScaledCut::operator()
result_type operator()(const reco::PhotonPtr &) const final
Definition:
PhoGenericRhoPtScaledCut.cc:57
EffectiveAreas::getEffectiveArea
const float getEffectiveArea(float eta) const
Definition:
EffectiveAreas.cc:44
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
PhoGenericRhoPtScaledCut(const edm::ParameterSet &c)
Definition:
PhoGenericRhoPtScaledCut.cc:36
EBEECutValuesT< double >
PhoGenericRhoPtScaledCut::lessThan_
bool lessThan_
Definition:
PhoGenericRhoPtScaledCut.cc:22
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
EffectiveAreas.h
PhoGenericRhoPtScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
PhoGenericRhoPtScaledCut.cc:53
edm::InputTag
Definition:
InputTag.h:15
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
Generated for CMSSW Reference Manual by
1.8.16