File indexing completed on 2025-04-02 23:19:45
0001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0002 #include "DataFormats/L1TMuonPhase2/interface/TrackerMuon.h"
0003 #include "L1Trigger/Phase2L1GMT/interface/TPS.h"
0004
0005 using namespace Phase2L1GMT;
0006
0007 TPS::TPS(const edm::ParameterSet& iConfig)
0008 : verbose_(iConfig.getParameter<int>("verbose")),
0009 tt_track_converter_(new TrackConverter(iConfig.getParameter<edm::ParameterSet>("trackConverter"))),
0010 tps_(new TPSAlgorithm(iConfig.getParameter<edm::ParameterSet>("trackMatching"))),
0011 isolation_(new Isolation(iConfig.getParameter<edm::ParameterSet>("isolation"))) {}
0012
0013 std::vector<l1t::TrackerMuon> TPS::processEvent(const std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> >& tracks,
0014 const l1t::MuonStubRefVector& muonStubs) {
0015
0016 std::array<std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> >, 9> loctracks;
0017 for (unsigned i = 0; i < 9; ++i)
0018 loctracks[i] = associateTracksWithNonant(tracks, i);
0019
0020
0021 std::array<std::vector<ConvertedTTTrack>, 9> convertedTracks;
0022 for (unsigned i = 0; i < 9; ++i)
0023 convertedTracks[i] = tt_track_converter_->convertTracks(loctracks.at(i));
0024
0025
0026 std::array<l1t::MuonStubRefVector, 9> stubs;
0027 for (int i = 0; i < 9; ++i)
0028 stubs[i] = associateStubsWithNonant(muonStubs, i);
0029
0030
0031 std::array<std::vector<PreTrackMatchedMuon>, 9> mus;
0032 for (int i = 0; i < 9; ++i)
0033 mus[i] = tps_->processNonant(convertedTracks.at(i), stubs.at(i));
0034
0035
0036 std::vector<std::vector<PreTrackMatchedMuon> > muCleaneds;
0037 muCleaneds.push_back(tps_->cleanNeighbor(mus[0], mus[8], mus[1], true));
0038 muCleaneds.push_back(tps_->cleanNeighbor(mus[1], mus[0], mus[2], false));
0039 muCleaneds.push_back(tps_->cleanNeighbor(mus[2], mus[1], mus[3], true));
0040 muCleaneds.push_back(tps_->cleanNeighbor(mus[3], mus[2], mus[4], false));
0041 muCleaneds.push_back(tps_->cleanNeighbor(mus[4], mus[3], mus[5], true));
0042 muCleaneds.push_back(tps_->cleanNeighbor(mus[5], mus[4], mus[6], false));
0043 muCleaneds.push_back(tps_->cleanNeighbor(mus[6], mus[5], mus[7], true));
0044 muCleaneds.push_back(tps_->cleanNeighbor(mus[7], mus[6], mus[8], false));
0045
0046 muCleaneds.push_back(tps_->cleanNeighbor(mus[8], mus[7], mus[0], false));
0047
0048
0049 std::vector<PreTrackMatchedMuon> mergedCleaned;
0050 for (auto&& v : muCleaneds) {
0051 mergedCleaned.insert(mergedCleaned.end(), v.begin(), v.end());
0052 }
0053
0054 std::vector<l1t::TrackerMuon> trackMatchedMuonsNoIso = tps_->convert(mergedCleaned, 32);
0055
0056
0057 std::vector<l1t::TrackerMuon> sortedTrackMuonsNoIso = tps_->sort(trackMatchedMuonsNoIso, 12);
0058
0059 tps_->SetQualityBits(sortedTrackMuonsNoIso);
0060
0061
0062 std::vector<ConvertedTTTrack> mergedconvertedTracks;
0063 for (auto&& v : convertedTracks) {
0064 mergedconvertedTracks.insert(mergedconvertedTracks.end(), v.begin(), v.end());
0065 }
0066
0067 isolation_->isolation_allmu_alltrk(sortedTrackMuonsNoIso, mergedconvertedTracks);
0068
0069
0070
0071
0072 tps_->outputGT(sortedTrackMuonsNoIso);
0073
0074 return sortedTrackMuonsNoIso;
0075 }
0076
0077 std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > TPS::associateTracksWithNonant(
0078 const std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> >& tracks, uint processor) {
0079 std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> > out;
0080 for (const auto& track : tracks) {
0081 if (track->phiSector() == processor) {
0082 out.push_back(track);
0083 }
0084 }
0085 return out;
0086 }
0087
0088 l1t::SAMuonRefVector TPS::associateMuonsWithNonant(const l1t::SAMuonRefVector& muons, uint processor) {
0089 l1t::SAMuonRefVector out;
0090
0091 ap_int<BITSPHI> center = ap_int<BITSPHI>(processor * 910);
0092
0093 for (const auto& s : muons) {
0094 ap_int<BITSSTUBCOORD> deltaPhi = s->hwPhi() - center;
0095 ap_uint<BITSPHI - 1> absDeltaPhi =
0096 (deltaPhi < 0) ? ap_uint<BITSPHI - 1>(-deltaPhi) : ap_uint<BITSPHI - 1>(deltaPhi);
0097 if (absDeltaPhi < 683)
0098 out.push_back(s);
0099 }
0100 return out;
0101 }
0102
0103 l1t::MuonStubRefVector TPS::associateStubsWithNonant(const l1t::MuonStubRefVector& allStubs, uint processor) {
0104 l1t::MuonStubRefVector out;
0105
0106 ap_int<BITSSTUBCOORD> center = ap_int<BITSSTUBCOORD>((processor * 910) / 8);
0107
0108 for (const auto& s : allStubs) {
0109 ap_int<BITSSTUBCOORD> phi = 0;
0110 if (s->quality() & 0x1)
0111 phi = s->coord1();
0112 else
0113 phi = s->coord2();
0114
0115 ap_int<BITSSTUBCOORD> deltaPhi = phi - center;
0116 ap_uint<BITSSTUBCOORD - 1> absDeltaPhi =
0117 (deltaPhi < 0) ? ap_uint<BITSSTUBCOORD - 1>(-deltaPhi) : ap_uint<BITSSTUBCOORD - 1>(deltaPhi);
0118 if (absDeltaPhi < 284)
0119 out.push_back(s);
0120 }
0121 return out;
0122 }