Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 //-------------------------------------------------
0002 //
0003 //   Class: DTuROSDigiToRaw
0004 //
0005 //   L1 DT uROS Raw-to-Digi
0006 //
0007 //
0008 //
0009 //   Author :
0010 //   J. Troconiz  - UAM
0011 //
0012 //
0013 //--------------------------------------------------
0014 
0015 #include "EventFilter/DTRawToDigi/plugins/DTuROSDigiToRaw.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 
0018 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0019 #include "EventFilter/DTRawToDigi/interface/DTROChainCoding.h"
0020 #include "EventFilter/Utilities/interface/DTCRC.h"
0021 
0022 #include <iostream>
0023 
0024 DTuROSDigiToRaw::DTuROSDigiToRaw(const edm::ParameterSet& pset) : eventNum(0) {
0025   produces<FEDRawDataCollection>();
0026 
0027   DTDigiInputTag_ = pset.getParameter<edm::InputTag>("digiColl");
0028 
0029   debug_ = pset.getUntrackedParameter<bool>("debug", false);
0030 
0031   for (int i = FEDNumbering::MINDTUROSFEDID; i <= FEDNumbering::MAXDTUROSFEDID; i++)
0032     feds_.push_back(i);
0033 
0034   nfeds_ = feds_.size();
0035 
0036   Raw_token = consumes<DTDigiCollection>(DTDigiInputTag_);
0037   mapping_token_ = esConsumes<DTReadOutMapping, DTReadOutMappingRcd>();
0038 }
0039 
0040 DTuROSDigiToRaw::~DTuROSDigiToRaw() {}
0041 
0042 void DTuROSDigiToRaw::produce(edm::Event& e, const edm::EventSetup& c) {
0043   FEDRawDataCollection data;
0044 
0045   if (!fillRawData(e, c, data))
0046     return;
0047 
0048   auto fed_product = std::make_unique<FEDRawDataCollection>(data);
0049 
0050   e.put(std::move(fed_product));
0051 }
0052 
0053 bool DTuROSDigiToRaw::fillRawData(edm::Event& e, const edm::EventSetup& c, FEDRawDataCollection& data) {
0054   eventNum = e.id().event();
0055 
0056   edm::Handle<DTDigiCollection> digis;
0057   e.getByToken(Raw_token, digis);
0058 
0059   edm::ESHandle<DTReadOutMapping> mapping = c.getHandle(mapping_token_);
0060 
0061   for (int w_i = 0; w_i < nfeds_; ++w_i) {
0062     process(feds_[w_i], digis, mapping, data);
0063   }
0064 
0065   return true;
0066 }
0067 
0068 void DTuROSDigiToRaw::process(int DTuROSFED,
0069                               edm::Handle<DTDigiCollection> digis,
0070                               edm::ESHandle<DTReadOutMapping> mapping,
0071                               FEDRawDataCollection& data) {
0072   clear();
0073 
0074   //--> DTDigi analysis
0075 
0076   DTDigiCollection::DigiRangeIterator dttax;
0077   for (dttax = digis->begin(); dttax != digis->end(); ++dttax) {
0078     const DTDigiCollection::Range& dttar = (*dttax).second;
0079     const DTLayerId dttal = (*dttax).first;
0080     for (DTDigiCollection::const_iterator ta = dttar.first; ta != dttar.second; ++ta) {
0081       int wheelId = dttal.wheel();
0082       int sectorId = dttal.sector();
0083       int stationId = dttal.station();
0084       int slId = dttal.superlayer();
0085       int layerId = dttal.layer();
0086       int cellId = (*ta).wire();
0087 
0088       int dduId, rosId, robId, tdcId, tdcChannel;
0089       if (!mapping->geometryToReadOut(
0090               wheelId, stationId, sectorId, slId, layerId, cellId, dduId, rosId, robId, tdcId, tdcChannel)) {
0091         int crate = theCRT(dduId, rosId);
0092 
0093         if (crate != DTuROSFED)
0094           continue;
0095 
0096         int slot = theSLT(dduId, rosId, robId);
0097         int link = theLNK(dduId, rosId, robId);
0098 
0099         int tdcTime = (*ta).countsTDC();
0100 
0101         bslts[slot - 1]++;
0102 
0103         int word = ((link & 0x7F) << 21) + ((tdcId & 0x03) << 19) + ((tdcChannel & 0x1F) << 14) + (tdcTime & 0x3FFF);
0104 
0105         wslts[slot - 1].push_back(word);
0106       }
0107     }
0108   }
0109 
0110   int lines = 4;
0111   int nslts = 0;
0112 
0113   for (int sltit = 0; sltit < DOCESLOTS; sltit++) {
0114     if (bslts[sltit] == 0)
0115       continue;
0116     nslts += 1;
0117     lines += 1;
0118 
0119     dslts[sltit] = ((bslts[sltit] + 1) / 2) + 5;
0120     lines += dslts[sltit];
0121   }
0122 
0123   FEDRawData& dttfdata = data.FEDData(DTuROSFED);
0124   dttfdata.resize(lines * 8);  // size in bytes
0125   unsigned char* lineFED = dttfdata.data();
0126 
0127   int dataWord1, dataWord2;
0128 
0129   //--> Header
0130 
0131   dataWord1 = 0x50000000 + (eventNum & 0xFFFFFF);
0132   dataWord2 = (DTuROSFED & 0xFFF) << 8;
0133 
0134   int newCRC = 0xFFFF;
0135   dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0136 
0137   *((int*)lineFED) = dataWord2;
0138   lineFED += 4;
0139   *((int*)lineFED) = dataWord1;
0140 
0141   //--> AMC sizes
0142 
0143   dataWord1 = (nslts & 0xF) << 20;
0144   dataWord2 = 0;
0145 
0146   dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0147 
0148   lineFED += 4;
0149   *((int*)lineFED) = dataWord2;
0150   lineFED += 4;
0151   *((int*)lineFED) = dataWord1;
0152 
0153   for (int sltit = 0; sltit < DOCESLOTS; sltit++) {
0154     if (bslts[sltit] == 0)
0155       continue;
0156 
0157     dataWord1 = (dslts[sltit] & 0xFFFFFF);
0158     dataWord2 = ((sltit + 1) & 0xF) << 16;
0159 
0160     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0161 
0162     lineFED += 4;
0163     *((int*)lineFED) = dataWord2;
0164     lineFED += 4;
0165     *((int*)lineFED) = dataWord1;
0166   }
0167 
0168   //--> AMC data
0169 
0170   for (int sltit = 0; sltit < DOCESLOTS; sltit++) {
0171     if (bslts[sltit] == 0)
0172       continue;
0173 
0174     dataWord1 = ((sltit + 1) & 0xF) << 24;
0175     dataWord2 = 0;
0176 
0177     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0178 
0179     lineFED += 4;
0180     *((int*)lineFED) = dataWord2;
0181     lineFED += 4;
0182     *((int*)lineFED) = dataWord1;
0183 
0184     dataWord1 = 0;
0185     dataWord2 = 0;
0186 
0187     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0188 
0189     lineFED += 4;
0190     *((int*)lineFED) = dataWord2;
0191     lineFED += 4;
0192     *((int*)lineFED) = dataWord1;
0193 
0194     for (int nhit = 0; nhit < bslts[sltit] / 2; nhit++) {
0195       dataWord1 = 0x20000000 + wslts[sltit].at(nhit * 2);
0196       dataWord2 = wslts[sltit].at(nhit * 2 + 1);
0197 
0198       dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0199 
0200       lineFED += 4;
0201       *((int*)lineFED) = dataWord2;
0202       lineFED += 4;
0203       *((int*)lineFED) = dataWord1;
0204     }
0205 
0206     if (bslts[sltit] % 2 == 1) {
0207       dataWord1 = 0x20000000 + wslts[sltit].at(bslts[sltit] - 1);
0208       dataWord2 = 0x1FFFFFFF;
0209 
0210       dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0211 
0212       lineFED += 4;
0213       *((int*)lineFED) = dataWord2;
0214       lineFED += 4;
0215       *((int*)lineFED) = dataWord1;
0216     }
0217 
0218     dataWord1 = 0x40000000;
0219     dataWord2 = 0;
0220 
0221     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0222 
0223     lineFED += 4;
0224     *((int*)lineFED) = dataWord2;
0225     lineFED += 4;
0226     *((int*)lineFED) = dataWord1;
0227 
0228     dataWord1 = 0x40000000;
0229     dataWord2 = 0;
0230 
0231     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0232 
0233     lineFED += 4;
0234     *((int*)lineFED) = dataWord2;
0235     lineFED += 4;
0236     *((int*)lineFED) = dataWord1;
0237 
0238     dataWord1 = 0;
0239     dataWord2 = (dslts[sltit] & 0xFFFFF);
0240 
0241     dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0242 
0243     lineFED += 4;
0244     *((int*)lineFED) = dataWord2;
0245     lineFED += 4;
0246     *((int*)lineFED) = dataWord1;
0247   }
0248 
0249   //--> Trailer - line 1
0250 
0251   dataWord1 = 0;
0252   dataWord2 = 0;
0253 
0254   dt_crc::calcCRC(dataWord1, dataWord2, newCRC);
0255 
0256   lineFED += 4;
0257   *((int*)lineFED) = dataWord2;
0258   lineFED += 4;
0259   *((int*)lineFED) = dataWord1;
0260 
0261   //--> Trailer - line 2
0262 
0263   dataWord1 = 0xA0000000 + (lines & 0xFFFFFF);
0264   dataWord2 = 0;
0265 
0266   dt_crc::calcCRC(dataWord1, dataWord2 & 0xFFFF, newCRC);
0267 
0268   dataWord2 += (newCRC & 0xFFFF) << 16;
0269 
0270   lineFED += 4;
0271   *((int*)lineFED) = dataWord2;
0272   lineFED += 4;
0273   *((int*)lineFED) = dataWord1;
0274 
0275   return;
0276 }
0277 
0278 void DTuROSDigiToRaw::clear() {
0279   for (int sltit = 0; sltit < DOCESLOTS; sltit++) {
0280     bslts[sltit] = 0;
0281     dslts[sltit] = 0;
0282     wslts[sltit].clear();
0283   }
0284 
0285   return;
0286 }
0287 
0288 int DTuROSDigiToRaw::theCRT(int ddu, int ros) {
0289   if (ros > 6 && ddu > 774)
0290     ddu = ddu - 5;
0291 
0292   if (ddu == 770)
0293     return FEDNumbering::MINDTUROSFEDID;
0294   else if (ddu == 771)
0295     return FEDNumbering::MINDTUROSFEDID;
0296   else if (ddu == 772)
0297     return FEDNumbering::MINDTUROSFEDID + 1;
0298   return FEDNumbering::MAXDTUROSFEDID;
0299 }
0300 
0301 int DTuROSDigiToRaw::theSLT(int ddu, int ros, int rob) {
0302   if (ros > 6 && ddu > 774)
0303     ddu = ddu - 5;
0304 
0305   int slot = ((ros - 1) / 3) + 1;
0306   if (rob == 23)
0307     slot = 5;
0308   if (ddu == 771)
0309     slot += 6;
0310   else if (ddu == 774)
0311     slot += 6;
0312   return slot;
0313 }
0314 
0315 int DTuROSDigiToRaw::theLNK(int ddu, int ros, int rob) {
0316   int link = rob;
0317   if (rob > 14)
0318     link = rob + 1;
0319   if (rob == 24)
0320     link = 15;
0321   link += ((ros - 1) % 3) * 24;
0322   if (rob == 23)
0323     link = ros - 1;
0324   return link;
0325 }
0326 
0327 #include "FWCore/Framework/interface/MakerMacros.h"
0328 DEFINE_FWK_MODULE(DTuROSDigiToRaw);