File indexing completed on 2024-07-03 04:18:09
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "SeedingRegionByL1.h"
0012
0013 #include <algorithm>
0014 #include <set>
0015 #include <vector>
0016
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018
0019 ticl::SeedingRegionByL1::SeedingRegionByL1(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
0020 : SeedingRegionAlgoBase(conf, sumes),
0021 l1GTCandsToken_(sumes.consumes<l1t::P2GTCandidateCollection>(conf.getParameter<edm::InputTag>("l1GTCandColl"))),
0022 algoVerbosity_(conf.getParameter<int>("algo_verbosity")),
0023 minPt_(conf.getParameter<double>("minPt")),
0024 minAbsEta_(conf.getParameter<double>("minAbsEta")),
0025 maxAbsEta_(conf.getParameter<double>("maxAbsEta")),
0026 quality_(conf.getParameter<int>("quality")),
0027 qualityIsMask_(conf.getParameter<bool>("qualityIsMask")),
0028 applyQuality_(conf.getParameter<bool>("applyQuality")) {}
0029
0030 void ticl::SeedingRegionByL1::makeRegions(const edm::Event &ev,
0031 const edm::EventSetup &es,
0032 std::vector<TICLSeedingRegion> &result) {
0033 auto l1GTCands = ev.getHandle(l1GTCandsToken_);
0034 edm::ProductID l1gtcandsId = l1GTCands.id();
0035
0036 for (size_t indx = 0; indx < (*l1GTCands).size(); indx++) {
0037 const auto &l1GTCand = (*l1GTCands)[indx];
0038 double offlinePt = l1GTCand.pt();
0039 bool passQuality(false);
0040
0041 if (applyQuality_) {
0042 if (qualityIsMask_) {
0043 passQuality = (l1GTCand.hwQualityFlags() & quality_);
0044 } else {
0045 passQuality = (l1GTCand.hwQualityFlags() == quality_);
0046 }
0047 } else {
0048 passQuality = true;
0049 }
0050
0051 if ((offlinePt < minPt_) || (std::abs(l1GTCand.eta()) < minAbsEta_) || (std::abs(l1GTCand.eta()) > maxAbsEta_) ||
0052 !passQuality) {
0053 continue;
0054 }
0055
0056 int iSide = int(l1GTCand.eta() > 0);
0057 result.emplace_back(GlobalPoint(l1GTCand.p4().X(), l1GTCand.p4().Y(), l1GTCand.p4().Z()),
0058 GlobalVector(l1GTCand.px(), l1GTCand.py(), l1GTCand.pz()),
0059 iSide,
0060 indx,
0061 l1gtcandsId);
0062 }
0063
0064 std::sort(result.begin(), result.end(), [](const TICLSeedingRegion &a, const TICLSeedingRegion &b) {
0065 return a.directionAtOrigin.perp2() > b.directionAtOrigin.perp2();
0066 });
0067 }
0068
0069 void ticl::SeedingRegionByL1::fillPSetDescription(edm::ParameterSetDescription &desc) {
0070 desc.add<edm::InputTag>("l1GTCandColl", edm::InputTag("L1TkPhotonsHGC", "EG"));
0071 desc.add<double>("minPt", 10);
0072 desc.add<double>("minAbsEta", 1.479);
0073 desc.add<double>("maxAbsEta", 4.0);
0074 desc.add<int>("quality", 5);
0075 desc.add<bool>("qualityIsMask", false);
0076 desc.add<bool>("applyQuality", false);
0077 SeedingRegionAlgoBase::fillPSetDescription(desc);
0078 }