File indexing completed on 2023-03-17 11:13:08
0001 #ifndef PHASE2GMT_TEMPORARY_ASSOCIATOR
0002 #define PHASE2GMT_TEMPORARY_ASSOCIATOR
0003
0004 #include "ap_int.h"
0005 #include "MuonROI.h"
0006 #include "DataFormats/L1TMuonPhase2/interface/MuonStub.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008
0009 namespace Phase2L1GMT {
0010
0011 class ROITempAssociator {
0012 public:
0013 ROITempAssociator(const edm::ParameterSet& iConfig) {}
0014 ~ROITempAssociator() {}
0015
0016 std::vector<MuonROI> associate(int bx,
0017 const l1t::ObjectRefBxCollection<l1t::RegionalMuonCand>& muons,
0018 const l1t::MuonStubRefVector& stubs) {
0019 std::vector<MuonROI> out;
0020 l1t::MuonStubRefVector usedStubs;
0021
0022 if (muons.size() > 0) {
0023 for (unsigned int i = 0; i < muons.size(bx); ++i) {
0024 const l1t::RegionalMuonCandRef& mu = muons.at(bx, i);
0025 uint pt = mu->hwPt();
0026 uint charge = mu->hwSign();
0027
0028 float eta = mu->hwEta() * 0.010875;
0029
0030 int globalPhi = 0;
0031 if (mu->trackFinderType() == l1t::bmtf) {
0032 globalPhi = mu->processor() * 48 + mu->hwPhi() - 24;
0033 } else {
0034 globalPhi = mu->processor() * 96 + mu->hwPhi() + 24;
0035 }
0036
0037 float phi = globalPhi * 2 * M_PI / 576.0;
0038 if (phi > (M_PI))
0039 phi = phi - 2 * M_PI;
0040 else if (phi < (-M_PI))
0041 phi = phi + 2 * M_PI;
0042
0043 MuonROI roi(bx, charge, pt, 1);
0044 roi.setMuonRef(mu);
0045 l1t::MuonStubRefVector cleanedStubs = clean(stubs, usedStubs);
0046
0047 for (unsigned int layer = 0; layer <= 4; ++layer) {
0048 l1t::MuonStubRefVector selectedStubs = selectLayerBX(cleanedStubs, bx, layer);
0049 int bestStubINT = -1;
0050 float dPhi = 1000.0;
0051
0052 for (uint i = 0; i < selectedStubs.size(); ++i) {
0053 const l1t::MuonStubRef& stub = selectedStubs[i];
0054 float deltaPhi = (stub->quality() & 0x1) ? stub->offline_coord1() - phi : stub->offline_coord2() - phi;
0055 if (deltaPhi > M_PI)
0056 deltaPhi = deltaPhi - 2 * M_PI;
0057 else if (deltaPhi < -M_PI)
0058 deltaPhi = deltaPhi + 2 * M_PI;
0059 deltaPhi = fabs(deltaPhi);
0060 float deltaEta = (stub->etaQuality() == 0 || (stub->etaQuality() & 0x1))
0061 ? fabs(stub->offline_eta1() - eta)
0062 : fabs(stub->offline_eta2() - eta);
0063 if (deltaPhi < (M_PI / 6.0) && deltaEta < 0.3 && deltaPhi < dPhi) {
0064 dPhi = deltaPhi;
0065 bestStubINT = i;
0066 }
0067 }
0068 if (bestStubINT >= 0) {
0069 roi.addStub(selectedStubs[bestStubINT]);
0070 usedStubs.push_back(selectedStubs[bestStubINT]);
0071 }
0072 }
0073 if (out.size() < 16 && !roi.stubs().empty())
0074 out.push_back(roi);
0075 }
0076 }
0077
0078
0079 l1t::MuonStubRefVector cleanedStubs = clean(stubs, usedStubs);
0080
0081 while (!cleanedStubs.empty()) {
0082 MuonROI roi(bx, 0, 0, 0);
0083 roi.addStub(cleanedStubs[0]);
0084 usedStubs.push_back(cleanedStubs[0]);
0085 for (unsigned int layer = 0; layer <= 4; ++layer) {
0086 if (layer == cleanedStubs[0]->tfLayer())
0087 continue;
0088 l1t::MuonStubRefVector selectedStubs = selectLayerBX(cleanedStubs, bx, layer);
0089 if (!selectedStubs.empty()) {
0090 roi.addStub(selectedStubs[0]);
0091 usedStubs.push_back(selectedStubs[0]);
0092 }
0093 }
0094 if (!roi.stubs().empty())
0095 if (out.size() < 16)
0096 out.push_back(roi);
0097 cleanedStubs = clean(cleanedStubs, usedStubs);
0098 }
0099 return out;
0100 }
0101
0102 private:
0103 l1t::MuonStubRefVector selectLayerBX(const l1t::MuonStubRefVector& all, int bx, uint layer) {
0104 l1t::MuonStubRefVector out;
0105 for (const auto& stub : all) {
0106 if (stub->bxNum() == bx && stub->tfLayer() == layer)
0107 out.push_back(stub);
0108 }
0109 return out;
0110 }
0111
0112 l1t::MuonStubRefVector clean(const l1t::MuonStubRefVector& all, const l1t::MuonStubRefVector& used) {
0113 l1t::MuonStubRefVector out;
0114 for (const auto& stub : all) {
0115 bool keep = true;
0116 for (const auto& st : used) {
0117 if (st == stub) {
0118 keep = false;
0119 break;
0120 }
0121 }
0122 if (keep)
0123 out.push_back(stub);
0124 }
0125 return out;
0126 }
0127 };
0128 }
0129
0130 #endif