Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:18:01

0001 #include "PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h"
0002 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0003 
0004 class PhoSCEtaMultiRangeCut : public CutApplicatorBase {
0005 public:
0006   PhoSCEtaMultiRangeCut(const edm::ParameterSet& c) : CutApplicatorBase(c), _absEta(c.getParameter<bool>("useAbsEta")) {
0007     const std::vector<edm::ParameterSet>& ranges = c.getParameterSetVector("allowedEtaRanges");
0008     for (const auto& range : ranges) {
0009       const double min = range.getParameter<double>("minEta");
0010       const double max = range.getParameter<double>("maxEta");
0011       _ranges.emplace_back(min, max);
0012     }
0013   }
0014 
0015   result_type operator()(const reco::PhotonPtr&) const final;
0016 
0017   double value(const reco::CandidatePtr& cand) const final;
0018 
0019   CandidateType candidateType() const final { return PHOTON; }
0020 
0021 private:
0022   const bool _absEta;
0023   std::vector<std::pair<double, double> > _ranges;
0024 };
0025 
0026 DEFINE_EDM_PLUGIN(CutApplicatorFactory, PhoSCEtaMultiRangeCut, "PhoSCEtaMultiRangeCut");
0027 
0028 CutApplicatorBase::result_type PhoSCEtaMultiRangeCut::operator()(const reco::PhotonPtr& cand) const {
0029   const reco::SuperClusterRef& scref = cand->superCluster();
0030   const double the_eta = (_absEta ? std::abs(scref->eta()) : scref->eta());
0031   bool result = false;
0032   for (const auto& range : _ranges) {
0033     if (the_eta >= range.first && the_eta < range.second) {
0034       result = true;
0035       break;
0036     }
0037   }
0038   return result;
0039 }
0040 
0041 double PhoSCEtaMultiRangeCut::value(const reco::CandidatePtr& cand) const {
0042   reco::PhotonPtr pho(cand);
0043   const reco::SuperClusterRef& scref = pho->superCluster();
0044   return (_absEta ? std::abs(scref->eta()) : scref->eta());
0045 }