Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:32

0001 //-------------------------------------------------
0002 //
0003 //   Class: DTTFFEDSim
0004 //
0005 //   L1 DT Track Finder Digi-to-Raw
0006 //
0007 //
0008 //
0009 //   Author :
0010 //   J. Troconiz  UAM Madrid
0011 //
0012 //--------------------------------------------------
0013 
0014 #include "EventFilter/DTTFRawToDigi/interface/DTTFFEDSim.h"
0015 
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 
0019 #include "EventFilter/Utilities/interface/DTCRC.h"
0020 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
0021 
0022 #include <iostream>
0023 
0024 using namespace std;
0025 
0026 DTTFFEDSim::DTTFFEDSim(const edm::ParameterSet &pset) : eventNum(0) {
0027   produces<FEDRawDataCollection>();
0028 
0029   DTDigiInputTag = pset.getParameter<edm::InputTag>("DTDigi_Source");
0030   DTPHTFInputTag = pset.getParameter<edm::InputTag>("DTTracks_Source");
0031 
0032   ChPh_tok = consumes<L1MuDTChambPhContainer>(DTDigiInputTag);
0033   ChTh_tok = consumes<L1MuDTChambThContainer>(DTDigiInputTag);
0034   Trk_tok = consumes<L1MuDTTrackContainer>(DTPHTFInputTag);
0035 }
0036 
0037 DTTFFEDSim::~DTTFFEDSim() {}
0038 
0039 void DTTFFEDSim::produce(edm::Event &e, const edm::EventSetup &c) {
0040   FEDRawDataCollection data;
0041 
0042   if (!fillRawData(e, data))
0043     return;
0044 
0045   unique_ptr<FEDRawDataCollection> fed_product(new FEDRawDataCollection(data));
0046 
0047   e.put(std::move(fed_product));
0048 }
0049 
0050 bool DTTFFEDSim::fillRawData(edm::Event &e, FEDRawDataCollection &data) {
0051   eventNum = e.id().event();
0052 
0053   int lines = 2;
0054 
0055   edm::Handle<L1MuDTChambPhContainer> phtrig;
0056   e.getByToken(ChPh_tok, phtrig);
0057   lines += phtrig->bxSize(-1, 1);
0058 
0059   edm::Handle<L1MuDTChambThContainer> thtrig;
0060   e.getByToken(ChTh_tok, thtrig);
0061   lines += thtrig->bxSize(-1, 1);
0062 
0063   edm::Handle<L1MuDTTrackContainer> trtrig;
0064   e.getByToken(Trk_tok, trtrig);
0065   lines += trtrig->bxSize(-1, 1) * 3;
0066 
0067   FEDRawData &dttfdata = data.FEDData(0x30C);
0068   dttfdata.resize(lines * 8);  // size in bytes
0069   unsigned char *LineFED = dttfdata.data();
0070 
0071   int *dataWord1 = new int;
0072   int *dataWord2 = new int;
0073 
0074   //--> Header
0075 
0076   *dataWord1 = 0x50000000 + (eventNum & 0xFFFFFF);
0077   *dataWord2 = 0x00030C00;
0078 
0079   int newCRC = 0xFFFF;
0080   dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0081 
0082   *((int *)LineFED) = *dataWord2;
0083   LineFED += 4;
0084   *((int *)LineFED) = *dataWord1;
0085 
0086   //--> DTTF data
0087 
0088   int TS1Id[4], TS2Id[4];  // word identifier for TS #1,#2 for stations
0089   TS1Id[0] = 0x0E;
0090   TS2Id[0] = 0x1E;
0091   TS1Id[1] = 0x2E;
0092   TS2Id[1] = 0x3E;
0093   TS1Id[3] = 0x4E;
0094   TS2Id[3] = 0x5E;
0095   TS1Id[2] = 0x8FFF8;
0096   TS2Id[2] = 0x9FFF8;
0097 
0098   // Input
0099   L1MuDTChambPhContainer::Phi_iterator tsphi;
0100 
0101   for (tsphi = phtrig->getContainer()->begin(); tsphi != phtrig->getContainer()->end(); tsphi++) {
0102     if (tsphi->code() != 7) {
0103       int wheelID = tsphi->whNum() + 1;
0104       if (wheelID <= 0)
0105         wheelID -= 2;
0106       int stationID = tsphi->stNum() - 1;
0107       int is2nd = tsphi->Ts2Tag();
0108 
0109       int channelNr = channel(wheelID, tsphi->scNum(), tsphi->bxNum() - is2nd);
0110       if (channelNr == 255)
0111         continue;
0112       int TSId = (is2nd == 0) ? TS1Id[stationID] : TS2Id[stationID];
0113 
0114       *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
0115 
0116       if (stationID != 2) {
0117         *dataWord2 = ((TSId & 0x0FF) << 24) + (~(tsphi->code() + 1) & 0x007) + ((~tsphi->phiB() & 0x3FF) << 3) +
0118                      ((~tsphi->phi() & 0xFFF) << 13);
0119       } else {
0120         *dataWord2 = ((TSId & 0xFFFFF) << 12) + (~(tsphi->code() + 1) & 0x00007) + ((~tsphi->phi() & 0x00FFF) << 3);
0121       }
0122 
0123       dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0124 
0125       LineFED += 4;
0126       *((int *)LineFED) = *dataWord2;
0127       LineFED += 4;
0128       *((int *)LineFED) = *dataWord1;
0129     }
0130   }
0131   // Input
0132 
0133   // Input
0134   L1MuDTChambThContainer::The_iterator tsthe;
0135 
0136   for (tsthe = thtrig->getContainer()->begin(); tsthe != thtrig->getContainer()->end(); tsthe++) {
0137     int wheelTh = tsthe->whNum();
0138     int sectorID = tsthe->scNum();
0139 
0140     int channelNr = channel(0, sectorID, tsthe->bxNum());
0141     if (channelNr == 255)
0142       continue;
0143     int TSId = wheelTh + 2;
0144 
0145     *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
0146 
0147     *dataWord2 = ((TSId & 0x07) << 28) + 0x0FFFFFFF;
0148 
0149     int stationID = tsthe->stNum() - 1;
0150     for (int bti = 0; bti < 7; bti++)
0151       if (wheelTh == -2 || wheelTh == -1 ||
0152           (wheelTh == 0 &&
0153            (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11)))
0154         *dataWord2 -= (tsthe->position(bti) & 0x1) << (stationID * 7 + bti);
0155       else
0156         *dataWord2 -= (tsthe->position(6 - bti) & 0x1) << (stationID * 7 + bti);
0157 
0158     dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0159 
0160     LineFED += 4;
0161     *((int *)LineFED) = *dataWord2;
0162     LineFED += 4;
0163     *((int *)LineFED) = *dataWord1;
0164   }
0165   // Input
0166 
0167   // Output
0168   L1MuDTTrackContainer::Trackiterator tstrk;
0169 
0170   for (tstrk = trtrig->getContainer()->begin(); tstrk != trtrig->getContainer()->end(); tstrk++) {
0171     int channelNr = channel(tstrk->whNum(), tstrk->scNum(), tstrk->bx());
0172     if (channelNr == 255)
0173       continue;
0174     int TSId = (tstrk->TrkTag() == 0) ? 0xAFFF : 0xBFFF;
0175 
0176     *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
0177 
0178     *dataWord2 = ((TSId & 0xFFFF) << 16) + (tstrk->stNum(4) & 0x0000F) + ((tstrk->stNum(3) & 0x0000F) << 4) +
0179                  ((tstrk->stNum(2) & 0x0000F) << 8) + ((tstrk->stNum(1) & 0x00003) << 12);
0180 
0181     dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0182 
0183     LineFED += 4;
0184     *((int *)LineFED) = *dataWord2;
0185     LineFED += 4;
0186     *((int *)LineFED) = *dataWord1;
0187 
0188     TSId = (tstrk->TrkTag() == 0) ? 0xCFFE : 0xDFFE;
0189 
0190     *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
0191 
0192     *dataWord2 = ((TSId & 0xFFFE) << 16) + (~tstrk->quality_packed() & 0x0007) + ((tstrk->phi_packed() & 0x00FF) << 3) +
0193                  ((~tstrk->charge_packed() & 0x0001) << 11) + ((~tstrk->pt_packed() & 0x001F) << 12);
0194 
0195     dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0196 
0197     LineFED += 4;
0198     *((int *)LineFED) = *dataWord2;
0199     LineFED += 4;
0200     *((int *)LineFED) = *dataWord1;
0201 
0202     channelNr = channel(0, tstrk->scNum(), tstrk->bx());
0203     if (channelNr == 255)
0204       continue;
0205     TSId = (tstrk->whNum() + 3) << 16;
0206     TSId += (tstrk->whNum() < 0) ? 0x8FFFC : 0x7FFFC;
0207 
0208     *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
0209 
0210     *dataWord2 = (TSId & 0xFFFFC) << 12;
0211 
0212     if (tstrk->TrkTag() == 0) {
0213       *dataWord2 += 0x3F80 + (tstrk->eta_packed() & 0x003F) + ((~tstrk->finehalo_packed() & 0x0001) << 6);
0214     } else {
0215       *dataWord2 += 0x007F + ((tstrk->eta_packed() & 0x003F) << 7) + ((~tstrk->finehalo_packed() & 0x0001) << 13);
0216     }
0217 
0218     dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
0219 
0220     LineFED += 4;
0221     *((int *)LineFED) = *dataWord2;
0222     LineFED += 4;
0223     *((int *)LineFED) = *dataWord1;
0224   }
0225   // Output
0226 
0227   //--> Trailer
0228 
0229   *dataWord1 = 0xA0000000 + (lines & 0xFFFFFF);
0230   *dataWord2 = 0;
0231 
0232   dt_crc::calcCRC(*dataWord1, *dataWord2 & 0xFFFF, newCRC);
0233 
0234   *dataWord2 += (newCRC & 0xFFFF) << 16;
0235 
0236   LineFED += 4;
0237   *((int *)LineFED) = *dataWord2;
0238   LineFED += 4;
0239   *((int *)LineFED) = *dataWord1;
0240 
0241   delete dataWord1;
0242   delete dataWord2;
0243   return true;
0244 }
0245 
0246 int DTTFFEDSim::channel(int wheel, int sector, int bx) {
0247   // wheel  :  -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
0248   //                           0 <=> ETTF
0249   // sector :  0 -> 11
0250   // bx     : -1 -> +1
0251 
0252   int myChannel = 255;
0253 
0254   if (abs(bx) > 1) {
0255     return myChannel;
0256   }
0257   if (sector < 0 || sector > 11) {
0258     return myChannel;
0259   }
0260   if (abs(wheel) > 3) {
0261     return myChannel;
0262   }
0263 
0264   myChannel = sector * 21 + wheel * 3 - bx + 10;
0265 
0266   if (myChannel > 125)
0267     myChannel += 2;
0268 
0269   return myChannel;
0270 }
0271 
0272 int DTTFFEDSim::bxNr(int channel) {
0273   int myChannel = channel;
0274 
0275   if (myChannel > 127)
0276     myChannel -= 2;
0277 
0278   if (myChannel < 0 || myChannel > 251) {
0279     return -999;
0280   }
0281 
0282   int myBx = 1 - (myChannel % 3);
0283 
0284   return myBx;
0285 }
0286 
0287 int DTTFFEDSim::sector(int channel) {
0288   int myChannel = channel;
0289 
0290   if (myChannel > 127)
0291     myChannel -= 2;
0292 
0293   if (myChannel < 0 || myChannel > 251) {
0294     return -999;
0295   }
0296 
0297   return myChannel / 21;
0298 }
0299 
0300 int DTTFFEDSim::wheel(int channel) {
0301   int myChannel = channel;
0302 
0303   if (myChannel > 127)
0304     myChannel -= 2;
0305 
0306   if (myChannel < 0 || myChannel > 251) {
0307     return -999;
0308   }
0309 
0310   int myWheel = ((myChannel % 21) / 3) - 3;
0311 
0312   return myWheel;
0313 }