File indexing completed on 2024-04-06 12:20:50
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMEtaProcessor.h"
0026
0027
0028
0029
0030
0031 #include <iostream>
0032 #include <iomanip>
0033 #include <bitset>
0034
0035
0036
0037
0038
0039 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMTFConfig.h"
0040 #include "L1Trigger/L1TMuonBarrel/src/L1MuBMSectorProcessor.h"
0041 #include "L1Trigger/L1TMuonBarrel/interface/L1MuBMTrackFinder.h"
0042 #include "CondFormats/L1TObjects/interface/L1MuDTEtaPattern.h"
0043 #include "L1Trigger/L1TMuonBarrel/interface/L1MuBMTEtaPatternLut.h"
0044 #include "CondFormats/DataRecord/interface/L1MuDTEtaPatternLutRcd.h"
0045 #include "L1Trigger/L1TMuonBarrel/interface/L1MuBMTEtaPatternLut.h"
0046 #include "CondFormats/DataRecord/interface/L1MuDTQualPatternLutRcd.h"
0047 #include "CondFormats/L1TObjects/interface/L1MuDTTFMasks.h"
0048 #include "CondFormats/DataRecord/interface/L1MuDTTFMasksRcd.h"
0049
0050 #include "DataFormats/L1TMuon/interface/L1MuBMTrack.h"
0051 #include "DataFormats/L1TMuon/interface/L1MuBMTrackSegEta.h"
0052 #include "DataFormats/L1TMuon/interface/BMTF/L1MuBMSecProcId.h"
0053
0054 using namespace std;
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064 L1MuBMEtaProcessor::L1MuBMEtaProcessor(const L1MuBMTrackFinder& tf, int id, edm::ConsumesCollector&& iC)
0065 : m_tf(tf),
0066 m_epid(id),
0067 m_foundPattern(0),
0068 m_tseta(15),
0069 m_bmtfParamsToken(iC.esConsumes()),
0070 m_DTDigiToken(iC.consumes<L1MuDTChambThContainer>(L1MuBMTFConfig::getBMThetaDigiInputTag())) {
0071 m_tseta.reserve(15);
0072 }
0073
0074
0075
0076
0077
0078 L1MuBMEtaProcessor::~L1MuBMEtaProcessor() {}
0079
0080
0081
0082
0083
0084
0085
0086
0087 void L1MuBMEtaProcessor::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
0088 auto const& params = c.getData(m_bmtfParamsToken);
0089 if (L1MuBMTFConfig::getEtaTF()) {
0090 receiveData(bx, e, params);
0091 runEtaTrackFinder(params);
0092 }
0093
0094 receiveAddresses();
0095 runEtaMatchingUnit(params);
0096
0097 assign();
0098 }
0099
0100
0101
0102
0103 void L1MuBMEtaProcessor::reset() {
0104 vector<const L1MuBMTrackSegEta*>::iterator iter = m_tseta.begin();
0105 while (iter != m_tseta.end()) {
0106 if (*iter) {
0107 delete *iter;
0108 *iter = nullptr;
0109 }
0110 iter++;
0111 }
0112
0113 m_tseta.clear();
0114
0115 for (int i = 0; i < 12; i++) {
0116 m_eta[i] = 99;
0117 m_fine[i] = false;
0118 m_pattern[i] = 0;
0119 m_address[i] = 0;
0120 m_TrackCand[i] = nullptr;
0121 m_TracKCand[i] = nullptr;
0122 }
0123
0124 m_foundPattern.clear();
0125
0126 m_mask = true;
0127 }
0128
0129
0130
0131
0132 void L1MuBMEtaProcessor::print() const {
0133 bool empty1 = true;
0134 for (int i = 0; i < 15; i++) {
0135 empty1 &= (m_tseta[i] == nullptr || m_tseta[i]->empty());
0136 }
0137
0138 bool empty2 = true;
0139 for (int i = 0; i < 12; i++) {
0140 empty2 &= (m_address[i] == 0);
0141 }
0142
0143 if (!empty1 || !empty2) {
0144 cout << "Eta processor " << m_epid << " : " << endl;
0145
0146
0147 if (!empty1) {
0148 cout << "Local pattern : " << endl;
0149 for (int i = 0; i < 15; i++) {
0150 if ((i + 5) % 5 == 0)
0151 cout << "station " << m_tseta[i]->station() << " : ";
0152 bitset<7> pos(m_tseta[i]->position());
0153 bitset<7> qua(m_tseta[i]->quality());
0154 for (int j = 6; j >= 0; j--) {
0155 cout << pos[j] + qua[j];
0156 }
0157 cout << " ";
0158 if ((i + 1) % 5 == 0)
0159 cout << endl;
0160 }
0161 cout << "Found patterns :" << endl;
0162 vector<int>::const_iterator iter;
0163 for (iter = m_foundPattern.begin(); iter != m_foundPattern.end(); iter++) {
0164 const L1MuDTEtaPattern p = theEtaPatternLUT.getPattern(*iter);
0165 int qualitycode = p.quality();
0166 cout << "ID = " << setw(4) << p.id() << " "
0167 << "eta = " << setw(3) << p.eta() << " "
0168 << "quality = " << setw(2) << qualitycode << " (" << quality(qualitycode, 1) << " "
0169 << quality(qualitycode, 2) << " " << quality(qualitycode, 3) << ")";
0170 for (int i = 0; i < 12; i++) {
0171 if (m_pattern[i] == p.id())
0172 cout << " <--";
0173 }
0174 cout << endl;
0175 }
0176 }
0177
0178 cout << "Received addresses : " << endl;
0179 for (int i = 0; i < 12; i++)
0180 cout << setw(3) << m_address[i] << " ";
0181 cout << endl;
0182
0183 if (!empty1) {
0184 cout << "Matched patterns : " << endl;
0185 for (int i = 0; i < 12; i++) {
0186 if (m_fine[i]) {
0187 const L1MuDTEtaPattern p = theEtaPatternLUT.getPattern(m_pattern[i]);
0188 int fineeta = p.eta();
0189 int coarseeta = theQualPatternLUT.getCoarseEta(i / 2 + 1, m_address[i]);
0190 cout << "Index = " << setw(2) << i << ", "
0191 << "address = " << setw(2) << m_address[i] << " --> "
0192 << "pattern = " << setw(4) << m_pattern[i] << " "
0193 << "eta (coarse) = " << setw(3) << coarseeta << " "
0194 << "eta (fine) = " << setw(3) << fineeta << " "
0195 << "quality = " << setw(2) << p.quality() << endl;
0196 }
0197 }
0198 }
0199
0200 cout << "Eta values and fine bits : " << endl;
0201 for (int i = 0; i < 12; i++)
0202 cout << setw(3) << m_eta[i] << " ";
0203 cout << endl;
0204 for (int i = 0; i < 12; i++)
0205 cout << setw(3) << m_fine[i] << " ";
0206 cout << endl;
0207 }
0208 }
0209
0210
0211
0212
0213 void L1MuBMEtaProcessor::receiveData(int bx, const edm::Event& e, const L1TMuonBarrelParams& bmtfParams) {
0214 msks = bmtfParams.l1mudttfmasks;
0215 theEtaPatternLUT.m_lut = bmtfParams.lutparams_.eta_lut_;
0216 theQualPatternLUT.m_lut = bmtfParams.lutparams_.qp_lut_;
0217
0218 edm::Handle<L1MuDTChambThContainer> dttrig;
0219 e.getByToken(m_DTDigiToken, dttrig);
0220
0221
0222 int bx_offset = 0;
0223 bx = bx + bx_offset;
0224
0225
0226
0227
0228 int sector = m_epid;
0229 for (int stat = 1; stat <= 3; stat++) {
0230 for (int wheel = -2; wheel <= 2; wheel++) {
0231 L1MuDTChambThDigi const* tseta = dttrig->chThetaSegm(wheel, stat, sector, bx);
0232 bitset<7> pos;
0233 bitset<7> qual;
0234
0235 int lwheel = wheel + 1;
0236 if (wheel < 0)
0237 lwheel = wheel - 1;
0238
0239 bool masked = false;
0240 if (stat == 1)
0241 masked = msks.get_etsoc_chdis_st1(lwheel, sector);
0242 if (stat == 2)
0243 masked = msks.get_etsoc_chdis_st2(lwheel, sector);
0244 if (stat == 3)
0245 masked = msks.get_etsoc_chdis_st3(lwheel, sector);
0246
0247 if (!masked)
0248 m_mask = false;
0249
0250 if (tseta && !masked) {
0251 if (wheel == -2 || wheel == -1 ||
0252 (wheel == 0 && (sector == 0 || sector == 3 || sector == 4 || sector == 7 || sector == 8 || sector == 11))) {
0253 for (int i = 0; i < 7; i++) {
0254 if (tseta->position(i))
0255 pos.set(6 - i);
0256 if (tseta->quality(i))
0257 qual.set(6 - i);
0258 }
0259 } else {
0260 for (int i = 0; i < 7; i++) {
0261 if (tseta->position(i))
0262 pos.set(i);
0263 if (tseta->quality(i))
0264 qual.set(i);
0265 }
0266 }
0267 }
0268
0269 const L1MuBMTrackSegEta* tmpts =
0270 new L1MuBMTrackSegEta(wheel, sector, stat, pos.to_ulong(), qual.to_ulong(), bx - bx_offset);
0271 m_tseta.push_back(tmpts);
0272 }
0273 }
0274 }
0275
0276
0277
0278
0279 void L1MuBMEtaProcessor::receiveAddresses() {
0280
0281
0282
0283 int sector = m_epid;
0284
0285 int i = 0;
0286 for (int wheel = -3; wheel <= 3; wheel++) {
0287 if (wheel == 0)
0288 continue;
0289 L1MuBMSecProcId tmpspid(wheel, sector);
0290 for (int number = 0; number < 2; number++) {
0291 const L1MuBMTrack* cand = m_tf.sp(tmpspid)->track(number);
0292 const L1MuBMTrack* canD = m_tf.sp(tmpspid)->tracK(number);
0293 if (cand) {
0294 m_address[i] = cand->address().trackAddressCode();
0295 if (!cand->empty()) {
0296 m_TrackCand[i] = const_cast<L1MuBMTrack*>(cand);
0297 m_TracKCand[i] = const_cast<L1MuBMTrack*>(canD);
0298 }
0299 i++;
0300 }
0301 }
0302 }
0303 }
0304
0305
0306
0307
0308 void L1MuBMEtaProcessor::runEtaTrackFinder(const L1TMuonBarrelParams& bmtfParams) {
0309 theEtaPatternLUT.m_lut = bmtfParams.lutparams_.eta_lut_;
0310
0311
0312 bool empty = true;
0313 for (int i = 0; i < 15; i++) {
0314 empty &= m_tseta[i]->empty();
0315 }
0316 if (empty)
0317 return;
0318
0319
0320
0321
0322 L1MuBMTEtaPatternLut::ETFLut_iter it = theEtaPatternLUT.begin();
0323 while (it != theEtaPatternLUT.end()) {
0324 const L1MuDTEtaPattern pattern = (*it).second;
0325 int qualitycode = pattern.quality();
0326
0327 bool good = true;
0328
0329 for (int station = 0; station < 3; station++) {
0330 int q = quality(qualitycode, station + 1);
0331 int wheel = pattern.wheel(station + 1);
0332 int bin = pattern.position(station + 1);
0333 if (bin == 0)
0334 continue;
0335 bitset<7> pos = m_tseta[wheel + 2 + station * 5]->position();
0336 bitset<7> qual = m_tseta[wheel + 2 + station * 5]->quality();
0337
0338 good &= pos.test(bin - 1);
0339
0340 if (q == 2)
0341 good &= qual.test(bin - 1);
0342 }
0343
0344 if (good)
0345 m_foundPattern.push_back(pattern.id());
0346
0347 it++;
0348 }
0349 }
0350
0351
0352
0353
0354 void L1MuBMEtaProcessor::runEtaMatchingUnit(const L1TMuonBarrelParams& bmtfParams) {
0355 theQualPatternLUT.m_lut = bmtfParams.lutparams_.qp_lut_;
0356
0357
0358 for (int i = 0; i < 12; i++) {
0359 int adr = m_address[i];
0360 if (adr == 0)
0361 continue;
0362 int sp = i / 2 + 1;
0363
0364
0365 if (!m_mask)
0366 m_eta[i] = theQualPatternLUT.getCoarseEta(sp, adr);
0367
0368 if (m_foundPattern.empty())
0369 continue;
0370
0371
0372
0373 const vector<short>& qualifiedPatterns = theQualPatternLUT.getQualifiedPatterns(sp, adr);
0374 vector<short>::const_iterator iter;
0375 vector<int>::const_iterator f_iter;
0376 for (iter = qualifiedPatterns.begin(); iter != qualifiedPatterns.end(); iter++) {
0377 f_iter = find(m_foundPattern.begin(), m_foundPattern.end(), (*iter));
0378
0379 if (f_iter != m_foundPattern.end()) {
0380 const L1MuDTEtaPattern p = theEtaPatternLUT.getPattern(*f_iter);
0381
0382 m_fine[i] = true;
0383 m_eta[i] = p.eta();
0384 ;
0385 m_pattern[i] = (*f_iter);
0386 break;
0387 }
0388 }
0389 }
0390
0391
0392
0393
0394
0395 for (int i = 0; i < 6; i++) {
0396 int idx1 = 2 * i;
0397 int idx2 = 2 * i + 1;
0398 int adr1 = m_address[idx1];
0399 int adr2 = m_address[idx2];
0400 if (adr1 == 0 || adr2 == 0)
0401 continue;
0402 if (adr1 == adr2 && !m_mask) {
0403
0404 m_eta[idx1] = theQualPatternLUT.getCoarseEta(i + 1, adr1);
0405 m_pattern[idx1] = 0;
0406 m_fine[idx1] = false;
0407 m_eta[idx2] = theQualPatternLUT.getCoarseEta(i + 1, adr2);
0408
0409 m_pattern[idx2] = 0;
0410 m_fine[idx2] = false;
0411 }
0412 }
0413 }
0414
0415
0416
0417
0418 void L1MuBMEtaProcessor::assign() {
0419 for (int i = 0; i < 12; i++) {
0420 if (m_TrackCand[i]) {
0421 if (m_eta[i] != 99) {
0422 m_TrackCand[i]->setEta(m_eta[i]);
0423 m_TracKCand[i]->setEta(m_eta[i]);
0424 } else {
0425
0426 }
0427 if (m_fine[i]) {
0428 m_TrackCand[i]->setFineEtaBit();
0429 m_TracKCand[i]->setFineEtaBit();
0430
0431 const L1MuDTEtaPattern p = theEtaPatternLUT.getPattern(m_pattern[i]);
0432 vector<const L1MuBMTrackSegEta*> TSeta;
0433 const L1MuBMTrackSegEta* ts = nullptr;
0434 for (int stat = 0; stat < 3; stat++) {
0435 int wh = p.wheel(stat + 1);
0436 int pos = p.position(stat + 1);
0437 if (pos == 0)
0438 continue;
0439 ts = m_tseta[wh + 2 + stat * 5];
0440 TSeta.push_back(ts);
0441 }
0442 m_TrackCand[i]->setTSeta(TSeta);
0443 m_TracKCand[i]->setTSeta(TSeta);
0444 }
0445 }
0446 }
0447 }
0448
0449
0450
0451
0452 int L1MuBMEtaProcessor::quality(int id, int stat) {
0453
0454
0455
0456
0457
0458
0459
0460
0461
0462 const int qualcode[27][3] = {{0, 0, 0}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {2, 0, 0}, {0, 2, 0}, {0, 0, 2},
0463 {1, 1, 0}, {1, 0, 1}, {0, 1, 1}, {2, 1, 0}, {1, 2, 0}, {2, 0, 1}, {1, 0, 2},
0464 {0, 2, 1}, {0, 1, 2}, {2, 2, 0}, {2, 0, 2}, {0, 2, 2}, {1, 1, 1}, {2, 1, 1},
0465 {1, 2, 1}, {1, 1, 2}, {2, 2, 1}, {2, 1, 2}, {1, 2, 2}, {2, 2, 2}};
0466
0467 return qualcode[id][stat - 1];
0468 }