PhysicsTools
SelectorUtils
plugins
MinPtCutInEtaRanges.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h
"
2
3
class
MinPtCutInEtaRanges
:
public
CutApplicatorBase
{
4
public
:
5
MinPtCutInEtaRanges
(
const
edm::ParameterSet
&
c
) :
CutApplicatorBase
(
c
),
_absEta
(
c
.getParameter<
bool
>(
"useAbsEta"
)) {
6
const
std::vector<edm::ParameterSet>&
ranges
=
c
.getParameterSetVector(
"allowedEtaRanges"
);
7
for
(
const
auto
&
range
:
ranges
) {
8
const
double
minEta
=
range
.getParameter<
double
>(
"minEta"
);
9
const
double
maxEta
=
range
.getParameter<
double
>(
"maxEta"
);
10
const
double
minPt
=
range
.getParameter<
double
>(
"minPt"
);
11
_ranges
.emplace_back(
minEta
,
maxEta
);
12
_minPt
.push_back(
minPt
);
13
}
14
}
15
16
double
value
(
const
reco::CandidatePtr
&
cand
)
const
final
{
return
cand
->pt(); }
17
18
result_type
asCandidate
(
const
argument_type
&)
const
final
;
19
20
private
:
21
const
bool
_absEta
;
22
std::vector<std::pair<double, double> >
_ranges
;
23
std::vector<double>
_minPt
;
// indexed as above
24
};
25
26
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
MinPtCutInEtaRanges
,
"MinPtCutInEtaRanges"
);
27
28
CutApplicatorBase::result_type
MinPtCutInEtaRanges::asCandidate
(
const
argument_type
&
cand
)
const
{
29
const
double
the_eta = (
_absEta
?
std::abs
(
cand
->eta()) :
cand
->eta());
30
bool
result
=
false
;
31
for
(
unsigned
i
= 0;
i
<
_ranges
.size(); ++
i
) {
32
const
auto
&
range
=
_ranges
[
i
];
33
if
(the_eta >=
range
.first && the_eta < range.second && cand->
pt
() >
_minPt
[
i
]) {
34
result
=
true
;
35
break
;
36
}
37
}
38
return
result
;
39
}
diffTwoXMLs.ranges
string ranges
Definition:
diffTwoXMLs.py:79
FastTimerService_cff.range
range
Definition:
FastTimerService_cff.py:34
MinPtCutInEtaRanges::_ranges
std::vector< std::pair< double, double > > _ranges
Definition:
MinPtCutInEtaRanges.cc:22
electrons_cff.bool
bool
Definition:
electrons_cff.py:393
mps_fire.i
i
Definition:
mps_fire.py:428
MinPtCutInEtaRanges::MinPtCutInEtaRanges
MinPtCutInEtaRanges(const edm::ParameterSet &c)
Definition:
MinPtCutInEtaRanges.cc:5
DiDispStaMuonMonitor_cfi.pt
pt
Definition:
DiDispStaMuonMonitor_cfi.py:39
MinPtCutInEtaRanges::asCandidate
result_type asCandidate(const argument_type &) const final
Definition:
MinPtCutInEtaRanges.cc:28
MinPtCutInEtaRanges::value
double value(const reco::CandidatePtr &cand) const final
Definition:
MinPtCutInEtaRanges.cc:16
candidate_functions::CandidateCut::result_type
bool result_type
Definition:
CandidateCut.h:11
MinPtCutInEtaRanges
Definition:
MinPtCutInEtaRanges.cc:3
maxEta
double maxEta
Definition:
PFJetBenchmarkAnalyzer.cc:76
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:124
edm::ParameterSet
Definition:
ParameterSet.h:47
edmplugin::PluginFactory
Definition:
PluginFactory.h:34
cand
Definition:
decayParser.h:32
beam_dqm_sourceclient-live_cfg.minPt
minPt
Definition:
beam_dqm_sourceclient-live_cfg.py:323
HltBtagPostValidation_cff.c
c
Definition:
HltBtagPostValidation_cff.py:31
edm::Ptr< Candidate >
CutApplicatorBase
Definition:
CutApplicatorBase.h:45
candidate_functions::CandidateCut::argument_type
reco::CandidatePtr argument_type
Definition:
CandidateCut.h:10
MinPtCutInEtaRanges::_minPt
std::vector< double > _minPt
Definition:
MinPtCutInEtaRanges.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
EgHLTOffEleSelection_cfi.minEta
minEta
Definition:
EgHLTOffEleSelection_cfi.py:11
MinPtCutInEtaRanges::_absEta
const bool _absEta
Definition:
MinPtCutInEtaRanges.cc:21
Generated for CMSSW Reference Manual by
1.8.16