CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/IORawData/DTCommissioning/src/DTROS8FileReader.cc

Go to the documentation of this file.
00001 
00008 #include <IORawData/DTCommissioning/src/DTROS8FileReader.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 DTROS8FileReader::DTROS8FileReader(const edm::ParameterSet& pset) : 
00033   runNum(1), eventNum(0) {
00034       
00035   const string & filename = pset.getParameter<string>("fileName");
00036 
00037   inputFile.open(filename.c_str());
00038   if( inputFile.fail() ) {
00039     throw cms::Exception("InputFileMissing") 
00040       << "DTROS8FileReader: the input file: " << filename <<" is not present";
00041   }
00042 }
00043 
00044 
00045 DTROS8FileReader::~DTROS8FileReader(){
00046   inputFile.close();
00047 }
00048 
00049 
00050 int DTROS8FileReader::fillRawData(EventID& eID,
00051                                   Timestamp& tstamp, 
00052                                   FEDRawDataCollection*& data){
00053   data = new FEDRawDataCollection();
00054 
00055   try {
00056     
00057 
00058     if( checkEndOfFile() ) throw 1; 
00059     
00060 
00061     // Get the total number of words from the 1st word in the payload
00062     int numberOfWords;
00063     int nread = 0;
00064     nread = inputFile.read(dataPointer<int>( &numberOfWords ), ros8WordLenght);
00065     if ( nread<=0 )  throw 1;
00066 
00067 
00068     // Get the event data (all words but the 1st)
00069     int* eventData = new int[numberOfWords];
00070     nread = inputFile.read(dataPointer<int>( eventData + 1 ), (numberOfWords-1) * ros8WordLenght );
00071     if ( nread<=0 ) throw 1;
00072     
00073 
00074     // Check that the event data size corresponds to the 1st word datum 
00075     if ( eventData[numberOfWords-1] != numberOfWords ) {
00076       cout << "[DTROS8FileReader]: word counter mismatch exception: "
00077            << numberOfWords << " " << eventData[numberOfWords-1] << endl;
00078       throw 99;
00079     }
00080 
00081     // The header added by the local DAQ occupies 8 words, starting from the 2nd 
00082     int* head = eventData + 1;
00083     
00084     /* 
00085       Header word 0: run number
00086       Header word 1: spill number
00087       Header word 2: event number
00088       Header word 3: reserved
00089       Header word 4: ROS data offset
00090       Header word 5: PU data offset
00091       Header word 6: reserved
00092       Header word 7: reserved
00093     */
00094 
00095     // WARNING: the event number is reset at a new spill
00096     eID = EventID( head[0], 1U, head[1]*head[2]);
00097 
00098     // The pointer to the ROS payload (the 1st word being the ROS words counter)
00099     int* rosData = eventData + head[4];
00100 
00101     // The ROS payload size
00102     int eventDataSize = *rosData * ros8WordLenght;
00103     // It has to be a multiple of 8 bytes. if not, adjust the size of the FED payload
00104     int adjustment = (eventDataSize/4)%2 == 1 ? 4 : 0;   
00105     
00106     //if ( (eventDataSize/4)%2 ) adjustment = 4;
00107 
00108 
00109     // The FED ID is always the first in the DT range
00110     FEDRawData& fedRawData = data->FEDData( FEDNumbering::MINDTFEDID );
00111     fedRawData.resize(eventDataSize+adjustment);
00112     
00113     // I pass only the ROS data to the Event
00114     copy(reinterpret_cast<unsigned char*>(rosData), 
00115          reinterpret_cast<unsigned char*>(rosData) + eventDataSize, fedRawData.data());
00116 
00117     // needed to get rid of memory leaks (?)
00118     delete[] eventData;
00119 
00120     return true;
00121   }
00122 
00123   catch( int i ) {
00124 
00125     if ( i == 1 ){
00126       cout << "[DTROS8FileReader]: END OF FILE REACHED. "
00127            << "No information read for the requested event" << endl;
00128       delete data; data=0;
00129       return false;
00130     }    
00131     else {
00132       cout << "[DTROS8FileReader]: PROBLEM WITH EVENT INFORMATION ON THE FILE. "
00133            << "EVENT DATA READING FAILED  code= " << i << endl;
00134       delete data; data=0;
00135       return false;
00136     }
00137 
00138   }
00139 
00140 }
00141 
00142 
00143 bool DTROS8FileReader::checkEndOfFile(){
00144 
00145   bool retval=false;
00146   if ( inputFile.eof() ) retval=true;
00147   return retval;
00148 
00149 }
00150 
00151 
00152