File indexing completed on 2023-03-17 11:10:31
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include <IORawData/DTCommissioning/plugins/DTNewROS8FileReader.h>
0012 #include <IORawData/DTCommissioning/plugins/DTFileReaderHelpers.h>
0013
0014 #include <DataFormats/FEDRawData/interface/FEDHeader.h>
0015 #include <DataFormats/FEDRawData/interface/FEDTrailer.h>
0016 #include <DataFormats/FEDRawData/interface/FEDNumbering.h>
0017
0018 #include "DataFormats/Provenance/interface/EventID.h"
0019 #include <DataFormats/Provenance/interface/Timestamp.h>
0020 #include <DataFormats/FEDRawData/interface/FEDRawData.h>
0021 #include <DataFormats/FEDRawData/interface/FEDRawDataCollection.h>
0022
0023 #include <FWCore/ParameterSet/interface/ParameterSet.h>
0024 #include <FWCore/Utilities/interface/Exception.h>
0025
0026 #include <string>
0027 #include <iosfwd>
0028 #include <iostream>
0029 #include <algorithm>
0030
0031 using namespace std;
0032 using namespace edm;
0033
0034 DTNewROS8FileReader::DTNewROS8FileReader(const edm::ParameterSet& pset) : runNumber(1), eventNum(0) {
0035 const string& filename = pset.getUntrackedParameter<string>("fileName");
0036
0037 inputFile.open(filename.c_str());
0038 if (inputFile.fail()) {
0039 throw cms::Exception("InputFileMissing")
0040 << "DTNewROS8FileReader: the input file: " << filename << " is not present";
0041 }
0042
0043 produces<FEDRawDataCollection>();
0044 }
0045
0046 DTNewROS8FileReader::~DTNewROS8FileReader() { inputFile.close(); }
0047
0048 int DTNewROS8FileReader::fillRawData(Event& e,
0049
0050 FEDRawDataCollection*& data) {
0051 EventID eID = e.id();
0052 data = new FEDRawDataCollection();
0053
0054 try {
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085 if (checkEndOfFile())
0086 throw 1;
0087
0088
0089 int numberEventHeadWords = 8;
0090
0091
0092 int numberOfWords = 0;
0093 int nread = 0;
0094 nread = inputFile.read(dataPointer<int>(&numberOfWords), ros8WordLenght);
0095 if (nread <= 0)
0096 throw 1;
0097
0098
0099
0100
0101 int datahead[numberEventHeadWords];
0102 for (int iih = 0; iih < numberEventHeadWords; iih++) {
0103 nread = inputFile.read(dataPointer<int>(&datahead[iih]), ros8WordLenght);
0104 }
0105
0106
0107
0108 int numberOfDataWords = numberOfWords - numberEventHeadWords - 2;
0109 int* eventData = new int[numberOfDataWords];
0110 nread = inputFile.read(dataPointer<int>(eventData), numberOfDataWords * ros8WordLenght);
0111 if (nread <= 0)
0112 throw 1;
0113
0114
0115
0116 int LastCounter = 0;
0117 nread = inputFile.read(dataPointer<int>(&LastCounter), ros8WordLenght);
0118 if (nread <= 0)
0119 throw 1;
0120
0121 if (LastCounter != numberOfWords) {
0122 cout << "[DTNewROS8FileReader]: word counter mismatch exception: " << numberOfWords << " " << LastCounter << endl;
0123 throw 99;
0124 }
0125
0126
0127 runNumber = datahead[0];
0128 cout << "[DTNewROS8FileReader]: Run Number: " << dec << runNumber << endl;
0129
0130
0131
0132 eventNum = datahead[1];
0133
0134 if (eventNum < 1)
0135 eventNum = 1;
0136
0137
0138 eID = EventID(runNumber, 1U, eventNum);
0139
0140
0141
0142
0143 int eventDataSize = numberOfDataWords * ros8WordLenght;
0144 int adjustment = (eventDataSize / 4) % 2 == 1 ? 4 : 0;
0145
0146
0147 FEDRawData& fedRawData = data->FEDData(FEDNumbering::MINDTFEDID);
0148 fedRawData.resize(eventDataSize + adjustment);
0149
0150
0151 copy(reinterpret_cast<unsigned char*>(eventData),
0152 reinterpret_cast<unsigned char*>(eventData) + eventDataSize,
0153 fedRawData.data());
0154
0155
0156 delete[] eventData;
0157
0158 return true;
0159
0160 } catch (int i) {
0161 if (i == 1) {
0162 cout << "[DTNewROS8FileReader]: END OF FILE REACHED. "
0163 << "No information read for the requested event" << endl;
0164 delete data;
0165 data = nullptr;
0166 return false;
0167 } else {
0168 cout << "[DTNewROS8FileReader]: PROBLEM WITH EVENT INFORMATION ON THE FILE. "
0169 << "EVENT DATA READING FAILED code= " << i << endl;
0170 delete data;
0171 data = nullptr;
0172 return false;
0173 }
0174 }
0175 }
0176
0177 void DTNewROS8FileReader::produce(Event& e, EventSetup const& es) {
0178 edm::Handle<FEDRawDataCollection> rawdata;
0179 FEDRawDataCollection* fedcoll = nullptr;
0180 fillRawData(e, fedcoll);
0181 std::unique_ptr<FEDRawDataCollection> bare_product(fedcoll);
0182 e.put(std::move(bare_product));
0183 }
0184
0185 bool DTNewROS8FileReader::checkEndOfFile() {
0186 bool retval = false;
0187 if (inputFile.eof())
0188 retval = true;
0189 return retval;
0190 }