File indexing completed on 2024-07-16 02:43:05
0001 #include "L1Trigger/Phase2L1GMT/interface/L1TPhase2GMTEndcapStubProcessor.h"
0002 #include "L1Trigger/L1TMuon/interface/MuonTriggerPrimitive.h"
0003 #include <cmath>
0004 #include <iostream>
0005 #include <string>
0006 #include <sstream>
0007
0008 L1TPhase2GMTEndcapStubProcessor::L1TPhase2GMTEndcapStubProcessor() : minBX_(-3), maxBX_(3) {}
0009
0010 L1TPhase2GMTEndcapStubProcessor::L1TPhase2GMTEndcapStubProcessor(const edm::ParameterSet& iConfig)
0011 : minBX_(iConfig.getParameter<int>("minBX")),
0012 maxBX_(iConfig.getParameter<int>("maxBX")),
0013 coord1LSB_(iConfig.getParameter<double>("coord1LSB")),
0014 coord2LSB_(iConfig.getParameter<double>("coord2LSB")),
0015 eta1LSB_(iConfig.getParameter<double>("eta1LSB")),
0016 eta2LSB_(iConfig.getParameter<double>("eta2LSB")),
0017 etaMatch_(iConfig.getParameter<double>("etaMatch")),
0018 phiMatch_(iConfig.getParameter<double>("phiMatch")),
0019 verbose_(iConfig.getParameter<unsigned int>("verbose")) {}
0020
0021 L1TPhase2GMTEndcapStubProcessor::~L1TPhase2GMTEndcapStubProcessor() {}
0022
0023 l1t::MuonStub L1TPhase2GMTEndcapStubProcessor::buildCSCOnlyStub(const CSCDetId& detid,
0024 const CSCCorrelatedLCTDigi& digi,
0025 const L1TMuon::GeometryTranslator* translator,
0026 unsigned int tag) {
0027 int endcap = detid.endcap();
0028 int station = detid.station();
0029 int chamber = detid.chamber();
0030 int ring = detid.ring();
0031
0032 L1TMuon::TriggerPrimitive primitive(detid, digi);
0033
0034 const GlobalPoint& gp = translator->getGlobalPoint(primitive);
0035
0036 int phi = int(gp.phi().value() / coord1LSB_);
0037 int eta1 = int(gp.eta() / eta1LSB_);
0038
0039 int wheel = 0;
0040 int sign = endcap == 1 ? -1 : 1;
0041
0042 if (ring == 3)
0043 wheel = sign * 3;
0044 else if (ring == 2)
0045 wheel = sign * 4;
0046 else if (ring == 1)
0047 wheel = sign * 5;
0048
0049 int sector = fabs(chamber);
0050
0051 int bx = digi.getBX() - 8;
0052 int quality = 1;
0053
0054 uint tfLayer = 0;
0055 if ((ring == 3 || ring == 2) && station == 1)
0056 tfLayer = 4;
0057 else if (ring == 1 && station == 1)
0058 tfLayer = 0;
0059 else if (station == 2)
0060 tfLayer = 2;
0061 else if (station == 3)
0062 tfLayer = 1;
0063 else if (station == 4)
0064 tfLayer = 3;
0065
0066 l1t::MuonStub stub(wheel, sector, station, tfLayer, phi, 0, tag, bx, quality, eta1, 0, 1, 0);
0067
0068 stub.setOfflineQuantities(gp.phi().value(), 0.0, gp.eta(), 0.0);
0069 return stub;
0070 }
0071
0072 l1t::MuonStub L1TPhase2GMTEndcapStubProcessor::buildRPCOnlyStub(const RPCDetId& detid,
0073 const RPCDigi& digi,
0074 const L1TMuon::GeometryTranslator* translator) {
0075 L1TMuon::TriggerPrimitive primitive(detid, digi);
0076 const GlobalPoint& gp = translator->getGlobalPoint(primitive);
0077
0078 int phi2 = int(gp.phi().value() / coord2LSB_);
0079 int eta2 = int(gp.eta() / eta2LSB_);
0080
0081 int wheel = -(6 - detid.ring()) * detid.region();
0082 int sector = (detid.sector() - 1) * 6 + detid.subsector();
0083 int station = detid.station();
0084 bool tag = detid.trIndex();
0085 int bx = digi.bx();
0086 int quality = 2;
0087
0088 int ring = detid.ring();
0089
0090 uint tfLayer = 0;
0091
0092 if ((ring == 3 || ring == 2) && station == 1)
0093 tfLayer = 4;
0094
0095
0096 else if (station == 2)
0097 tfLayer = 2;
0098 else if (station == 3)
0099 tfLayer = 1;
0100 else if (station == 4)
0101 tfLayer = 3;
0102
0103 l1t::MuonStub stub(wheel, sector, station, tfLayer, 0, phi2, tag, bx, quality, 0, eta2, 2, 0);
0104 stub.setOfflineQuantities(gp.phi().value(), gp.phi().value(), gp.eta(), gp.eta());
0105 return stub;
0106 }
0107
0108 l1t::MuonStubCollection L1TPhase2GMTEndcapStubProcessor::combineStubs(const l1t::MuonStubCollection& cscStubs,
0109 const l1t::MuonStubCollection& rpcStubs) {
0110 l1t::MuonStubCollection out;
0111 l1t::MuonStubCollection usedRPC;
0112 l1t::MuonStubCollection cleanedRPC;
0113
0114 l1t::MuonStubCollection cleanedCSC;
0115
0116
0117 l1t::MuonStubCollection allCSC = cscStubs;
0118
0119 while (!allCSC.empty()) {
0120 l1t::MuonStub stub = allCSC[0];
0121 l1t::MuonStubCollection freeCSC;
0122 for (uint i = 1; i < allCSC.size(); ++i) {
0123 if ((stub.etaRegion() == allCSC[i].etaRegion()) && (stub.depthRegion() == allCSC[i].depthRegion()) &&
0124 (fabs(deltaPhi(stub.offline_coord1(), allCSC[i].offline_coord1())) < 0.001)) {
0125 if (fabs(stub.offline_eta1() - allCSC[i].offline_eta1()) > 0.001) {
0126 stub.setEta(stub.eta1(), allCSC[i].eta1(), 3);
0127 stub.setOfflineQuantities(stub.offline_coord1(), 0.0, stub.offline_eta1(), allCSC[i].offline_eta1());
0128 }
0129 } else {
0130 freeCSC.push_back(allCSC[i]);
0131 }
0132 }
0133 cleanedCSC.push_back(stub);
0134 allCSC = freeCSC;
0135 }
0136
0137 for (const auto& csc : cleanedCSC) {
0138 int nRPC = 0;
0139 float phiF = 0.0;
0140 float etaF = 0.0;
0141 int phi = 0;
0142 int eta = 0;
0143 for (const auto& rpc : rpcStubs) {
0144 if (csc.depthRegion() != rpc.depthRegion())
0145 continue;
0146 if (fabs(deltaPhi(csc.offline_coord1(), rpc.offline_coord2())) < phiMatch_ &&
0147 fabs(csc.offline_eta1() - rpc.offline_eta2()) < etaMatch_ && csc.bxNum() == rpc.bxNum()) {
0148 phiF += rpc.offline_coord2();
0149 etaF += rpc.offline_eta2();
0150 phi += rpc.coord2();
0151 eta += rpc.eta2();
0152 nRPC++;
0153 usedRPC.push_back(rpc);
0154 }
0155 }
0156
0157 int finalRPCPhi = 0;
0158 int finalRPCEta = 0;
0159 double offline_finalRPCPhi = 0;
0160 double offline_finalRPCEta = 0;
0161 if (nRPC != 0) {
0162 finalRPCPhi = phi / nRPC;
0163 finalRPCEta = eta / nRPC;
0164 offline_finalRPCPhi = phiF / nRPC;
0165 offline_finalRPCEta = etaF / nRPC;
0166 l1t::MuonStub stub(csc.etaRegion(),
0167 csc.phiRegion(),
0168 csc.depthRegion(),
0169 csc.tfLayer(),
0170 csc.coord1(),
0171 finalRPCPhi,
0172 0,
0173 csc.bxNum(),
0174 3,
0175 csc.eta1(),
0176 finalRPCEta,
0177 3,
0178 0);
0179 stub.setOfflineQuantities(csc.offline_coord1(), offline_finalRPCPhi, csc.offline_eta1(), offline_finalRPCEta);
0180 out.push_back(stub);
0181 } else {
0182 out.push_back(csc);
0183 }
0184 }
0185
0186
0187
0188 for (const auto& rpc : rpcStubs) {
0189 bool keep = true;
0190 for (const auto& rpc2 : usedRPC) {
0191 if (rpc == rpc2) {
0192 keep = false;
0193 break;
0194 }
0195 }
0196 if (keep)
0197 cleanedRPC.push_back(rpc);
0198 }
0199
0200 while (!cleanedRPC.empty()) {
0201 l1t::MuonStubCollection freeRPC;
0202
0203 int nRPC = 1;
0204 float phiF = cleanedRPC[0].offline_coord2();
0205 float etaF = cleanedRPC[0].offline_eta2();
0206 int phi = cleanedRPC[0].coord2();
0207 int eta = cleanedRPC[0].eta2();
0208
0209 for (unsigned i = 1; i < cleanedRPC.size(); ++i) {
0210 if (fabs(deltaPhi(cleanedRPC[0].offline_coord2(), cleanedRPC[i].offline_coord2())) < phiMatch_ &&
0211 cleanedRPC[0].depthRegion() == cleanedRPC[i].depthRegion() &&
0212 fabs(cleanedRPC[0].offline_eta2() - cleanedRPC[i].offline_eta2()) < etaMatch_ &&
0213 cleanedRPC[0].bxNum() == cleanedRPC[i].bxNum()) {
0214 phiF += cleanedRPC[i].offline_coord2();
0215 etaF += cleanedRPC[i].offline_eta2();
0216 phi += cleanedRPC[i].coord2();
0217 eta += cleanedRPC[i].eta2();
0218 nRPC++;
0219 } else {
0220 freeRPC.push_back(cleanedRPC[i]);
0221 }
0222 }
0223 l1t::MuonStub stub(cleanedRPC[0].etaRegion(),
0224 cleanedRPC[0].phiRegion(),
0225 cleanedRPC[0].depthRegion(),
0226 cleanedRPC[0].tfLayer(),
0227 0,
0228 phi / nRPC,
0229 0,
0230 cleanedRPC[0].bxNum(),
0231 2,
0232 0,
0233 eta / nRPC,
0234 2,
0235 0);
0236 stub.setOfflineQuantities(phiF / nRPC, phiF / nRPC, etaF / nRPC, etaF / nRPC);
0237 out.push_back(stub);
0238 cleanedRPC = freeRPC;
0239 };
0240 return out;
0241 }
0242
0243 l1t::MuonStubCollection L1TPhase2GMTEndcapStubProcessor::makeStubs(
0244 const MuonDigiCollection<CSCDetId, CSCCorrelatedLCTDigi>& csc,
0245 const MuonDigiCollection<RPCDetId, RPCDigi>& cleaned,
0246 const L1TMuon::GeometryTranslator* t,
0247 const edm::EventSetup& iSetup) {
0248 l1t::MuonStubCollection cscStubs;
0249 auto chamber = csc.begin();
0250 auto chend = csc.end();
0251 for (; chamber != chend; ++chamber) {
0252 auto digi = (*chamber).second.first;
0253 auto dend = (*chamber).second.second;
0254 unsigned int tag = 0;
0255 for (; digi != dend; ++digi) {
0256 l1t::MuonStub stub = buildCSCOnlyStub((*chamber).first, *digi, t, tag);
0257 tag = tag + 1;
0258 if (stub.bxNum() >= minBX_ && stub.bxNum() <= maxBX_)
0259 cscStubs.push_back(stub);
0260 }
0261 }
0262
0263 l1t::MuonStubCollection rpcStubs;
0264
0265 auto rpcchamber = cleaned.begin();
0266 auto rpcchend = cleaned.end();
0267 for (; rpcchamber != rpcchend; ++rpcchamber) {
0268 if ((*rpcchamber).first.region() == 0)
0269 continue;
0270 auto digi = (*rpcchamber).second.first;
0271 auto dend = (*rpcchamber).second.second;
0272 for (; digi != dend; ++digi) {
0273 l1t::MuonStub stub = buildRPCOnlyStub((*rpcchamber).first, *digi, t);
0274 if (stub.bxNum() >= minBX_ && stub.bxNum() <= maxBX_)
0275 rpcStubs.push_back(stub);
0276 }
0277 }
0278
0279 l1t::MuonStubCollection combinedStubs = combineStubs(cscStubs, rpcStubs);
0280
0281 if (verbose_) {
0282 edm::LogInfo("EndcapStub") << "CSC Stubs";
0283 for (const auto& stub : cscStubs)
0284 edm::LogInfo("EndcapStub") << "CSC Stub bx=" << stub.bxNum() << " TF=" << stub.tfLayer()
0285 << " etaRegion=" << stub.etaRegion() << " phiRegion=" << stub.phiRegion()
0286 << " depthRegion=" << stub.depthRegion() << " coord1=" << stub.offline_coord1() << ","
0287 << stub.coord1() << " coord2=" << stub.offline_coord2() << "," << stub.coord2()
0288 << " eta1=" << stub.offline_eta1() << "," << stub.eta1()
0289 << " eta2=" << stub.offline_eta2() << "," << stub.eta2()
0290 << " quality=" << stub.quality() << " etaQuality=" << stub.etaQuality();
0291
0292 edm::LogInfo("EndcapStub") << "RPC Stubs";
0293 for (const auto& stub : rpcStubs)
0294 edm::LogInfo("EndcapStub") << "RPC Stub bx=" << stub.bxNum() << " TF=" << stub.tfLayer()
0295 << " etaRegion=" << stub.etaRegion() << " phiRegion=" << stub.phiRegion()
0296 << " depthRegion=" << stub.depthRegion() << " coord1=" << stub.offline_coord1() << ","
0297 << stub.coord1() << " coord2=" << stub.offline_coord2() << "," << stub.coord2()
0298 << " eta1=" << stub.offline_eta1() << "," << stub.eta1()
0299 << " eta2=" << stub.offline_eta2() << "," << stub.eta2()
0300 << " quality=" << stub.quality() << " etaQuality=" << stub.etaQuality();
0301
0302 for (const auto& stub : combinedStubs)
0303 edm::LogInfo("EndcapStub") << "Combined Stub bx=" << stub.bxNum() << " TF=" << stub.tfLayer()
0304 << " etaRegion=" << stub.etaRegion() << " phiRegion=" << stub.phiRegion()
0305 << " depthRegion=" << stub.depthRegion() << " coord1=" << stub.offline_coord1() << ","
0306 << stub.coord1() << " coord2=" << stub.offline_coord2() << "," << stub.coord2()
0307 << " eta1=" << stub.offline_eta1() << "," << stub.eta1()
0308 << " eta2=" << stub.offline_eta2() << "," << stub.eta2()
0309 << " quality=" << stub.quality() << " etaQuality=" << stub.etaQuality();
0310 }
0311
0312 return combinedStubs;
0313 }