Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \file
0002  *
0003  *  $Date: 2010/03/12 10:04:19 $
0004  *  $Revision: 1.22 $
0005  *  \author M. Zanetti
0006  */
0007 
0008 #include <IORawData/DTCommissioning/plugins/DTDDUFileReader.h>
0009 #include <IORawData/DTCommissioning/plugins/DTFileReaderHelpers.h>
0010 
0011 #include <DataFormats/FEDRawData/interface/FEDHeader.h>
0012 #include <DataFormats/FEDRawData/interface/FEDTrailer.h>
0013 #include <DataFormats/FEDRawData/interface/FEDNumbering.h>
0014 
0015 #include "DataFormats/Provenance/interface/EventID.h"
0016 #include <DataFormats/Provenance/interface/Timestamp.h>
0017 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
0018 #include <DataFormats/FEDRawData/interface/FEDRawDataCollection.h>
0019 
0020 #include <FWCore/ParameterSet/interface/ParameterSet.h>
0021 #include <FWCore/Utilities/interface/Exception.h>
0022 
0023 #include <string>
0024 #include <iosfwd>
0025 #include <iostream>
0026 #include <algorithm>
0027 
0028 using namespace std;
0029 using namespace edm;
0030 
0031 DTDDUFileReader::DTDDUFileReader(const edm::ParameterSet& pset) : runNumber(1), eventNumber(1) {
0032   const string& filename = pset.getUntrackedParameter<string>("fileName");
0033 
0034   readFromDMA = pset.getUntrackedParameter<bool>("readFromDMA", true);
0035   numberOfHeaderWords = pset.getUntrackedParameter<int>("numberOfHeaderWords", 10);
0036   skipEvents = pset.getUntrackedParameter<int>("skipEvents", 0);
0037 
0038   inputFile.open(filename.c_str());
0039   if (inputFile.fail()) {
0040     throw cms::Exception("InputFileMissing") << "[DTDDUFileReader]: the input file: " << filename << " is not present";
0041   } else {
0042     cout << "DTDDUFileReader: DaqSource file '" << filename << "' was succesfully opened" << endl;
0043   }
0044 
0045   uint32_t runNumber_tmp = 1;
0046   inputFile.read(dataPointer<uint32_t>(&runNumber_tmp), 4);
0047   runNumber = runNumber_tmp;
0048 
0049   inputFile.ignore(4 * (numberOfHeaderWords - 1));
0050 
0051   if (skipEvents) {
0052     cout << "" << endl;
0053     cout << "   Dear user, please be patient, " << skipEvents << " are being skipped .." << endl;
0054     cout << "" << endl;
0055   }
0056 
0057   produces<FEDRawDataCollection>();
0058 }
0059 
0060 DTDDUFileReader::~DTDDUFileReader() { inputFile.close(); }
0061 
0062 int DTDDUFileReader::fillRawData(Event& e,
0063                                  //              Timestamp& tstamp,
0064                                  FEDRawDataCollection*& data) {
0065   EventID eID = e.id();
0066   data = new FEDRawDataCollection();
0067 
0068   vector<uint64_t> eventData;
0069   size_t estimatedEventDimension = 102400;  // dimensione hardcoded
0070   eventData.reserve(estimatedEventDimension);
0071   uint64_t word = 0;
0072 
0073   bool haederTag = false;
0074   bool dataTag = true;
0075   bool headerAlreadyFound = false;
0076 
0077   int wordCount = 0;
0078 
0079   // getting the data word by word from the file
0080   // do it until you get the DDU trailer
0081   while (!isTrailer(word, dataTag, wordCount)) {
0082     //while ( !isTrailer(word) ) {
0083 
0084     if (readFromDMA) {
0085       int nread;
0086       word = dmaUnpack(dataTag, nread);
0087       if (nread <= 0) {
0088         cout << "[DTDDUFileReader]: ERROR! no more words and failed to get the trailer" << endl;
0089         delete data;
0090         data = nullptr;
0091         return false;
0092       }
0093     }
0094 
0095     else {
0096       int nread = inputFile.read(dataPointer<uint64_t>(&word), dduWordLength);
0097       dataTag = false;
0098       if (nread <= 0) {
0099         cout << "[DTDDUFileReader]: ERROR! failed to get the trailer" << endl;
0100         delete data;
0101         data = nullptr;
0102         return false;
0103       }
0104     }
0105 
0106     // get the DDU header
0107     if (!headerAlreadyFound)
0108       if (isHeader(word, dataTag)) {
0109         headerAlreadyFound = true;
0110         haederTag = true;
0111       }
0112 
0113     // from now on fill the eventData with the ROS data
0114     if (haederTag) {
0115       if (readFromDMA) {
0116         // swapping only the 8 byte words
0117         if (dataTag) {
0118           swap(word);
0119         }  // WARNING also the ddu status words have been swapped!
0120            // Control the correct interpretation in DDUUnpacker
0121       }
0122 
0123       eventData.push_back(word);
0124       wordCount++;
0125     }
0126   }
0127 
0128   //     FEDTrailer candidate(reinterpret_cast<const unsigned char*>(&word));
0129   //     cout<<"EventSize: pushed back "<<eventData.size()
0130   //    <<";  From trailer "<<candidate.fragmentLength()<<endl;
0131 
0132   // next event reading will start with meaningless trailer+header from DTLocalDAQ
0133   // those will be skipped automatically when seeking for the DDU header
0134   //if (eventData.size() > estimatedEventDimension) throw 2;
0135 
0136   // Eventually skipping events
0137   if ((int)eventNumber >= skipEvents) {
0138     // Setting the Event ID
0139     eID = EventID(runNumber, 1U, eventNumber);
0140 
0141     // eventDataSize = (Number Of Words)* (Word Size)
0142     int eventDataSize = eventData.size() * dduWordLength;
0143 
0144     if (dduID < 770 || dduID > 775) {
0145       cout << "[DTDDUFileReader]: ERROR. DDU ID out of range. DDU id=" << dduID << endl;
0146       // The FED ID is always the first in the DT range
0147       dduID = FEDNumbering::MINDTFEDID;
0148     }
0149     FEDRawData& fedRawData = data->FEDData(dduID);
0150     fedRawData.resize(eventDataSize);
0151 
0152     copy(reinterpret_cast<unsigned char*>(&eventData[0]),
0153          reinterpret_cast<unsigned char*>(&eventData[0]) + eventDataSize,
0154          fedRawData.data());
0155   }
0156 
0157   return true;
0158 }
0159 
0160 void DTDDUFileReader::produce(Event& e, EventSetup const& es) {
0161   edm::Handle<FEDRawDataCollection> rawdata;
0162   FEDRawDataCollection* fedcoll = nullptr;
0163   fillRawData(e, fedcoll);
0164   std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0165   e.put(std::move(bare_product));
0166 }
0167 
0168 void DTDDUFileReader::swap(uint64_t& word) {
0169   twoNibble64* newWorld = reinterpret_cast<twoNibble64*>(&word);
0170 
0171   uint32_t msBits_tmp = newWorld->msBits;
0172   newWorld->msBits = newWorld->lsBits;
0173   newWorld->lsBits = msBits_tmp;
0174 }
0175 
0176 uint64_t DTDDUFileReader::dmaUnpack(bool& isData, int& nread) {
0177   uint64_t dduWord = 0;
0178 
0179   uint32_t td[4];
0180   // read 4 32-bits word from the file;
0181   nread = inputFile.read(dataPointer<uint32_t>(&td[0]), 4);
0182   nread += inputFile.read(dataPointer<uint32_t>(&td[1]), 4);
0183   nread += inputFile.read(dataPointer<uint32_t>(&td[2]), 4);
0184   nread += inputFile.read(dataPointer<uint32_t>(&td[3]), 4);
0185 
0186   uint32_t data[2] = {0, 0};
0187   // adjust 4 32-bits words  into 2 32-bits words
0188   data[0] |= td[3] & 0x3ffff;
0189   data[0] |= (td[2] << 18) & 0xfffc0000;
0190   data[1] |= (td[2] >> 14) & 0x0f;
0191   data[1] |= (td[1] << 4) & 0x3ffff0;
0192   data[1] |= (td[0] << 22) & 0xffc00000;
0193 
0194   isData = (td[0] >> 10) & 0x01;
0195 
0196   // push_back to a 64 word
0197   dduWord = (uint64_t(data[1]) << 32) | data[0];
0198 
0199   return dduWord;
0200 }
0201 
0202 bool DTDDUFileReader::isHeader(uint64_t word, bool dataTag) {
0203   bool it_is = false;
0204   FEDHeader candidate(reinterpret_cast<const unsigned char*>(&word));
0205   if (candidate.check()) {
0206     // if ( candidate.check() && !dataTag) {
0207     it_is = true;
0208     dduID = candidate.sourceID();
0209     eventNumber++;
0210   }
0211 
0212   return it_is;
0213 }
0214 
0215 bool DTDDUFileReader::isTrailer(uint64_t word, bool dataTag, unsigned int wordCount) {
0216   bool it_is = false;
0217   FEDTrailer candidate(reinterpret_cast<const unsigned char*>(&word));
0218   if (candidate.check()) {
0219     //  if ( candidate.check() && !dataTag) {
0220     //cout<<"[DTDDUFileReader] "<<wordCount<<" - "<<candidate.fragmentLength()<<endl;
0221     if (wordCount == candidate.fragmentLength())
0222       it_is = true;
0223   }
0224 
0225   return it_is;
0226 }
0227 
0228 bool DTDDUFileReader::checkEndOfFile() {
0229   bool retval = false;
0230   if (inputFile.eof())
0231     retval = true;
0232   return retval;
0233 }