RecoEgamma
PhotonIdentification
plugins
cuts
PhoSCEtaMultiRangeCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/Photon.h
"
3
4
class
PhoSCEtaMultiRangeCut
:
public
CutApplicatorBase
{
5
public
:
6
PhoSCEtaMultiRangeCut
(
const
edm::ParameterSet
&
c
) :
CutApplicatorBase
(
c
),
_absEta
(
c
.getParameter<
bool
>(
"useAbsEta"
)) {
7
const
std::vector<edm::ParameterSet>&
ranges
=
c
.getParameterSetVector(
"allowedEtaRanges"
);
8
for
(
const
auto
&
range
:
ranges
) {
9
const
double
min
=
range
.getParameter<
double
>(
"minEta"
);
10
const
double
max
=
range
.getParameter<
double
>(
"maxEta"
);
11
_ranges
.emplace_back(
min
,
max
);
12
}
13
}
14
15
result_type
operator()
(
const
reco::PhotonPtr
&)
const
final
;
16
17
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
;
18
19
CandidateType
candidateType
()
const
final {
return
PHOTON
; }
20
21
private
:
22
const
bool
_absEta
;
23
std::vector<std::pair<double, double> >
_ranges
;
24
};
25
26
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
PhoSCEtaMultiRangeCut
,
"PhoSCEtaMultiRangeCut"
);
27
28
CutApplicatorBase::result_type
PhoSCEtaMultiRangeCut::operator()
(
const
reco::PhotonPtr
&
cand
)
const
{
29
const
reco::SuperClusterRef
& scref =
cand
->superCluster();
30
const
double
the_eta = (
_absEta
?
std::abs
(scref->eta()) : scref->eta());
31
bool
result
=
false
;
32
for
(
const
auto
&
range
:
_ranges
) {
33
if
(the_eta >=
range
.first && the_eta <
range
.second) {
34
result
=
true
;
35
break
;
36
}
37
}
38
return
result
;
39
}
40
41
double
PhoSCEtaMultiRangeCut::value
(
const
reco::CandidatePtr
&
cand
)
const
{
42
reco::PhotonPtr
pho(
cand
);
43
const
reco::SuperClusterRef
& scref = pho->superCluster();
44
return
(
_absEta
?
std::abs
(scref->eta()) : scref->eta());
45
}
diffTwoXMLs.ranges
string ranges
Definition:
diffTwoXMLs.py:79
FastTimerService_cff.range
range
Definition:
FastTimerService_cff.py:34
PhoSCEtaMultiRangeCut::operator()
result_type operator()(const reco::PhotonPtr &) const final
Definition:
PhoSCEtaMultiRangeCut.cc:28
electrons_cff.bool
bool
Definition:
electrons_cff.py:366
CutApplicatorBase::PHOTON
Definition:
CutApplicatorBase.h:47
min
T min(T a, T b)
Definition:
MathUtil.h:58
watchdog.const
const
Definition:
watchdog.py:83
edm::Ref< SuperClusterCollection >
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
PhoSCEtaMultiRangeCut::PhoSCEtaMultiRangeCut
PhoSCEtaMultiRangeCut(const edm::ParameterSet &c)
Definition:
PhoSCEtaMultiRangeCut.cc:6
Photon.h
PhoSCEtaMultiRangeCut::value
double value(const reco::CandidatePtr &cand) const final
Definition:
PhoSCEtaMultiRangeCut.cc:41
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
edm::ParameterSet
Definition:
ParameterSet.h:47
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
SiStripPI::max
Definition:
SiStripPayloadInspectorHelper.h:169
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
cand
Definition:
decayParser.h:32
edm::Ptr
Definition:
AssociationVector.h:31
PhoSCEtaMultiRangeCut::_absEta
const bool _absEta
Definition:
PhoSCEtaMultiRangeCut.cc:22
CutApplicatorBase
Definition:
CutApplicatorBase.h:45
PhoSCEtaMultiRangeCut
Definition:
PhoSCEtaMultiRangeCut.cc:4
PhoSCEtaMultiRangeCut::candidateType
CandidateType candidateType() const final
Definition:
PhoSCEtaMultiRangeCut.cc:19
PhoSCEtaMultiRangeCut::_ranges
std::vector< std::pair< double, double > > _ranges
Definition:
PhoSCEtaMultiRangeCut.cc:23
mps_fire.result
result
Definition:
mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition:
Abs.h:22
CutApplicatorBase.h
c
auto & c
Definition:
CAHitNtupletGeneratorKernelsImpl.h:56
Generated for CMSSW Reference Manual by
1.8.16