Main Page
Namespaces
Classes
Package Documentation
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/ElectronIdentification/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 {
19
return
PHOTON
;
20
}
21
22
private
:
23
ThreadSafeStringCut<StringObjectFunction<reco::Photon>
,
reco::Photon
>
varFunc_
;
24
bool
lessThan_
;
25
//cut value is constTerm + linearRhoTerm_*rho + linearPtTerm*pt + quadraticPtTerm*pt*pt
26
//note EBEECutValues & Effective areas are conceptually the same thing, both are eta
27
//binned constants, just Effective areas have arbitary rather than barrel/endcap binnng
28
EBEECutValues
constTerm_
;
29
EffectiveAreas
linearRhoTerm_
;
30
EBEECutValues
linearPtTerm_
;
31
EBEECutValues
quadraticPtTerm_
;
32
33
edm::Handle<double>
rhoHandle_
;
34
};
35
36
37
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
38
PhoGenericRhoPtScaledCut
,
39
"PhoGenericRhoPtScaledCut"
);
40
41
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
(
const
edm::ParameterSet
& params) :
42
CutApplicatorWithEventContentBase
(params),
43
varFunc_
(params.getParameter<
std
::
string
>(
"cutVariable"
)),
44
lessThan_
(params.getParameter<
bool
>(
"lessThan"
)),
45
constTerm_
(params,
"constTerm"
),
46
linearRhoTerm_
(params.getParameter<
edm
::FileInPath>(
"effAreasConfigFile"
).fullPath()),
47
linearPtTerm_
(params,
"linearPtTerm"
),
48
quadraticPtTerm_
(params,
"quadPtTerm"
)
49
{
50
edm::InputTag
rhoTag
= params.
getParameter
<
edm::InputTag
>(
"rho"
);
51
contentTags_
.emplace(
"rho"
,rhoTag);
52
}
53
54
void
PhoGenericRhoPtScaledCut::setConsumes
(
edm::ConsumesCollector
& cc) {
55
auto
rho
= cc.
consumes
<
double
>(
contentTags_
[
"rho"
]);
56
contentTokens_
.emplace(
"rho"
,
rho
);
57
}
58
59
void
PhoGenericRhoPtScaledCut::getEventContent
(
const
edm::EventBase
&
ev
) {
60
ev.
getByLabel
(
contentTags_
[
"rho"
],
rhoHandle_
);
61
}
62
63
CutApplicatorBase::result_type
64
PhoGenericRhoPtScaledCut::
65
operator()
(
const
reco::PhotonPtr
& pho)
const
{
66
const
double
rho
= (*rhoHandle_);
67
68
const
float
var
=
varFunc_
(*pho);
69
70
const
float
et
= pho->et();
71
const
float
absEta =
std::abs
(pho->superCluster()->eta());
72
const
float
cutValue =
constTerm_
(pho) +
linearRhoTerm_
.
getEffectiveArea
(absEta)*rho +
73
linearPtTerm_
(pho)*et +
quadraticPtTerm_
(pho)*et*
et
;
74
if
(
lessThan_
)
return
var < cutValue;
75
else
return
var>=cutValue;
76
}
77
78
double
PhoGenericRhoPtScaledCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
79
reco::PhotonPtr
pho(cand);
80
return
varFunc_
(*pho);
81
}
edm::ConsumesCollector::consumes
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition:
ConsumesCollector.h:49
DDAxes::rho
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
JetChargeProducer_cfi.var
var
Definition:
JetChargeProducer_cfi.py:4
CutApplicatorWithEventContentBase
Definition:
CutApplicatorWithEventContentBase.h:19
PhoGenericRhoPtScaledCut::lessThan_
bool lessThan_
Definition:
PhoGenericRhoPtScaledCut.cc:24
PhoGenericRhoPtScaledCut::getEventContent
void getEventContent(const edm::EventBase &) final
Definition:
PhoGenericRhoPtScaledCut.cc:59
regressionModifier_cfi.rhoTag
rhoTag
Definition:
regressionModifier_cfi.py:5
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
CutApplicatorWithEventContentBase::contentTokens_
std::unordered_map< std::string, edm::EDGetToken > contentTokens_
Definition:
CutApplicatorWithEventContentBase.h:44
PhoGenericRhoPtScaledCut::linearPtTerm_
EBEECutValues linearPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:30
EnergyCorrector.c
c
Definition:
EnergyCorrector.py:44
edm::Handle< double >
Photon.h
EffectiveAreas.h
std
Definition:
JetResolutionObject.h:80
ev
bool ev
Definition:
Hydjet2Hadronizer.cc:95
PhoGenericRhoPtScaledCut::PhoGenericRhoPtScaledCut
PhoGenericRhoPtScaledCut(const edm::ParameterSet &c)
Definition:
PhoGenericRhoPtScaledCut.cc:41
EffectiveAreas::getEffectiveArea
const float getEffectiveArea(float eta) const
Definition:
EffectiveAreas.cc:47
PhoGenericRhoPtScaledCut::rhoHandle_
edm::Handle< double > rhoHandle_
Definition:
PhoGenericRhoPtScaledCut.cc:33
reco::Photon
Definition:
Photon.h:22
EBEECutValuesT< double >
PhoGenericRhoPtScaledCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
PhoGenericRhoPtScaledCut.cc:78
CutApplicatorWithEventContentBase::contentTags_
std::unordered_map< std::string, edm::InputTag > contentTags_
Definition:
CutApplicatorWithEventContentBase.h:39
PhoGenericRhoPtScaledCut::candidateType
CandidateType candidateType() const final
Definition:
PhoGenericRhoPtScaledCut.cc:18
edmplugin::PluginFactory
Definition:
PluginFactory.h:33
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
PhoGenericRhoPtScaledCut::varFunc_
ThreadSafeStringCut< StringObjectFunction< reco::Photon >, reco::Photon > varFunc_
Definition:
PhoGenericRhoPtScaledCut.cc:23
CutApplicatorWithEventContentBase.h
edm::Ptr
Definition:
AssociationVector.h:30
PhoGenericRhoPtScaledCut
Definition:
PhoGenericRhoPtScaledCut.cc:7
PhoGenericRhoPtScaledCut::setConsumes
void setConsumes(edm::ConsumesCollector &) final
Definition:
PhoGenericRhoPtScaledCut.cc:54
edm::EventBase
Definition:
EventBase.h:46
electrons_cff.bool
bool
Definition:
electrons_cff.py:397
EBEECutValues.h
PhoGenericRhoPtScaledCut::quadraticPtTerm_
EBEECutValues quadraticPtTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:31
stringResolutionProvider_cfi.et
et
define resolution functions of each parameter
Definition:
stringResolutionProvider_cfi.py:13
PhoGenericRhoPtScaledCut::linearRhoTerm_
EffectiveAreas linearRhoTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:29
edm::EventBase::getByLabel
bool getByLabel(InputTag const &, Handle< T > &) const
Definition:
EventBase.h:92
edm
HLT enums.
Definition:
AlignableModifier.h:17
edm::InputTag
Definition:
InputTag.h:15
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
ThreadSafeStringCut
Definition:
ThreadSafeStringCut.h:17
edm::ParameterSet
Definition:
ParameterSet.h:36
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:121
cand
Definition:
decayParser.h:34
CutApplicatorBase::PHOTON
Definition:
CutApplicatorBase.h:48
EffectiveAreas
Definition:
EffectiveAreas.h:8
PhoGenericRhoPtScaledCut::constTerm_
EBEECutValues constTerm_
Definition:
PhoGenericRhoPtScaledCut.cc:28
ThreadSafeStringCut.h
edm::ConsumesCollector
Definition:
ConsumesCollector.h:39
PhoGenericRhoPtScaledCut::operator()
result_type operator()(const reco::PhotonPtr &) const final
Definition:
PhoGenericRhoPtScaledCut.cc:65
Generated for CMSSW Reference Manual by
1.8.11