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.16 $
0005  *  \author M. Zanetti
0006  */
0007 
0008 #include <IORawData/DTCommissioning/plugins/DTROS8FileReader.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 DTROS8FileReader::DTROS8FileReader(const edm::ParameterSet& pset) : runNum(1), eventNum(0) {
0032   //const string & filename = pset.getParameter<string>("fileName");  // Don't work, it produces an error saying that the DTNewROS8FileReader plugin doesn't exist!!!
0033 
0034   const string& filename = pset.getUntrackedParameter<string>("fileName");
0035 
0036   inputFile.open(filename.c_str());
0037   if (inputFile.fail()) {
0038     throw cms::Exception("InputFileMissing") << "DTROS8FileReader: the input file: " << filename << " is not present";
0039   }
0040 
0041   produces<FEDRawDataCollection>();
0042 }
0043 
0044 DTROS8FileReader::~DTROS8FileReader() { inputFile.close(); }
0045 
0046 int DTROS8FileReader::fillRawData(Event& e,
0047                                   //                  Timestamp& tstamp,
0048                                   FEDRawDataCollection*& data) {
0049   EventID eID = e.id();
0050   data = new FEDRawDataCollection();
0051 
0052   try {
0053     if (checkEndOfFile())
0054       throw 1;
0055 
0056     // Get the total number of words from the 1st word in the payload
0057     int numberOfWords = 0;
0058     int nread = 0;
0059     nread = inputFile.read(dataPointer<int>(&numberOfWords), ros8WordLenght);
0060     if (nread <= 0)
0061       throw 1;
0062 
0063     // Get the event data (all words but the 1st)
0064     int* eventData = new int[numberOfWords];
0065     nread = inputFile.read(dataPointer<int>(eventData + 1), (numberOfWords - 1) * ros8WordLenght);
0066     if (nread <= 0)
0067       throw 1;
0068 
0069     // Check that the event data size corresponds to the 1st word datum
0070     if (numberOfWords <= 0 || eventData[numberOfWords - 1] != numberOfWords) {
0071       cout << "[DTROS8FileReader]: word counter mismatch exception: " << numberOfWords << " "
0072            << eventData[numberOfWords - 1] << endl;
0073       throw 99;
0074     }
0075 
0076     // The header added by the local DAQ occupies 8 words, starting from the 2nd
0077     int* head = eventData + 1;
0078 
0079     /* 
0080       Header word 0: run number
0081       Header word 1: spill number
0082       Header word 2: event number
0083       Header word 3: reserved
0084       Header word 4: ROS data offset
0085       Header word 5: PU data offset
0086       Header word 6: reserved
0087       Header word 7: reserved
0088     */
0089 
0090     // WARNING: the event number is reset at a new spill
0091     eID = EventID(head[0], 1U, head[1] * head[2]);
0092 
0093     // The pointer to the ROS payload (the 1st word being the ROS words counter)
0094     int* rosData = eventData + head[4];
0095 
0096     // The ROS payload size
0097     int eventDataSize = *rosData * ros8WordLenght;
0098     // It has to be a multiple of 8 bytes. if not, adjust the size of the FED payload
0099     int adjustment = (eventDataSize / 4) % 2 == 1 ? 4 : 0;
0100 
0101     //if ( (eventDataSize/4)%2 ) adjustment = 4;
0102 
0103     // The FED ID is always the first in the DT range
0104     FEDRawData& fedRawData = data->FEDData(FEDNumbering::MINDTFEDID);
0105     fedRawData.resize(eventDataSize + adjustment);
0106 
0107     // I pass only the ROS data to the Event
0108     copy(reinterpret_cast<unsigned char*>(rosData),
0109          reinterpret_cast<unsigned char*>(rosData) + eventDataSize,
0110          fedRawData.data());
0111 
0112     // needed to get rid of memory leaks (?)
0113     delete[] eventData;
0114 
0115     return true;
0116   }
0117 
0118   catch (int i) {
0119     if (i == 1) {
0120       cout << "[DTROS8FileReader]: END OF FILE REACHED. "
0121            << "No information read for the requested event" << endl;
0122       delete data;
0123       data = nullptr;
0124       return false;
0125     } else {
0126       cout << "[DTROS8FileReader]: PROBLEM WITH EVENT INFORMATION ON THE FILE. "
0127            << "EVENT DATA READING FAILED  code= " << i << endl;
0128       delete data;
0129       data = nullptr;
0130       return false;
0131     }
0132   }
0133 }
0134 
0135 void DTROS8FileReader::produce(Event& e, EventSetup const& es) {
0136   edm::Handle<FEDRawDataCollection> rawdata;
0137   FEDRawDataCollection* fedcoll = nullptr;
0138   fillRawData(e, fedcoll);
0139   std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0140   e.put(std::move(bare_product));
0141 }
0142 
0143 bool DTROS8FileReader::checkEndOfFile() {
0144   bool retval = false;
0145   if (inputFile.eof())
0146     retval = true;
0147   return retval;
0148 }