Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:36:08

0001 //-------------------------------------------------
0002 //
0003 //   Class: DTTFFEDReader
0004 //
0005 //   L1 DT Track Finder Raw-to-Digi
0006 //
0007 //
0008 //
0009 //   Author :
0010 //   J. Troconiz  UAM Madrid
0011 //
0012 //--------------------------------------------------
0013 
0014 #include "EventFilter/DTTFRawToDigi/interface/DTTFFEDReader.h"
0015 
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 #include "EventFilter/Utilities/interface/DTCRC.h"
0021 
0022 #include <iostream>
0023 
0024 using namespace std;
0025 
0026 DTTFFEDReader::DTTFFEDReader(const edm::ParameterSet &pset) {
0027   produces<L1MuDTChambPhContainer>();
0028   produces<L1MuDTChambThContainer>();
0029   produces<L1MuDTTrackContainer>("DATA");
0030 
0031   DTTFInputTag = pset.getParameter<edm::InputTag>("DTTF_FED_Source");
0032 
0033   verbose_ = pset.getUntrackedParameter<bool>("verbose", false);
0034 
0035   Raw_token = consumes<FEDRawDataCollection>(DTTFInputTag);
0036 }
0037 
0038 DTTFFEDReader::~DTTFFEDReader() {}
0039 
0040 void DTTFFEDReader::produce(edm::Event &e, const edm::EventSetup &c) {
0041   unique_ptr<L1MuDTChambPhContainer> phi_product(new L1MuDTChambPhContainer);
0042   unique_ptr<L1MuDTChambThContainer> the_product(new L1MuDTChambThContainer);
0043   unique_ptr<L1MuDTTrackContainer> tra_product(new L1MuDTTrackContainer);
0044 
0045   L1MuDTChambPhContainer::Phi_Container phi_data;
0046   L1MuDTChambThContainer::The_Container the_data;
0047   L1MuDTTrackContainer::TrackContainer tra_data;
0048 
0049   if (!fillRawData(e, phi_data, the_data, tra_data))
0050     return;
0051 
0052   phi_product->setContainer(phi_data);
0053   the_product->setContainer(the_data);
0054   tra_product->setContainer(tra_data);
0055 
0056   e.put(std::move(phi_product));
0057   e.put(std::move(the_product));
0058   e.put(std::move(tra_product), "DATA");
0059 }
0060 
0061 bool DTTFFEDReader::fillRawData(edm::Event &e,
0062                                 L1MuDTChambPhContainer::Phi_Container &phi_data,
0063                                 L1MuDTChambThContainer::The_Container &the_data,
0064                                 L1MuDTTrackContainer::TrackContainer &tra_data) {
0065   analyse(e);
0066 
0067   phi_data = p_data();
0068   the_data = t_data();
0069   tra_data = k_data();
0070 
0071   return true;
0072 }
0073 
0074 //--------------
0075 // Operations --
0076 //--------------
0077 void DTTFFEDReader::analyse(edm::Event &e) {
0078   clear();
0079   process(e);
0080   match();
0081   return;
0082 }
0083 
0084 // process data
0085 void DTTFFEDReader::process(edm::Event &e) {
0086   // Container
0087   vector<int> DTTFWordContainer;
0088   vector<int>::iterator DTTFiterator;
0089 
0090   // Header constituents
0091   int BOEevTy, DTTFId;
0092 
0093   // DTTF Payload constituents
0094   int DTTFWord;
0095   int DTTFChan, bitsID;
0096   int addr1[2] = {3, 3};
0097   int addr2[2] = {15, 15};
0098   int addr3[2] = {15, 15};
0099   int addr4[2] = {15, 15};
0100 
0101   // Trailer constituents
0102   int evtLgth, CRC;
0103 
0104   //--> Header
0105 
0106   edm::Handle<FEDRawDataCollection> data;
0107   e.getByToken(Raw_token, data);
0108   FEDRawData dttfdata = data->FEDData(0x030C);
0109   if (dttfdata.size() == 0)
0110     return;
0111 
0112   int *dataWord1 = new int;
0113   int *dataWord2 = new int;
0114   unsigned char *LineFED = dttfdata.data();
0115   *dataWord2 = *((int *)LineFED);
0116   LineFED += 4;
0117   *dataWord1 = *((int *)LineFED);
0118   int lines = 1;  // already counting header
0119 
0120   BOEevTy = ((*dataWord1) & 0xFF000000) >> 24;  // positions 57 ->64
0121   DTTFId = ((*dataWord2) & 0x000FFF00) >> 8;    // positions 9 ->20
0122 
0123   if ((BOEevTy != 0x50) || (DTTFId != 0x030C)) {
0124     if (verbose_)
0125       edm::LogWarning("dttf_unpacker") << "Not a DTTF header " << hex << *dataWord1;
0126     delete dataWord1;
0127     delete dataWord2;
0128     return;
0129   }
0130 
0131   int newCRC = 0xFFFF;
0132   dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0133 
0134   //--> DTTF data
0135 
0136   LineFED += 4;
0137   *dataWord2 = *((int *)LineFED);
0138   LineFED += 4;
0139   *dataWord1 = *((int *)LineFED);
0140   int chkEOE = ((*dataWord1) & 0xFFF00000) >> 20;
0141   lines++;
0142 
0143   while (chkEOE != 0xA00) {
0144     dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0145 
0146     DTTFWord = *dataWord1;
0147     DTTFWordContainer.push_back(DTTFWord);
0148     DTTFWord = *dataWord2;
0149     DTTFWordContainer.push_back(DTTFWord);
0150 
0151     LineFED += 4;
0152     *dataWord2 = *((int *)LineFED);
0153     LineFED += 4;
0154     *dataWord1 = *((int *)LineFED);
0155     chkEOE = ((*dataWord1) & 0xFFF00000) >> 20;
0156     lines++;
0157 
0158     if (lines > 3026) {
0159       if (verbose_)
0160         edm::LogWarning("dttf_unpacker") << "Warning : number of DTTF lines > 3026 ";  // 3026 = 1(header) +
0161                                                                                        // 3024(max # PHTF-ETTF
0162                                                                                        // 64 bits words) +
0163                                                                                        // 1(trailer)
0164       delete dataWord1;
0165       delete dataWord2;
0166       return;
0167     }
0168 
0169   }  // end while-Data loop
0170 
0171   //--> Trailer
0172 
0173   evtLgth = ((*dataWord1) & 0x00FFFFFF);    // positions 33 ->56
0174   CRC = ((*dataWord2) & 0xFFFF0000) >> 16;  // positions 17 ->32
0175 
0176   dt_crc::calcCRC(*dataWord1, (*dataWord2) & 0xFFFF, newCRC);
0177 
0178   if (newCRC != CRC) {
0179     if (verbose_)
0180       edm::LogWarning("dttf_unpacker") << "Calculated CRC " << hex << newCRC << " differs from CRC in trailer " << hex
0181                                        << CRC;
0182     delete dataWord1;
0183     delete dataWord2;
0184     return;
0185   }
0186 
0187   if (lines != evtLgth) {
0188     if (verbose_)
0189       edm::LogWarning("dttf_unpacker") << "Number of words read != event length " << dec << lines << " " << evtLgth;
0190     delete dataWord1;
0191     delete dataWord2;
0192     return;
0193   }
0194 
0195   // --> analyse event
0196 
0197   for (DTTFiterator = DTTFWordContainer.begin(); DTTFiterator != DTTFWordContainer.end(); DTTFiterator++) {
0198     DTTFChan = ((*DTTFiterator) & 0xFF000000) >> 24;
0199     DTTFiterator++;
0200     bitsID = ((*DTTFiterator) & 0xF0000000) >> 28;
0201 
0202     int bxID = bxNr(DTTFChan);
0203     if (bxID == -999)
0204       continue;
0205     int wheelID = wheel(DTTFChan);
0206     if (wheelID == -999)
0207       continue;
0208     int sectorID = sector(DTTFChan);
0209     if (sectorID == -999)
0210       continue;
0211 
0212     // Input
0213     if (wheelID != 0 && bitsID <= 0x9) {
0214       int wheelPh = (abs(wheelID) - 1) * wheelID / abs(wheelID);
0215       int stationID = 0;
0216       int ra = 0;
0217       int ba = 0;
0218       int tsqual = 0;
0219       int ts2tag = 0;
0220 
0221       if ((bitsID >> 1) == 0) {
0222         stationID = 1;
0223       }
0224       if ((bitsID >> 1) == 1) {
0225         stationID = 2;
0226       }
0227       if ((bitsID >> 1) == 4) {
0228         stationID = 3;
0229       }
0230       if ((bitsID >> 1) == 2) {
0231         stationID = 4;
0232       }
0233 
0234       if (stationID != 3) {
0235         ts2tag = (bitsID) & 0x1;
0236         tsqual = (~(*DTTFiterator) & 0x07) - 1;
0237         ba = (~(*DTTFiterator) & 0x1FF8) >> 3;
0238         if (ba > 0x1FF)
0239           ba -= 0x400;
0240         ra = (~(*DTTFiterator) & 0x1FFE000) >> 13;
0241         if (ra > 0x7FF)
0242           ra -= 0x1000;
0243       } else {
0244         ts2tag = (bitsID) & 0x1;
0245         tsqual = (~(*DTTFiterator) & 0x07) - 1;
0246         ra = (~(*DTTFiterator) & 0x7FF8) >> 3;
0247         if (ra > 0x7FF)
0248           ra -= 0x1000;
0249       }
0250 
0251       if (tsqual != 7 && wheelID != -1) {
0252         phiSegments.push_back(
0253             L1MuDTChambPhDigi(bxID + ts2tag, wheelPh, sectorID, stationID, ra, ba, tsqual, ts2tag, 0));
0254       }
0255     }
0256     // Input
0257 
0258     // Input
0259     if (wheelID == 0 && bitsID <= 0x4) {
0260       int wheelTh = bitsID - 2;
0261 
0262       int posALL, posBTI[7];
0263 
0264       if (wheelTh == -2 || wheelTh == -1 ||
0265           (wheelTh == 0 &&
0266            (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11))) {
0267         posALL = ~(*DTTFiterator) & 0x7F;
0268         posBTI[0] = ~(*DTTFiterator) & 0x01;
0269         posBTI[1] = (~(*DTTFiterator) & 0x02) >> 1;
0270         posBTI[2] = (~(*DTTFiterator) & 0x04) >> 2;
0271         posBTI[3] = (~(*DTTFiterator) & 0x08) >> 3;
0272         posBTI[4] = (~(*DTTFiterator) & 0x10) >> 4;
0273         posBTI[5] = (~(*DTTFiterator) & 0x20) >> 5;
0274         posBTI[6] = (~(*DTTFiterator) & 0x40) >> 6;
0275 
0276         if (posALL) {
0277           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 1, posBTI));
0278         }
0279 
0280         posALL = ~(*DTTFiterator) & 0x3F80;
0281         posBTI[0] = (~(*DTTFiterator) & 0x0080) >> 7;
0282         posBTI[1] = (~(*DTTFiterator) & 0x0100) >> 8;
0283         posBTI[2] = (~(*DTTFiterator) & 0x0200) >> 9;
0284         posBTI[3] = (~(*DTTFiterator) & 0x0400) >> 10;
0285         posBTI[4] = (~(*DTTFiterator) & 0x0800) >> 11;
0286         posBTI[5] = (~(*DTTFiterator) & 0x1000) >> 12;
0287         posBTI[6] = (~(*DTTFiterator) & 0x2000) >> 13;
0288 
0289         if (posALL) {
0290           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 2, posBTI));
0291         }
0292 
0293         posALL = ~(*DTTFiterator) & 0x1FC000;
0294         posBTI[0] = (~(*DTTFiterator) & 0x004000) >> 14;
0295         posBTI[1] = (~(*DTTFiterator) & 0x008000) >> 15;
0296         posBTI[2] = (~(*DTTFiterator) & 0x010000) >> 16;
0297         posBTI[3] = (~(*DTTFiterator) & 0x020000) >> 17;
0298         posBTI[4] = (~(*DTTFiterator) & 0x040000) >> 18;
0299         posBTI[5] = (~(*DTTFiterator) & 0x080000) >> 19;
0300         posBTI[6] = (~(*DTTFiterator) & 0x100000) >> 20;
0301 
0302         if (posALL) {
0303           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 3, posBTI));
0304         }
0305       }
0306 
0307       else {
0308         posALL = ~(*DTTFiterator) & 0x7F;
0309         posBTI[6] = ~(*DTTFiterator) & 0x01;
0310         posBTI[5] = (~(*DTTFiterator) & 0x02) >> 1;
0311         posBTI[4] = (~(*DTTFiterator) & 0x04) >> 2;
0312         posBTI[3] = (~(*DTTFiterator) & 0x08) >> 3;
0313         posBTI[2] = (~(*DTTFiterator) & 0x10) >> 4;
0314         posBTI[1] = (~(*DTTFiterator) & 0x20) >> 5;
0315         posBTI[0] = (~(*DTTFiterator) & 0x40) >> 6;
0316 
0317         if (posALL) {
0318           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 1, posBTI));
0319         }
0320 
0321         posALL = ~(*DTTFiterator) & 0x3F80;
0322         posBTI[6] = (~(*DTTFiterator) & 0x0080) >> 7;
0323         posBTI[5] = (~(*DTTFiterator) & 0x0100) >> 8;
0324         posBTI[4] = (~(*DTTFiterator) & 0x0200) >> 9;
0325         posBTI[3] = (~(*DTTFiterator) & 0x0400) >> 10;
0326         posBTI[2] = (~(*DTTFiterator) & 0x0800) >> 11;
0327         posBTI[1] = (~(*DTTFiterator) & 0x1000) >> 12;
0328         posBTI[0] = (~(*DTTFiterator) & 0x2000) >> 13;
0329 
0330         if (posALL) {
0331           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 2, posBTI));
0332         }
0333 
0334         posALL = ~(*DTTFiterator) & 0x1FC000;
0335         posBTI[6] = (~(*DTTFiterator) & 0x004000) >> 14;
0336         posBTI[5] = (~(*DTTFiterator) & 0x008000) >> 15;
0337         posBTI[4] = (~(*DTTFiterator) & 0x010000) >> 16;
0338         posBTI[3] = (~(*DTTFiterator) & 0x020000) >> 17;
0339         posBTI[2] = (~(*DTTFiterator) & 0x040000) >> 18;
0340         posBTI[1] = (~(*DTTFiterator) & 0x080000) >> 19;
0341         posBTI[0] = (~(*DTTFiterator) & 0x100000) >> 20;
0342 
0343         if (posALL) {
0344           theSegments.push_back(L1MuDTChambThDigi(bxID, wheelTh, sectorID, 3, posBTI));
0345         }
0346       }
0347     }
0348     // Input
0349 
0350     // Addresses
0351     if (wheelID != 0 && bitsID >= 0xA && bitsID <= 0xB) {
0352       int candID = bitsID - 0xA;
0353 
0354       addr4[candID] = ((*DTTFiterator) & 0x0F);
0355       addr3[candID] = ((*DTTFiterator) & 0xF0) >> 4;
0356       addr2[candID] = ((*DTTFiterator) & 0xF00) >> 8;
0357       addr1[candID] = ((*DTTFiterator) & 0x3000) >> 12;
0358     }
0359     // Addresses
0360 
0361     // Output
0362     if (wheelID != 0 && bitsID >= 0xC) {
0363       int muonID = 0;
0364       int pt = 0;
0365       int ch = 0;
0366       int phi = 0;
0367       int qual = 0;
0368 
0369       muonID = (bitsID & 0x1);
0370       qual = (~(*DTTFiterator) & 0x07);
0371       phi = ((*DTTFiterator) & 0x7F8) >> 3;
0372       ch = (~(*DTTFiterator) & 0x800) >> 11;
0373       pt = (~(*DTTFiterator) & 0x1F000) >> 12;
0374 
0375       if (qual != 0) {
0376         dtTracks.push_back(L1MuDTTrackCand(0,
0377                                            phi,
0378                                            0,
0379                                            pt,
0380                                            ch,
0381                                            1,
0382                                            0,
0383                                            qual,
0384                                            bxID,
0385                                            wheelID,
0386                                            sectorID,
0387                                            muonID,
0388                                            addr1[muonID],
0389                                            addr2[muonID],
0390                                            addr3[muonID],
0391                                            addr4[muonID]));
0392       }
0393     }
0394     // Output
0395 
0396     // Output
0397     if (wheelID == 0 && bitsID >= 0x8) {
0398       int wheelTh = bitsID & 0x7;
0399 
0400       int etaALL;
0401 
0402       etaALL = ~(*DTTFiterator) & 0x007F;
0403       if (etaALL) {
0404         etTrack[bxID + 1][sectorID][wheelTh][0] = (*DTTFiterator) & 0x003F;
0405         efTrack[bxID + 1][sectorID][wheelTh][0] = (~(*DTTFiterator) & 0x0040) >> 6;
0406       }
0407 
0408       etaALL = (~(*DTTFiterator) & 0x3F80) >> 7;
0409       if (etaALL) {
0410         etTrack[bxID + 1][sectorID][wheelTh][1] = ((*DTTFiterator) & 0x1F80) >> 7;
0411         efTrack[bxID + 1][sectorID][wheelTh][1] = (~(*DTTFiterator) & 0x2000) >> 13;
0412       }
0413     }
0414     // Output
0415 
0416   }  // end for-loop container content
0417 
0418   delete dataWord1;
0419   delete dataWord2;
0420   return;
0421 }
0422 
0423 void DTTFFEDReader::match() {
0424   for (L1MuDTTrackContainer::TrackIterator i = dtTracks.begin(); i != dtTracks.end(); i++) {
0425     int bxTh = i->bx() + 1;
0426     int sectorTh = i->scNum();
0427     int wheelTh = i->whNum() + 3;
0428     if (wheelTh > 3)
0429       wheelTh -= 1;
0430     int muonTh = i->TrkTag();
0431 
0432     i->setEtaPacked(etTrack[bxTh][sectorTh][wheelTh][muonTh]);
0433     i->setFineHaloPacked(efTrack[bxTh][sectorTh][wheelTh][muonTh]);
0434   }
0435 
0436   return;
0437 }
0438 
0439 // access data
0440 const L1MuDTChambPhContainer::Phi_Container &DTTFFEDReader::p_data() { return phiSegments; }
0441 
0442 const L1MuDTChambThContainer::The_Container &DTTFFEDReader::t_data() { return theSegments; }
0443 
0444 const L1MuDTTrackContainer::TrackContainer &DTTFFEDReader::k_data() { return dtTracks; }
0445 
0446 void DTTFFEDReader::clear() {
0447   phiSegments.clear();
0448   theSegments.clear();
0449   dtTracks.clear();
0450 
0451   for (int i = 0; i < 3; i++) {
0452     for (int j = 0; j < 12; j++) {
0453       for (int k = 0; k < 6; k++) {
0454         for (int l = 0; l < 2; l++) {
0455           etTrack[i][j][k][l] = 0;
0456           efTrack[i][j][k][l] = 0;
0457         }
0458       }
0459     }
0460   }
0461 
0462   return;
0463 }
0464 
0465 int DTTFFEDReader::channel(int wheel, int sector, int bx) {
0466   // wheel  :  -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
0467   //                           0 <=> ETTF
0468   // sector :  0 -> 11
0469   // bx     : -1 -> +1
0470 
0471   int myChannel = 255;
0472 
0473   if (abs(bx) > 1) {
0474     return myChannel;
0475   }
0476   if (sector < 0 || sector > 11) {
0477     return myChannel;
0478   }
0479   if (abs(wheel) > 3) {
0480     return myChannel;
0481   }
0482 
0483   myChannel = sector * 21 + wheel * 3 - bx + 10;
0484 
0485   if (myChannel > 125)
0486     myChannel += 2;
0487 
0488   return myChannel;
0489 }
0490 
0491 int DTTFFEDReader::bxNr(int channel) {
0492   int myChannel = channel;
0493 
0494   if (myChannel > 127)
0495     myChannel -= 2;
0496 
0497   if (myChannel < 0 || myChannel > 251) {
0498     return -999;
0499   }
0500 
0501   int myBx = 1 - (myChannel % 3);
0502 
0503   return myBx;
0504 }
0505 
0506 int DTTFFEDReader::sector(int channel) {
0507   int myChannel = channel;
0508 
0509   if (myChannel > 127)
0510     myChannel -= 2;
0511 
0512   if (myChannel < 0 || myChannel > 251) {
0513     return -999;
0514   }
0515 
0516   return myChannel / 21;
0517 }
0518 
0519 int DTTFFEDReader::wheel(int channel) {
0520   int myChannel = channel;
0521 
0522   if (myChannel > 127)
0523     myChannel -= 2;
0524 
0525   if (myChannel < 0 || myChannel > 251) {
0526     return -999;
0527   }
0528 
0529   int myWheel = ((myChannel % 21) / 3) - 3;
0530 
0531   return myWheel;
0532 }