src
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 "
CommonTools/Egamma/interface/EffectiveAreas.h
"
4
#include "
CommonTools/Utils/interface/StringObjectFunction.h
"
5
#include "
CommonTools/Utils/interface/ThreadSafeFunctor.h
"
6
#include "
RecoEgamma/EgammaTools/interface/EBEECutValues.h
"
7
8
class
PhoGenericRhoPtScaledCut
:
public
CutApplicatorWithEventContentBase
{
9
public
:
10
PhoGenericRhoPtScaledCut
(
const
edm::ParameterSet
&
c
);
11
12
result_type
operator()
(
const
reco::PhotonPtr
&)
const
final
;
13
14
void
setConsumes
(
edm::ConsumesCollector
&)
final
;
15
void
getEventContent
(
const
edm::EventBase
&)
final
;
16
17
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
18
19
CandidateType
candidateType
()
const
final {
return
PHOTON
; }
20
21
private
:
22
ThreadSafeFunctor<StringObjectFunction<reco::Photon>
>
varFunc_
;
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
27
EBEECutValues
constTerm_
;
28
EffectiveAreas
linearRhoTerm_
;
29
EBEECutValues
linearPtTerm_
;
30
EBEECutValues
quadraticPtTerm_
;
31
32
edm::Handle<double>
rhoHandle_
;
33
};
34
35
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
PhoGenericRhoPtScaledCut
,
"PhoGenericRhoPtScaledCut"
);
36
37
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
(
const
edm::ParameterSet
&
params
)
38
:
CutApplicatorWithEventContentBase
(
params
),
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
49
void
PhoGenericRhoPtScaledCut::setConsumes
(
edm::ConsumesCollector
&
cc
) {
50
auto
rho
=
cc
.consumes<
double
>(
contentTags_
[
"rho"
]);
51
contentTokens_
.emplace(
"rho"
,
rho
);
52
}
53
54
void
PhoGenericRhoPtScaledCut::getEventContent
(
const
edm::EventBase
&
ev
) {
55
ev
.getByLabel(
contentTags_
[
"rho"
],
rhoHandle_
);
56
}
57
58
CutApplicatorBase::result_type
PhoGenericRhoPtScaledCut::operator()
(
const
reco::PhotonPtr
& pho)
const
{
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
73
double
PhoGenericRhoPtScaledCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
74
reco::PhotonPtr
pho(
cand
);
75
return
varFunc_
(*pho);
76
}
DDAxes::rho
makeMEIFBenchmarkPlots.ev
ev
Definition:
makeMEIFBenchmarkPlots.py:55
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
EgHLTOffHistBins_cfi.et
et
Definition:
EgHLTOffHistBins_cfi.py:8
PhoGenericRhoPtScaledCut::lessThan_
bool lessThan_
Definition:
PhoGenericRhoPtScaledCut.cc:23
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:35
PhoGenericRhoPtScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
PhoGenericRhoPtScaledCut.cc:54
runTheMatrix.const
const
Definition:
runTheMatrix.py:341
trigObjTnPSource_cfi.var
var
Definition:
trigObjTnPSource_cfi.py:21
gpuPixelDoublets::cc
uint32_t cc[maxCellsPerHit]
Definition:
gpuFishbone.h:49
PhoGenericRhoPtScaledCut::linearPtTerm_
EBEECutValues linearPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:29
edm::Handle< double >
Photon.h
EffectiveAreas.h
std
Definition:
JetResolutionObject.h:76
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
PhoGenericRhoPtScaledCut(const edm::ParameterSet &c)
Definition:
PhoGenericRhoPtScaledCut.cc:37
PhoGenericRhoPtScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
PhoGenericRhoPtScaledCut.cc:32
EBEECutValuesT< double >
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
StringObjectFunction.h
PhoGenericRhoPtScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
PhoGenericRhoPtScaledCut.cc:73
ThreadSafeFunctor
Definition:
ThreadSafeFunctor.h:12
edmplugin::PluginFactory
Definition:
PluginFactory.h:35
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
PhoGenericRhoPtScaledCut::candidateType
CandidateType candidateType() const final
Definition:
PhoGenericRhoPtScaledCut.cc:19
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
HLT_2023v12_cff.rhoTag
rhoTag
Definition:
HLT_2023v12_cff.py:17942
CutApplicatorWithEventContentBase.h
edm::Ptr
Definition:
AssociationVector.h:31
PhoGenericRhoPtScaledCut
Definition:
PhoGenericRhoPtScaledCut.cc:8
PhoGenericRhoPtScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
PhoGenericRhoPtScaledCut.cc:49
edm::EventBase
Definition:
EventBase.h:47
submitPVValidationJobs.params
def params
Definition:
submitPVValidationJobs.py:482
EffectiveAreas::getEffectiveArea
const float getEffectiveArea(float eta) const
Definition:
EffectiveAreas.cc:76
PhoGenericRhoPtScaledCut::operator()
result_type operator()(const reco::PhotonPtr &) const final
Definition:
PhoGenericRhoPtScaledCut.cc:58
ThreadSafeFunctor.h
electrons_cff.bool
bool
Definition:
electrons_cff.py:315
EBEECutValues.h
PhoGenericRhoPtScaledCut::quadraticPtTerm_
EBEECutValues quadraticPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:30
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:40
PhoGenericRhoPtScaledCut::varFunc_
ThreadSafeFunctor< StringObjectFunction< reco::Photon > > varFunc_
Definition:
PhoGenericRhoPtScaledCut.cc:22
PhoGenericRhoPtScaledCut::linearRhoTerm_
EffectiveAreas linearRhoTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:28
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::InputTag
Definition:
InputTag.h:15
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edm::ParameterSet
Definition:
ParameterSet.h:48
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:123
cand
Definition:
decayParser.h:32
CutApplicatorBase::PHOTON
Definition:
CutApplicatorBase.h:47
EffectiveAreas
Definition:
EffectiveAreas.h:7
PhoGenericRhoPtScaledCut::constTerm_
EBEECutValues constTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:27
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
contentValuesFiles.fullPath
fullPath
Definition:
contentValuesFiles.py:64
Generated for CMSSW Reference Manual by
1.8.14