CMS 3D CMS Logo

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