CMS 3D CMS Logo

DTDDUFileReader.cc

Go to the documentation of this file.
00001 
00008 #include <IORawData/DTCommissioning/src/DTDDUFileReader.h>
00009 #include <IORawData/DTCommissioning/src/DTFileReaderHelpers.h>
00010 
00011 #include <DataFormats/FEDRawData/interface/FEDHeader.h>
00012 #include <DataFormats/FEDRawData/interface/FEDTrailer.h>
00013 #include <DataFormats/FEDRawData/interface/FEDNumbering.h>
00014 
00015 #include "DataFormats/Provenance/interface/EventID.h"
00016 #include <DataFormats/Provenance/interface/Timestamp.h>
00017 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
00018 #include <DataFormats/FEDRawData/interface/FEDRawDataCollection.h>
00019 
00020 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00021 
00022 #include <string>
00023 #include <iosfwd>
00024 #include <iostream>
00025 #include <algorithm>
00026    
00027 using namespace std;
00028 using namespace edm;
00029 
00030 
00031 DTDDUFileReader::DTDDUFileReader(const edm::ParameterSet& pset) : 
00032   runNumber(1), eventNumber(1) {
00033 
00034   const string & filename = pset.getUntrackedParameter<string>("fileName");
00035 
00036   readFromDMA = pset.getUntrackedParameter<bool>("readFromDMA",true);
00037   numberOfHeaderWords = pset.getUntrackedParameter<int>("numberOfHeaderWords",10);
00038   skipEvents = pset.getUntrackedParameter<int>("skipEvents",0);  
00039 
00040   inputFile.open(filename.c_str());
00041   if( inputFile.fail() ) {
00042     throw cms::Exception("InputFileMissing") 
00043       << "[DTDDUFileReader]: the input file: " << filename <<" is not present";
00044   } else {
00045     cout << "DTDDUFileReader: DaqSource file '" << filename << "' was succesfully opened" << endl;
00046   }
00047   
00048   inputFile.ignore(4*numberOfHeaderWords);
00049   
00050   if (skipEvents) { 
00051     cout<<""<<endl;
00052     cout<<"   Dear user, pleas be patient, "<<skipEvents<<" are being skipped .."<<endl;
00053     cout<<""<<endl;
00054   }
00055 
00056 }
00057 
00058 
00059 DTDDUFileReader::~DTDDUFileReader(){
00060       inputFile.close();
00061 }
00062 
00063 bool DTDDUFileReader::fillRawData(EventID& eID,
00064                                   Timestamp& tstamp, 
00065                                   FEDRawDataCollection*& data){
00066   data = new FEDRawDataCollection();
00067 
00068   vector<uint64_t> eventData;
00069   size_t estimatedEventDimension = 102400; // dimensione hardcoded
00070   eventData.reserve(estimatedEventDimension); 
00071   uint64_t word = 0;
00072 
00073   bool haederTag = false;
00074   bool dataTag = true;
00075   
00076   
00077   int wordCount = 0;
00078   
00079   // getting the data word by word from the file
00080   // do it until you get the DDU trailer
00081   while ( !isTrailer(word, dataTag, wordCount) ) {
00082     //while ( !isTrailer(word) ) { 
00083     
00084     if (readFromDMA) {
00085       int nread;
00086       word = dmaUnpack(dataTag,nread);
00087       if ( nread<=0 ) {
00088         cout<<"[DTDDUFileReader]: ERROR! no more words and failed to get the trailer"<<endl;
00089         delete data; data=0;
00090         return false;
00091       }
00092     }
00093     
00094     else {
00095       int nread = inputFile.read(dataPointer<uint64_t>( &word ), dduWordLength);
00096       dataTag = false;
00097       if ( nread<=0 ) {
00098         cout<<"[DTDDUFileReader]: ERROR! failed to get the trailer"<<endl;
00099         delete data; data=0;
00100         return false;
00101       }
00102     }
00103     
00104     // get the DDU header
00105     if (isHeader(word,dataTag)) haederTag=true;
00106     
00107     // from now on fill the eventData with the ROS data
00108     if (haederTag) {
00109       
00110       if (readFromDMA) {
00111         // swapping only the 8 byte words
00112         if (dataTag) {
00113           swap(word);
00114         } // WARNING also the ddu status words have been swapped!
00115         // Control the correct interpretation in DDUUnpacker
00116       }
00117       
00118       eventData.push_back(word);
00119       wordCount++;
00120     }
00121     
00122   } 
00123   
00124   //     FEDTrailer candidate(reinterpret_cast<const unsigned char*>(&word));
00125   //     cout<<"EventSize: pushed back "<<eventData.size()
00126   //    <<";  From trailer "<<candidate.lenght()<<endl;
00127   
00128   // next event reading will start with meaningless trailer+header from DTLocalDAQ
00129   // those will be skipped automatically when seeking for the DDU header
00130   //if (eventData.size() > estimatedEventDimension) throw 2;
00131   
00132   // Eventually skipping events
00133   if ((int)eventNumber >= skipEvents) {
00134 
00135     // Setting the Event ID
00136     eID = EventID( runNumber, eventNumber);
00137     
00138     // eventDataSize = (Number Of Words)* (Word Size)
00139     int eventDataSize = eventData.size()*dduWordLength;
00140     
00141     // The FED ID is always the first in the DT range
00142     FEDRawData& fedRawData = data->FEDData( FEDNumbering::getDTFEDIds().first );
00143     fedRawData.resize(eventDataSize);
00144     
00145     copy(reinterpret_cast<unsigned char*>(&eventData[0]),
00146          reinterpret_cast<unsigned char*>(&eventData[0]) + eventDataSize, fedRawData.data());
00147 
00148   }
00149 
00150   return true;
00151     
00152 }
00153 
00154 void DTDDUFileReader::swap(uint64_t & word) {
00155   
00156   twoNibble64* newWorld = reinterpret_cast<twoNibble64*>(&word);
00157 
00158   uint32_t msBits_tmp = newWorld->msBits;
00159   newWorld->msBits = newWorld->lsBits;
00160   newWorld->lsBits = msBits_tmp;
00161 }
00162 
00163 
00164 uint64_t DTDDUFileReader::dmaUnpack(bool & isData, int & nread) {
00165   
00166   uint64_t dduWord = 0;
00167 
00168   uint32_t td[4];
00169   // read 4 32-bits word from the file;
00170   nread = inputFile.read(dataPointer<uint32_t>( &td[0] ), 4);
00171   nread += inputFile.read(dataPointer<uint32_t>( &td[1] ), 4);
00172   nread += inputFile.read(dataPointer<uint32_t>( &td[2] ), 4);
00173   nread += inputFile.read(dataPointer<uint32_t>( &td[3] ), 4);
00174 
00175   uint32_t data[2];
00176   // adjust 4 32-bits words  into 2 32-bits words
00177   data[0] |= td[3] & 0x3ffff;
00178   data[0] |= (td[2] << 18 ) & 0xfffc0000;
00179   data[1] |= (td[2] >> 14 ) & 0x0f;
00180   data[1] |= (td[1] << 4 ) & 0x3ffff0;
00181   data[1] |= (td[0] << 22 ) & 0xffc00000;
00182 
00183   isData = ( td[0] >> 10 ) & 0x01;
00184 
00185   // push_back to a 64 word
00186   dduWord = (uint64_t(data[1]) << 32) | data[0];
00187 
00188   return dduWord;
00189 }
00190 
00191 bool DTDDUFileReader::isHeader(uint64_t word, bool dataTag) {
00192 
00193   bool it_is = false;
00194   FEDHeader candidate(reinterpret_cast<const unsigned char*>(&word));
00195   if ( candidate.check() ) {
00196     // if ( candidate.check() && !dataTag) {
00197     it_is = true;
00198    eventNumber++;
00199   }
00200  
00201   return it_is;
00202 }
00203 
00204 
00205 bool DTDDUFileReader::isTrailer(uint64_t word, bool dataTag, int wordCount) {
00206 
00207   bool it_is = false;
00208   FEDTrailer candidate(reinterpret_cast<const unsigned char*>(&word));
00209   if ( candidate.check() ) {
00210     //  if ( candidate.check() && !dataTag) {
00211     //cout<<"[DTDDUFileReader] "<<wordCount<<" - "<<candidate.lenght()<<endl;
00212     if ( wordCount == candidate.lenght())
00213       it_is = true;
00214   }
00215  
00216   return it_is;
00217 }
00218 
00219 
00220 bool DTDDUFileReader::checkEndOfFile(){
00221 
00222   bool retval=false;
00223   if ( inputFile.eof() ) retval=true;
00224   return retval;
00225 
00226 }
00227 
00228 
00229 

Generated on Tue Jun 9 17:39:23 2009 for CMSSW by  doxygen 1.5.4