Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:42

0001 //-------------------------------------------------
0002 //
0003 //   Class: L1MuDTEtaProcessor
0004 //
0005 //   Description: Eta Processor
0006 //
0007 //                An Eta Processor consists of :
0008 //                a receiver unit,
0009 //                one Eta Track Finder (ETF) and
0010 //                one Eta Matching Unit (EMU)
0011 //
0012 //
0013 //   Author :
0014 //   N. Neumeister            CERN EP
0015 //   J. Troconiz              UAM Madrid
0016 //
0017 //--------------------------------------------------
0018 
0019 //-----------------------
0020 // This Class's Header --
0021 //-----------------------
0022 
0023 #include "L1Trigger/DTTrackFinder/src/L1MuDTEtaProcessor.h"
0024 
0025 //---------------
0026 // C++ Headers --
0027 //---------------
0028 
0029 #include <iostream>
0030 #include <iomanip>
0031 #include <bitset>
0032 
0033 //-------------------------------
0034 // Collaborating Class Headers --
0035 //-------------------------------
0036 
0037 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTFConfig.h"
0038 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTrackSegEta.h"
0039 #include "L1Trigger/DTTrackFinder/interface/L1MuDTSecProcId.h"
0040 #include "L1Trigger/DTTrackFinder/src/L1MuDTSectorProcessor.h"
0041 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTrackFinder.h"
0042 #include "L1Trigger/DTTrackFinder/interface/L1MuDTTrack.h"
0043 #include "CondFormats/L1TObjects/interface/L1MuDTEtaPattern.h"
0044 #include "CondFormats/L1TObjects/interface/L1MuDTEtaPatternLut.h"
0045 #include "CondFormats/DataRecord/interface/L1MuDTEtaPatternLutRcd.h"
0046 #include "CondFormats/L1TObjects/interface/L1MuDTQualPatternLut.h"
0047 #include "CondFormats/DataRecord/interface/L1MuDTQualPatternLutRcd.h"
0048 #include "CondFormats/L1TObjects/interface/L1MuDTTFMasks.h"
0049 #include "CondFormats/DataRecord/interface/L1MuDTTFMasksRcd.h"
0050 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThDigi.h"
0051 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
0052 
0053 using namespace std;
0054 
0055 // --------------------------------
0056 //       class L1MuDTEtaProcessor
0057 //---------------------------------
0058 
0059 //----------------
0060 // Constructors --
0061 //----------------
0062 
0063 L1MuDTEtaProcessor::L1MuDTEtaProcessor(const L1MuDTTrackFinder& tf, int id, edm::ConsumesCollector iC)
0064     : m_tf(tf),
0065       m_epid(id),
0066       m_foundPattern(0),
0067       m_tseta(15),
0068       m_DTDigiToken(iC.consumes<L1MuDTChambThContainer>(m_tf.config()->getDTDigiInputTag())),
0069       theEtaToken(iC.esConsumes()),
0070       theQualToken(iC.esConsumes()),
0071       theMsksToken(iC.esConsumes()) {
0072   m_tseta.reserve(15);
0073 }
0074 
0075 //--------------
0076 // Destructor --
0077 //--------------
0078 
0079 L1MuDTEtaProcessor::~L1MuDTEtaProcessor() {}
0080 
0081 //--------------
0082 // Operations --
0083 //--------------
0084 
0085 //
0086 // run Eta Processor
0087 //
0088 void L1MuDTEtaProcessor::run(int bx, const edm::Event& e, const edm::EventSetup& c) {
0089   if (m_tf.config()->getEtaTF()) {
0090     receiveData(bx, e, c);
0091     runEtaTrackFinder(c);
0092   }
0093 
0094   receiveAddresses();
0095   runEtaMatchingUnit(c);
0096 
0097   assign();
0098 }
0099 
0100 //
0101 // reset Eta Processor
0102 //
0103 void L1MuDTEtaProcessor::reset() {
0104   vector<const L1MuDTTrackSegEta*>::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 // print track candidates found in Eta Processor
0131 //
0132 void L1MuDTEtaProcessor::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     // print local pattern
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 // receive data ( 15*3 DTBX eta trigger primitives )
0212 //
0213 void L1MuDTEtaProcessor::receiveData(int bx, const edm::Event& e, const edm::EventSetup& c) {
0214   msks = c.getHandle(theMsksToken);
0215 
0216   edm::Handle<L1MuDTChambThContainer> dttrig;
0217   e.getByToken(m_DTDigiToken, dttrig);
0218 
0219   // const int bx_offset = dttrig->correctBX();
0220   int bx_offset = 0;
0221   bx = bx + bx_offset;
0222 
0223   //
0224   // get 5*3 eta track segments
0225   //
0226   int sector = m_epid;
0227   for (int stat = 1; stat <= 3; stat++) {
0228     for (int wheel = -2; wheel <= 2; wheel++) {
0229       L1MuDTChambThDigi const* tseta = dttrig->chThetaSegm(wheel, stat, sector, bx);
0230       bitset<7> pos;
0231       bitset<7> qual;
0232 
0233       int lwheel = wheel + 1;
0234       if (wheel < 0)
0235         lwheel = wheel - 1;
0236 
0237       bool masked = false;
0238       if (stat == 1)
0239         masked = msks->get_etsoc_chdis_st1(lwheel, sector);
0240       if (stat == 2)
0241         masked = msks->get_etsoc_chdis_st2(lwheel, sector);
0242       if (stat == 3)
0243         masked = msks->get_etsoc_chdis_st3(lwheel, sector);
0244 
0245       if (!masked)
0246         m_mask = false;
0247 
0248       if (tseta && !masked) {
0249         if (wheel == -2 || wheel == -1 ||
0250             (wheel == 0 && (sector == 0 || sector == 3 || sector == 4 || sector == 7 || sector == 8 || sector == 11))) {
0251           for (int i = 0; i < 7; i++) {
0252             if (tseta->position(i))
0253               pos.set(6 - i);
0254             if (tseta->quality(i))
0255               qual.set(6 - i);
0256           }
0257         } else {
0258           for (int i = 0; i < 7; i++) {
0259             if (tseta->position(i))
0260               pos.set(i);
0261             if (tseta->quality(i))
0262               qual.set(i);
0263           }
0264         }
0265       }
0266 
0267       const L1MuDTTrackSegEta* tmpts =
0268           new L1MuDTTrackSegEta(wheel, sector, stat, pos.to_ulong(), qual.to_ulong(), bx - bx_offset);
0269       m_tseta.push_back(tmpts);
0270     }
0271   }
0272 }
0273 
0274 //
0275 // receive track addresses from 6 Sector Processors
0276 //
0277 void L1MuDTEtaProcessor::receiveAddresses() {
0278   // get track address code of all track segments
0279   // 6*2 times 5 bits; valid range [1-22]
0280 
0281   int sector = m_epid;
0282 
0283   int i = 0;
0284   for (int wheel = -3; wheel <= 3; wheel++) {
0285     if (wheel == 0)
0286       continue;
0287     L1MuDTSecProcId tmpspid(wheel, sector);
0288     for (int number = 0; number < 2; number++) {
0289       const L1MuDTTrack* cand = m_tf.sp(tmpspid)->track(number);
0290       const L1MuDTTrack* canD = m_tf.sp(tmpspid)->tracK(number);
0291       if (cand) {
0292         m_address[i] = cand->address().trackAddressCode();
0293         if (!cand->empty()) {
0294           m_TrackCand[i] = const_cast<L1MuDTTrack*>(cand);
0295           m_TracKCand[i] = const_cast<L1MuDTTrack*>(canD);
0296         }
0297         i++;
0298       }
0299     }
0300   }
0301 }
0302 
0303 //
0304 // run Eta Track Finder (ETF)
0305 //
0306 void L1MuDTEtaProcessor::runEtaTrackFinder(const edm::EventSetup& c) {
0307   theEtaPatternLUT = c.getHandle(theEtaToken);
0308 
0309   // check if there are any data
0310   bool empty = true;
0311   for (int i = 0; i < 15; i++) {
0312     empty &= m_tseta[i]->empty();
0313   }
0314   if (empty)
0315     return;
0316 
0317   // Pattern comparator:
0318   // loop over all patterns and compare with local chamber pattern
0319   // result : list of valid pattern IDs ( m_foundPattern )
0320   L1MuDTEtaPatternLut::ETFLut_iter it = theEtaPatternLUT->begin();
0321   while (it != theEtaPatternLUT->end()) {
0322     const L1MuDTEtaPattern pattern = (*it).second;
0323     int qualitycode = pattern.quality();
0324 
0325     bool good = true;
0326 
0327     for (int station = 0; station < 3; station++) {
0328       int q = quality(qualitycode, station + 1);
0329       int wheel = pattern.wheel(station + 1);
0330       int bin = pattern.position(station + 1);
0331       if (bin == 0)
0332         continue;
0333       bitset<7> pos = m_tseta[wheel + 2 + station * 5]->position();
0334       bitset<7> qual = m_tseta[wheel + 2 + station * 5]->quality();
0335       // compare position
0336       good &= pos.test(bin - 1);
0337       // compare quality
0338       if (q == 2)
0339         good &= qual.test(bin - 1);
0340     }
0341 
0342     if (good)
0343       m_foundPattern.push_back(pattern.id());
0344 
0345     it++;
0346   }
0347 }
0348 
0349 //
0350 // run Eta Matching Unit (EMU)
0351 //
0352 void L1MuDTEtaProcessor::runEtaMatchingUnit(const edm::EventSetup& c) {
0353   theQualPatternLUT = c.getHandle(theQualToken);
0354 
0355   // loop over all addresses
0356   for (int i = 0; i < 12; i++) {
0357     int adr = m_address[i];
0358     if (adr == 0)
0359       continue;
0360     int sp = i / 2 + 1;  //sector processor [1,6]
0361 
0362     // assign coarse eta value
0363     if (!m_mask)
0364       m_eta[i] = theQualPatternLUT->getCoarseEta(sp, adr);
0365     if (m_eta[i] == 99)
0366       m_eta[i] = 32;
0367     if (m_eta[i] > 31)
0368       m_eta[i] -= 64;
0369     m_eta[i] += 32;
0370 
0371     if (m_foundPattern.empty())
0372       continue;
0373 
0374     // get list of qualified patterns ordered by quality
0375     // and compare with found patterns
0376     const vector<short>& qualifiedPatterns = theQualPatternLUT->getQualifiedPatterns(sp, adr);
0377     vector<short>::const_iterator iter;
0378     vector<int>::const_iterator f_iter;
0379     for (iter = qualifiedPatterns.begin(); iter != qualifiedPatterns.end(); iter++) {
0380       f_iter = find(m_foundPattern.begin(), m_foundPattern.end(), (*iter));
0381       // found
0382       if (f_iter != m_foundPattern.end()) {
0383         const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(*f_iter);
0384         // assign fine eta value
0385         m_fine[i] = true;
0386         m_eta[i] = p.eta();  // improved eta
0387         if (m_eta[i] == 99)
0388           m_eta[i] = 32;
0389         if (m_eta[i] > 31)
0390           m_eta[i] -= 64;
0391         m_eta[i] += 32;
0392         m_pattern[i] = (*f_iter);
0393         break;
0394       }
0395     }
0396   }
0397 
0398   // if both tracks from one sector processor deliver the same track address
0399   // both tracks get only a coarse eta value!
0400 
0401   // loop over sector processors
0402   for (int i = 0; i < 6; i++) {
0403     int idx1 = 2 * i;
0404     int idx2 = 2 * i + 1;
0405     int adr1 = m_address[idx1];
0406     int adr2 = m_address[idx2];
0407     if (adr1 == 0 || adr2 == 0)
0408       continue;
0409     if (adr1 == adr2 && !m_mask) {
0410       // both tracks get coarse (default) eta value
0411       m_eta[idx1] = theQualPatternLUT->getCoarseEta(i + 1, adr1);
0412       if (m_eta[idx1] == 99)
0413         m_eta[idx1] = 32;
0414       if (m_eta[idx1] > 31)
0415         m_eta[idx1] -= 64;
0416       m_eta[idx1] += 32;
0417       m_pattern[idx1] = 0;
0418       m_fine[idx1] = false;
0419       m_eta[idx2] = theQualPatternLUT->getCoarseEta(i + 1, adr2);
0420       if (m_eta[idx2] == 99)
0421         m_eta[idx2] = 32;
0422       if (m_eta[idx2] > 31)
0423         m_eta[idx2] -= 64;
0424       m_eta[idx2] += 32;
0425       m_pattern[idx2] = 0;
0426       m_fine[idx2] = false;
0427     }
0428   }
0429 }
0430 
0431 //
0432 // assign values to track candidates
0433 //
0434 void L1MuDTEtaProcessor::assign() {
0435   for (int i = 0; i < 12; i++) {
0436     if (m_TrackCand[i]) {
0437       if (m_eta[i] != 99) {
0438         m_TrackCand[i]->setEta(m_eta[i]);
0439         m_TracKCand[i]->setEta(m_eta[i]);
0440       } else {
0441         //        if ( i/2 != 2 ) cerr << "L1MuDTEtaProcessor: assign invalid eta" << " " << m_address[i] << endl;
0442       }
0443       if (m_fine[i]) {
0444         m_TrackCand[i]->setFineEtaBit();
0445         m_TracKCand[i]->setFineEtaBit();
0446         // find all contributing track segments
0447         const L1MuDTEtaPattern p = theEtaPatternLUT->getPattern(m_pattern[i]);
0448         vector<const L1MuDTTrackSegEta*> TSeta;
0449         const L1MuDTTrackSegEta* ts = nullptr;
0450         for (int stat = 0; stat < 3; stat++) {
0451           int wh = p.wheel(stat + 1);
0452           int pos = p.position(stat + 1);
0453           if (pos == 0)
0454             continue;
0455           ts = m_tseta[wh + 2 + stat * 5];
0456           TSeta.push_back(ts);
0457         }
0458         m_TrackCand[i]->setTSeta(TSeta);
0459         m_TracKCand[i]->setTSeta(TSeta);
0460       }
0461     }
0462   }
0463 }
0464 
0465 //
0466 // get quality:  id [0,26], stat [1,3]
0467 //
0468 int L1MuDTEtaProcessor::quality(int id, int stat) {
0469   // quality codes as defined in CMS Note 2001/027
0470   // This QualityPatterns are used to have a defined Quality-Identifier for
0471   // all possible found tracks.
0472   // Therefore three integer numbers ( Station 1, 2, 3 from left to right )
0473   // can have numbers like 0, 1 or 2
0474   //    0 ... no hit is given
0475   //    1 ... a LTRG is given
0476   //    2 ... a HTRG is given
0477 
0478   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},
0479                                {1, 1, 0}, {1, 0, 1}, {0, 1, 1}, {2, 1, 0}, {1, 2, 0}, {2, 0, 1}, {1, 0, 2},
0480                                {0, 2, 1}, {0, 1, 2}, {2, 2, 0}, {2, 0, 2}, {0, 2, 2}, {1, 1, 1}, {2, 1, 1},
0481                                {1, 2, 1}, {1, 1, 2}, {2, 2, 1}, {2, 1, 2}, {1, 2, 2}, {2, 2, 2}};
0482 
0483   return qualcode[id][stat - 1];
0484 }