CMS 3D CMS Logo

GsfEleSCEtaMultiRangeCut.cc
Go to the documentation of this file.
3 
5 public:
7  : CutApplicatorBase(c), _absEta(c.getParameter<bool>("useAbsEta")) {
8  const std::vector<edm::ParameterSet>& ranges = c.getParameterSetVector("allowedEtaRanges");
9  for (const auto& range : ranges) {
10  const double min = range.getParameter<double>("minEta");
11  const double max = range.getParameter<double>("maxEta");
12  _ranges.emplace_back(min, max);
13  }
14  }
15 
16  result_type operator()(const reco::GsfElectronPtr&) const final;
17 
18  double value(const reco::CandidatePtr& cand) const final;
19 
21 
22 private:
23  const bool _absEta;
24  std::vector<std::pair<double, double> > _ranges;
25 };
26 
28 
30  const reco::SuperClusterRef& scref = cand->superCluster();
31  const double the_eta = (_absEta ? std::abs(scref->eta()) : scref->eta());
32  bool result = false;
33  for (const auto& range : _ranges) {
34  if (the_eta >= range.first && the_eta < range.second) {
35  result = true;
36  break;
37  }
38  }
39  return result;
40 }
41 
44  const reco::SuperClusterRef& scref = ele->superCluster();
45  return (_absEta ? std::abs(scref->eta()) : scref->eta());
46 }
diffTwoXMLs.ranges
string ranges
Definition: diffTwoXMLs.py:79
FastTimerService_cff.range
range
Definition: FastTimerService_cff.py:34
electrons_cff.bool
bool
Definition: electrons_cff.py:393
min
T min(T a, T b)
Definition: MathUtil.h:58
GsfEleSCEtaMultiRangeCut::_absEta
const bool _absEta
Definition: GsfEleSCEtaMultiRangeCut.cc:23
GsfEleSCEtaMultiRangeCut::value
double value(const reco::CandidatePtr &cand) const final
Definition: GsfEleSCEtaMultiRangeCut.cc:42
GsfEleSCEtaMultiRangeCut::candidateType
CandidateType candidateType() const final
Definition: GsfEleSCEtaMultiRangeCut.cc:20
watchdog.const
const
Definition: watchdog.py:83
edm::Ref< SuperClusterCollection >
candidate_functions::CandidateCut::result_type
bool result_type
Definition: CandidateCut.h:11
GsfElectron.h
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
GsfEleSCEtaMultiRangeCut
Definition: GsfEleSCEtaMultiRangeCut.cc:4
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
GsfEleSCEtaMultiRangeCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const final
Definition: GsfEleSCEtaMultiRangeCut.cc:29
cand
Definition: decayParser.h:32
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
GsfEleSCEtaMultiRangeCut::GsfEleSCEtaMultiRangeCut
GsfEleSCEtaMultiRangeCut(const edm::ParameterSet &c)
Definition: GsfEleSCEtaMultiRangeCut.cc:6
edm::Ptr< reco::GsfElectron >
CutApplicatorBase
Definition: CutApplicatorBase.h:45
GsfEleSCEtaMultiRangeCut::_ranges
std::vector< std::pair< double, double > > _ranges
Definition: GsfEleSCEtaMultiRangeCut.cc:24
reco::GsfElectron::superCluster
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:163
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CutApplicatorBase.h
CutApplicatorBase::ELECTRON
Definition: CutApplicatorBase.h:47