File indexing completed on 2023-10-25 09:55:29
0001 #include "L1Trigger/L1TMuonEndCap/interface/MicroGMTConverter.h"
0002
0003 MicroGMTConverter::MicroGMTConverter() {}
0004
0005 MicroGMTConverter::~MicroGMTConverter() {}
0006
0007 void MicroGMTConverter::convert(const int global_event_BX,
0008 const EMTFTrack& in_track,
0009 l1t::RegionalMuonCand& out_cand) const {
0010 l1t::tftype tftype = (in_track.Endcap() == 1) ? l1t::tftype::emtf_pos : l1t::tftype::emtf_neg;
0011 int sector = in_track.Sector() - 1;
0012
0013 out_cand.setHwPt(in_track.GMT_pt());
0014 out_cand.setHwPtUnconstrained(in_track.GMT_pt_dxy());
0015 out_cand.setHwDXY(in_track.GMT_dxy());
0016 out_cand.setHwPhi(in_track.GMT_phi());
0017 out_cand.setHwEta(in_track.GMT_eta());
0018 out_cand.setHwSign(in_track.GMT_charge());
0019 out_cand.setHwSignValid(in_track.GMT_charge_valid());
0020 out_cand.setHwQual(in_track.GMT_quality());
0021 out_cand.setHwHF(false);
0022 out_cand.setTFIdentifiers(sector, tftype);
0023
0024 const EMTFPtLUT& ptlut_data = in_track.PtLUT();
0025
0026
0027 int me1_ch_id =
0028 (ptlut_data.bt_vi[0] == 0 && ptlut_data.bt_vi[1] != 0) ? ptlut_data.bt_ci[1] + 16 : ptlut_data.bt_ci[0];
0029 int me2_ch_id = ptlut_data.bt_ci[2];
0030 int me3_ch_id = ptlut_data.bt_ci[3];
0031 int me4_ch_id = ptlut_data.bt_ci[4];
0032
0033 int me1_seg_id = (ptlut_data.bt_vi[0] == 0 && ptlut_data.bt_vi[1] != 0) ? ptlut_data.bt_si[1] : ptlut_data.bt_si[0];
0034 int me2_seg_id = ptlut_data.bt_si[2];
0035 int me3_seg_id = ptlut_data.bt_si[3];
0036 int me4_seg_id = ptlut_data.bt_si[4];
0037
0038 auto get_gmt_chamber_me1 = [](int ch) {
0039 int gmt_ch = 0;
0040 if (ch == 10)
0041 gmt_ch = 1;
0042 else if (ch == 11)
0043 gmt_ch = 2;
0044 else if (ch == 12)
0045 gmt_ch = 3;
0046 else if (ch == 3 + 16)
0047 gmt_ch = 4;
0048 else if (ch == 6 + 16)
0049 gmt_ch = 5;
0050 else if (ch == 9 + 16)
0051 gmt_ch = 6;
0052 return gmt_ch;
0053 };
0054
0055 auto get_gmt_chamber = [](int ch) {
0056 int gmt_ch = 0;
0057 if (ch == 10)
0058 gmt_ch = 1;
0059 else if (ch == 11)
0060 gmt_ch = 2;
0061 else if (ch == 3)
0062 gmt_ch = 3;
0063 else if (ch == 9)
0064 gmt_ch = 4;
0065 return gmt_ch;
0066 };
0067
0068
0069 int EMTF_kBX = (abs(global_event_BX) % 2048) - 25 + in_track.BX();
0070 if (EMTF_kBX < 0)
0071 EMTF_kBX += 2048;
0072
0073 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME1Seg, me1_seg_id);
0074 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME1Ch, get_gmt_chamber_me1(me1_ch_id));
0075 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME2Seg, me2_seg_id);
0076 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME2Ch, get_gmt_chamber(me2_ch_id));
0077 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME3Seg, me3_seg_id);
0078 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME3Ch, get_gmt_chamber(me3_ch_id));
0079 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME4Seg, me4_seg_id);
0080 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kME4Ch, get_gmt_chamber(me4_ch_id));
0081 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kTrkNum, in_track.Track_num());
0082 out_cand.setTrackSubAddress(l1t::RegionalMuonCand::kBX, EMTF_kBX);
0083 }
0084
0085 void MicroGMTConverter::convert_all(const edm::Event& iEvent,
0086 const EMTFTrackCollection& in_tracks,
0087 l1t::RegionalMuonCandBxCollection& out_cands) const {
0088 int gmtMinBX = -2;
0089 int gmtMaxBX = +2;
0090
0091 out_cands.clear();
0092 out_cands.setBXRange(gmtMinBX, gmtMaxBX);
0093
0094 for (const auto& in_track : in_tracks) {
0095 int bx = in_track.BX();
0096
0097 if (gmtMinBX <= bx && bx <= gmtMaxBX) {
0098 l1t::RegionalMuonCand out_cand;
0099
0100 convert(iEvent.bunchCrossing(), in_track, out_cand);
0101 out_cands.push_back(bx, out_cand);
0102 }
0103 }
0104
0105
0106 emtf::sort_uGMT_muons(out_cands);
0107 }
0108
0109 namespace emtf {
0110
0111 void sort_uGMT_muons(l1t::RegionalMuonCandBxCollection& cands) {
0112 int minBX = cands.getFirstBX();
0113 int maxBX = cands.getLastBX();
0114 int emtfMinProc = 0;
0115 int emtfMaxProc = 11;
0116
0117
0118 auto sortedCands = std::make_unique<l1t::RegionalMuonCandBxCollection>();
0119 sortedCands->clear();
0120 sortedCands->setBXRange(minBX, maxBX);
0121 for (int iBX = minBX; iBX <= maxBX; ++iBX) {
0122 for (int proc = emtfMinProc; proc <= emtfMaxProc; proc++) {
0123 for (l1t::RegionalMuonCandBxCollection::const_iterator cand = cands.begin(iBX); cand != cands.end(iBX);
0124 ++cand) {
0125 int cand_proc = cand->processor();
0126 if (cand->trackFinderType() == l1t::tftype::emtf_neg)
0127 cand_proc += 6;
0128 if (cand_proc != proc)
0129 continue;
0130 sortedCands->push_back(iBX, *cand);
0131 }
0132 }
0133 }
0134
0135
0136 std::swap(cands, *sortedCands);
0137 sortedCands.reset();
0138 }
0139
0140 }