Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-04-29 02:41:11

0001 #include "EventFilter/Utilities/interface/DAQSource.h"
0002 #include "EventFilter/Utilities/interface/DAQSourceModelsFRD.h"
0003 
0004 #include <iostream>
0005 #include <sstream>
0006 #include <sys/types.h>
0007 #include <sys/file.h>
0008 #include <sys/time.h>
0009 #include <unistd.h>
0010 #include <vector>
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
0015 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
0016 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
0017 
0018 #include "DataFormats/TCDS/interface/TCDSRaw.h"
0019 
0020 #include "FWCore/Framework/interface/Event.h"
0021 #include "EventFilter/Utilities/interface/GlobalEventNumber.h"
0022 #include "EventFilter/Utilities/interface/DAQSourceModels.h"
0023 #include "EventFilter/Utilities/interface/DAQSource.h"
0024 
0025 #include "EventFilter/Utilities/interface/AuxiliaryMakers.h"
0026 
0027 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
0028 #include "DataFormats/Provenance/interface/EventID.h"
0029 #include "DataFormats/Provenance/interface/Timestamp.h"
0030 #include "EventFilter/Utilities/interface/crc32c.h"
0031 
0032 using namespace edm::streamer;
0033 
0034 void DataModeFRD::readEvent(edm::EventPrincipal& eventPrincipal) {
0035   std::unique_ptr<FEDRawDataCollection> rawData(new FEDRawDataCollection);
0036   bool tcdsInRange;
0037   unsigned char* tcds_pointer = nullptr;
0038   edm::Timestamp tstamp = fillFEDRawDataCollection(*rawData, tcdsInRange, tcds_pointer);
0039 
0040   if (daqSource_->useL1EventID()) {
0041     uint32_t L1EventID = event_->event();
0042     edm::EventID eventID = edm::EventID(daqSource_->eventRunNumber(), daqSource_->currentLumiSection(), L1EventID);
0043     edm::EventAuxiliary aux(
0044         eventID, daqSource_->processGUID(), tstamp, event_->isRealData(), edm::EventAuxiliary::PhysicsTrigger);
0045     aux.setProcessHistoryID(daqSource_->processHistoryID());
0046     daqSource_->makeEventWrapper(eventPrincipal, aux);
0047   } else if (tcds_pointer == nullptr) {
0048     uint32_t L1EventID = event_->event();
0049     throw cms::Exception("DAQSource::read") << "No TCDS FED in event with FEDHeader EID -: " << L1EventID;
0050   } else {
0051     const FEDHeader fedHeader(tcds_pointer);
0052     tcds::Raw_v1 const* tcds = reinterpret_cast<tcds::Raw_v1 const*>(tcds_pointer + FEDHeader::length);
0053     edm::EventAuxiliary aux =
0054         evf::evtn::makeEventAuxiliary(tcds,
0055                                       daqSource_->eventRunNumber(),      //TODO_ eventRunNumber_
0056                                       daqSource_->currentLumiSection(),  //currentLumiSection_
0057                                       event_->isRealData(),
0058                                       static_cast<edm::EventAuxiliary::ExperimentType>(fedHeader.triggerType()),
0059                                       daqSource_->processGUID(),
0060                                       !daqSource_->fileListLoopMode(),
0061                                       !tcdsInRange);
0062     aux.setProcessHistoryID(daqSource_->processHistoryID());
0063     daqSource_->makeEventWrapper(eventPrincipal, aux);
0064   }
0065 
0066   std::unique_ptr<edm::WrapperBase> edp(new edm::Wrapper<FEDRawDataCollection>(std::move(rawData)));
0067   eventPrincipal.put(
0068       daqProvenanceHelpers_[0]->productDescription(), std::move(edp), daqProvenanceHelpers_[0]->dummyProvenance());
0069 }
0070 
0071 edm::Timestamp DataModeFRD::fillFEDRawDataCollection(FEDRawDataCollection& rawData,
0072                                                      bool& tcdsInRange,
0073                                                      unsigned char*& tcds_pointer) {
0074   edm::TimeValue_t time;
0075   timeval stv;
0076   gettimeofday(&stv, nullptr);
0077   time = stv.tv_sec;
0078   time = (time << 32) + stv.tv_usec;
0079   edm::Timestamp tstamp(time);
0080 
0081   uint32_t eventSize = event_->eventSize();
0082   unsigned char* event = (unsigned char*)event_->payload();
0083   tcds_pointer = nullptr;
0084   tcdsInRange = false;
0085   uint16_t selectedTCDSFed = 0;
0086   unsigned int fedsInEvent = 0;
0087   while (eventSize > 0) {
0088     assert(eventSize >= FEDTrailer::length);
0089     eventSize -= FEDTrailer::length;
0090     const FEDTrailer fedTrailer(event + eventSize);
0091     const uint32_t fedSize = fedTrailer.fragmentLength() << 3;  //trailer length counts in 8 bytes
0092     assert(eventSize >= fedSize - FEDHeader::length);
0093     eventSize -= (fedSize - FEDHeader::length);
0094     const FEDHeader fedHeader(event + eventSize);
0095     const uint16_t fedId = fedHeader.sourceID();
0096     if (fedId > FEDNumbering::MAXFEDID) {
0097       throw cms::Exception("DAQSource::fillFEDRawDataCollection") << "Out of range FED ID : " << fedId;
0098     } else if (fedId >= MINTCDSuTCAFEDID_ && fedId <= MAXTCDSuTCAFEDID_) {
0099       if (!selectedTCDSFed) {
0100         selectedTCDSFed = fedId;
0101         tcds_pointer = event + eventSize;
0102         if (fedId >= FEDNumbering::MINTCDSuTCAFEDID && fedId <= FEDNumbering::MAXTCDSuTCAFEDID) {
0103           tcdsInRange = true;
0104         }
0105       } else
0106         throw cms::Exception("DAQSource::fillFEDRawDataCollection")
0107             << "Second TCDS FED ID " << fedId << " found. First ID: " << selectedTCDSFed;
0108     }
0109     //take event ID from GTPE FED
0110     FEDRawData& fedData = rawData.FEDData(fedId);
0111     fedData.resize(fedSize);
0112     memcpy(fedData.data(), event + eventSize, fedSize);
0113 
0114     fedsInEvent++;
0115     if (verifyFEDs_ || !expectedFedsInEvent_) {
0116       if (fedIdSet_.find(fedId) == fedIdSet_.end()) {
0117         if (expectedFedsInEvent_)
0118           throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0119               << "FED Id: " << fedId << " was not found in previous events";
0120         else
0121           fedIdSet_.insert(fedId);
0122       }
0123     }
0124   }
0125   assert(eventSize == 0);
0126 
0127   if (!fedsInEvent)
0128     throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0129         << "Event " << event_->event() << " does not contain any FEDs";
0130   else if (!expectedFedsInEvent_) {
0131     expectedFedsInEvent_ = fedsInEvent;
0132     if (fedIdSet_.size() != fedsInEvent)
0133       throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0134           << "First received event: " << event_->event() << " contains duplicate FEDs";
0135   } else if (fedsInEvent != expectedFedsInEvent_)
0136     throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0137         << "Event " << event_->event() << " does not contain same number of FEDs as previous: " << fedsInEvent << "/"
0138         << expectedFedsInEvent_;
0139 
0140   return tstamp;
0141 }
0142 
0143 std::vector<std::shared_ptr<const edm::DaqProvenanceHelper>>& DataModeFRD::makeDaqProvenanceHelpers() {
0144   //set FRD data collection
0145   daqProvenanceHelpers_.clear();
0146   daqProvenanceHelpers_.emplace_back(std::make_shared<const edm::DaqProvenanceHelper>(
0147       edm::TypeID(typeid(FEDRawDataCollection)), "FEDRawDataCollection", "FEDRawDataCollection", "DAQSource"));
0148   return daqProvenanceHelpers_;
0149 }
0150 
0151 void DataModeFRD::makeDataBlockView(unsigned char* addr, RawInputFile* rawFile) {
0152   dataBlockAddr_ = addr;
0153   dataBlockMax_ = rawFile->currentChunkSize();
0154   eventCached_ = false;
0155   nextEventView(rawFile);
0156   eventCached_ = true;
0157 }
0158 
0159 bool DataModeFRD::nextEventView(RawInputFile*) {
0160   if (eventCached_)
0161     return true;
0162   event_ = std::make_unique<FRDEventMsgView>(dataBlockAddr_);
0163   if (event_->size() > dataBlockMax_) {
0164     throw cms::Exception("DAQSource::getNextEvent")
0165         << " event id:" << event_->event() << " lumi:" << event_->lumi() << " run:" << event_->run()
0166         << " of size:" << event_->size() << " bytes does not fit into a chunk of size:" << dataBlockMax_ << " bytes";
0167   }
0168   if (event_->version() < 5)
0169     throw cms::Exception("DAQSource::getNextEvent")
0170         << "Unsupported FRD version " << event_->version() << ". Minimum supported is v5.";
0171   return true;
0172 }
0173 
0174 bool DataModeFRD::checksumValid() {
0175   crc_ = 0;
0176   crc_ = crc32c(crc_, (const unsigned char*)event_->payload(), event_->eventSize());
0177   if (crc_ != event_->crc32c())
0178     return false;
0179   else
0180     return true;
0181 }
0182 
0183 std::string DataModeFRD::getChecksumError() const {
0184   std::stringstream ss;
0185   ss << "Found a wrong crc32c checksum: expected 0x" << std::hex << event_->crc32c() << " but calculated 0x" << crc_;
0186   return ss.str();
0187 }
0188 
0189 /*
0190  * FRD preRead
0191  */
0192 
0193 void DataModeFRDPreUnpack::unpackEvent(edm::streamer::FRDEventMsgView* eview,
0194                                        UnpackedRawEventWrapper* ec,
0195                                        unsigned int ls) {
0196   //TODO: also walk the file and build checksum
0197   FEDRawDataCollection* rawData = new FEDRawDataCollection;
0198   bool tcdsInRange;
0199   unsigned char* tcds_pointer = nullptr;
0200   std::string errmsg;
0201   bool err = false;
0202   edm::Timestamp tstamp = fillFEDRawDataCollection(eview, *rawData, tcdsInRange, tcds_pointer, err, errmsg);
0203   ec->setRawData(rawData);
0204 
0205   uint32_t L1EventID = eview->event();
0206   if (err) {
0207     ec->setError(errmsg);
0208   } else if (daqSource_->useL1EventID()) {
0209     //filelist mode run override not available with this model currently (source sets it too late)
0210     edm::EventID eventID = edm::EventID(eview->run(), ls, L1EventID);
0211     ec->setAux(new edm::EventAuxiliary(
0212         eventID, daqSource_->processGUID(), tstamp, eview->isRealData(), edm::EventAuxiliary::PhysicsTrigger));
0213     ec->aux()->setProcessHistoryID(daqSource_->processHistoryID());
0214   } else if (tcds_pointer == nullptr) {
0215     std::stringstream ss;
0216     ss << "No TCDS FED in event with FEDHeader EID -: " << L1EventID;
0217     ec->setError(ss.str());
0218   } else {
0219     const FEDHeader fedHeader(tcds_pointer);
0220     tcds::Raw_v1 const* tcds = reinterpret_cast<tcds::Raw_v1 const*>(tcds_pointer + FEDHeader::length);
0221     edm::EventAuxiliary* aux = new edm::EventAuxiliary();  //allocate empty aux
0222     *aux = evf::evtn::makeEventAuxiliary(tcds,
0223                                          eview->run(),
0224                                          ls,
0225                                          eview->isRealData(),
0226                                          static_cast<edm::EventAuxiliary::ExperimentType>(fedHeader.triggerType()),
0227                                          daqSource_->processGUID(),
0228                                          !daqSource_->fileListLoopMode(),
0229                                          !tcdsInRange);
0230     ec->setAux(aux);
0231     ec->aux()->setProcessHistoryID(daqSource_->processHistoryID());
0232     ec->setRun(eview->run());
0233   }
0234 }
0235 
0236 void DataModeFRDPreUnpack::readEvent(edm::EventPrincipal& eventPrincipal) {
0237   if (ec_->error())
0238     throw cms::Exception("DAQSource::read") << ec_->errmsg();
0239 
0240   daqSource_->makeEventWrapper(eventPrincipal, *ec_->aux());
0241 
0242   std::unique_ptr<edm::WrapperBase> edp(new edm::Wrapper<FEDRawDataCollection>(std::move(ec_->rawDataRef())));
0243   eventPrincipal.put(
0244       daqProvenanceHelpers_[0]->productDescription(), std::move(edp), daqProvenanceHelpers_[0]->dummyProvenance());
0245 }
0246 
0247 void DataModeFRDPreUnpack::unpackFile(RawInputFile* currentFile) {
0248   const uint64_t fileSize = currentFile->fileSize_;
0249   const unsigned rawHeaderSize = currentFile->rawHeaderSize_;
0250 
0251   //TODO: set threadError for issues in this function
0252   if (rawHeaderSize > 0) {
0253     assert(fileSize >= rawHeaderSize);
0254   }
0255   assert(fileSize >= headerSize());
0256 
0257   uint64_t bufpos = rawHeaderSize;
0258 
0259   while (bufpos < fileSize) {  //loop while there is file/events to read
0260 
0261     assert(bufpos + headerSize() <= fileSize);
0262 
0263     //fit to buffer model
0264     auto dataBlockAddr = (unsigned char*)currentFile->chunks_[0]->buf_ + bufpos;
0265 
0266     //first view for header only, check if it fits
0267     auto eview = std::make_unique<FRDEventMsgView>(dataBlockAddr);
0268 
0269     assert(bufpos + eview->size() <= fileSize);
0270     bufpos += eview->size();
0271 
0272     //create event wrapper
0273     //we will store this per each event queued to fwk
0274     UnpackedRawEventWrapper* ec = new UnpackedRawEventWrapper();
0275 
0276     assert(eview->version() >= 5);
0277 
0278     //crc check
0279     uint32_t crc = crc32c(0, (const unsigned char*)eview->payload(), eview->eventSize());
0280     if (crc != eview->crc32c()) {
0281       std::stringstream ss;
0282       ss << "Found a wrong crc32c checksum: expected 0x" << std::hex << eview->crc32c() << " but calculated 0x" << crc;
0283       ec->setChecksumError(ss.str());
0284       //unpackEvent(eview.get(), ec);
0285     } else
0286       unpackEvent(eview.get(), ec, currentFile->lumi_);
0287     currentFile->queue(ec);
0288   }
0289 }
0290 
0291 edm::Timestamp DataModeFRDPreUnpack::fillFEDRawDataCollection(edm::streamer::FRDEventMsgView* eview,
0292                                                               FEDRawDataCollection& rawData,
0293                                                               bool& tcdsInRange,
0294                                                               unsigned char*& tcds_pointer,
0295                                                               bool& err,
0296                                                               std::string& errmsg) {
0297   edm::TimeValue_t time;
0298   timeval stv;
0299   gettimeofday(&stv, nullptr);
0300   time = stv.tv_sec;
0301   time = (time << 32) + stv.tv_usec;
0302   edm::Timestamp tstamp(time);
0303 
0304   try {
0305     uint32_t eventSize = eview->eventSize();
0306     unsigned char* event = (unsigned char*)eview->payload();
0307     tcds_pointer = nullptr;
0308     tcdsInRange = false;
0309     uint16_t selectedTCDSFed = 0;
0310     unsigned int fedsInEvent = 0;
0311     while (eventSize > 0) {
0312       assert(eventSize >= FEDTrailer::length);
0313       eventSize -= FEDTrailer::length;
0314       const FEDTrailer fedTrailer(event + eventSize);
0315       const uint32_t fedSize = fedTrailer.fragmentLength() << 3;  //trailer length counts in 8 bytes
0316       assert(eventSize >= fedSize - FEDHeader::length);
0317       eventSize -= (fedSize - FEDHeader::length);
0318       const FEDHeader fedHeader(event + eventSize);
0319       const uint16_t fedId = fedHeader.sourceID();
0320       if (fedId > FEDNumbering::MAXFEDID)
0321         throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection") << "Out of range FED ID : " << fedId;
0322       else if (fedId >= MINTCDSuTCAFEDID_ && fedId <= MAXTCDSuTCAFEDID_) {
0323         if (!selectedTCDSFed) {
0324           selectedTCDSFed = fedId;
0325           tcds_pointer = event + eventSize;
0326           if (fedId >= FEDNumbering::MINTCDSuTCAFEDID && fedId <= FEDNumbering::MAXTCDSuTCAFEDID) {
0327             tcdsInRange = true;
0328           }
0329         } else
0330           throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0331               << "Second TCDS FED ID " << fedId << " found. First ID: " << selectedTCDSFed;
0332       }
0333       //take event ID from GTPE FED
0334       FEDRawData& fedData = rawData.FEDData(fedId);
0335       fedData.resize(fedSize);
0336       memcpy(fedData.data(), event + eventSize, fedSize);
0337 
0338       fedsInEvent++;
0339       if (verifyFEDs_ || !expectedFedsInEvent_) {
0340         if (fedIdSet_.find(fedId) == fedIdSet_.end()) {
0341           if (expectedFedsInEvent_)
0342             throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0343                 << "FEDID " << fedId << " was not found in previous events";
0344           else
0345             fedIdSet_.insert(fedId);
0346         }
0347       }
0348     }
0349     assert(eventSize == 0);
0350 
0351     if (!fedsInEvent)
0352       throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0353           << "Event " << event_->event() << " does not contain any FEDs";
0354     else if (!expectedFedsInEvent_) {
0355       expectedFedsInEvent_ = fedsInEvent;
0356       if (fedIdSet_.size() != fedsInEvent)
0357         throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0358             << "First received event: " << event_->event() << " contains duplicate FEDs";
0359     } else if (fedsInEvent != expectedFedsInEvent_)
0360       throw cms::Exception("DataModeFRDPreUnpack:::fillFRDCollection")
0361           << "Event " << event_->event() << " does not contain same number of FEDs as previous: " << fedsInEvent << "/"
0362           << expectedFedsInEvent_;
0363 
0364   } catch (cms::Exception& e) {
0365     err = true;
0366     errmsg = e.what();
0367   }
0368   return tstamp;
0369 }
0370 
0371 std::vector<std::shared_ptr<const edm::DaqProvenanceHelper>>& DataModeFRDPreUnpack::makeDaqProvenanceHelpers() {
0372   //set FRD data collection
0373   daqProvenanceHelpers_.clear();
0374   daqProvenanceHelpers_.emplace_back(std::make_shared<const edm::DaqProvenanceHelper>(
0375       edm::TypeID(typeid(FEDRawDataCollection)), "FEDRawDataCollection", "FEDRawDataCollection", "DAQSource"));
0376   return daqProvenanceHelpers_;
0377 }
0378 
0379 void DataModeFRDPreUnpack::makeDataBlockView(unsigned char* addr, RawInputFile* rawFile) {
0380   dataBlockAddr_ = addr;
0381   dataBlockMax_ = rawFile->currentChunkSize();
0382   eventCached_ = false;
0383   nextEventView(rawFile);
0384   eventCached_ = true;
0385 }
0386 
0387 bool DataModeFRDPreUnpack::nextEventView(RawInputFile* currentFile) {
0388   if (eventCached_)
0389     return true;
0390   event_ = std::make_unique<FRDEventMsgView>(dataBlockAddr_);
0391   if (event_->size() > dataBlockMax_) {
0392     throw cms::Exception("DAQSource::getNextEvent")
0393         << " event id:" << event_->event() << " lumi:" << event_->lumi() << " run:" << event_->run()
0394         << " of size:" << event_->size() << " bytes does not fit into a chunk of size:" << dataBlockMax_ << " bytes";
0395   }
0396 
0397   if (event_->version() < 5)
0398     throw cms::Exception("DAQSource::getNextEvent")
0399         << "Unsupported FRD version " << event_->version() << ". Minimum supported is v5.";
0400 
0401   currentFile->popQueue(ec_);
0402   return true;
0403 }
0404 
0405 bool DataModeFRDPreUnpack::checksumValid() { return !ec_->checksumError(); }
0406 
0407 std::string DataModeFRDPreUnpack::getChecksumError() const { return ec_->errmsg(); }
0408 
0409 /*
0410  * FRD Multi Source
0411  */
0412 
0413 void DataModeFRDStriped::makeDirectoryEntries(std::vector<std::string> const& baseDirs,
0414                                               std::vector<int> const& numSources,
0415                                               std::vector<int> const& sourceIDs,
0416                                               std::string const& sourceIdentifier,
0417                                               std::string const& runDir) {
0418   std::filesystem::path runDirP(runDir);
0419   for (auto& baseDir : baseDirs) {
0420     std::filesystem::path baseDirP(baseDir);
0421     buPaths_.emplace_back(baseDirP / runDirP);
0422   }
0423 }
0424 
0425 void DataModeFRDStriped::readEvent(edm::EventPrincipal& eventPrincipal) {
0426   assert(!events_.empty());
0427   std::unique_ptr<FEDRawDataCollection> rawData(new FEDRawDataCollection);
0428   bool tcdsInRange;
0429   unsigned char* tcds_pointer = nullptr;
0430   edm::Timestamp tstamp = fillFRDCollection(*rawData, tcdsInRange, tcds_pointer);
0431 
0432   auto const& event = events_[0];
0433   if (daqSource_->useL1EventID()) {
0434     uint32_t L1EventID = event->event();
0435     edm::EventID eventID = edm::EventID(daqSource_->eventRunNumber(), daqSource_->currentLumiSection(), L1EventID);
0436     edm::EventAuxiliary aux(
0437         eventID, daqSource_->processGUID(), tstamp, event->isRealData(), edm::EventAuxiliary::PhysicsTrigger);
0438     aux.setProcessHistoryID(daqSource_->processHistoryID());
0439     daqSource_->makeEventWrapper(eventPrincipal, aux);
0440   } else if (tcds_pointer == nullptr) {
0441     uint32_t L1EventID = event->event();
0442     throw cms::Exception("DAQSource::read") << "No TCDS FED in event with FEDHeader EID -: " << L1EventID;
0443   } else {
0444     const FEDHeader fedHeader(tcds_pointer);
0445     tcds::Raw_v1 const* tcds = reinterpret_cast<tcds::Raw_v1 const*>(tcds_pointer + FEDHeader::length);
0446     edm::EventAuxiliary aux =
0447         evf::evtn::makeEventAuxiliary(tcds,
0448                                       daqSource_->eventRunNumber(),
0449                                       daqSource_->currentLumiSection(),
0450                                       event->isRealData(),
0451                                       static_cast<edm::EventAuxiliary::ExperimentType>(fedHeader.triggerType()),
0452                                       daqSource_->processGUID(),
0453                                       !daqSource_->fileListLoopMode(),
0454                                       !tcdsInRange);
0455     aux.setProcessHistoryID(daqSource_->processHistoryID());
0456     daqSource_->makeEventWrapper(eventPrincipal, aux);
0457   }
0458   std::unique_ptr<edm::WrapperBase> edp(new edm::Wrapper<FEDRawDataCollection>(std::move(rawData)));
0459   eventPrincipal.put(
0460       daqProvenanceHelpers_[0]->productDescription(), std::move(edp), daqProvenanceHelpers_[0]->dummyProvenance());
0461   eventCached_ = false;
0462 }
0463 
0464 edm::Timestamp DataModeFRDStriped::fillFRDCollection(FEDRawDataCollection& rawData,
0465                                                      bool& tcdsInRange,
0466                                                      unsigned char*& tcds_pointer) {
0467   edm::TimeValue_t time;
0468   timeval stv;
0469   gettimeofday(&stv, nullptr);
0470   time = stv.tv_sec;
0471   time = (time << 32) + stv.tv_usec;
0472   edm::Timestamp tstamp(time);
0473 
0474   tcds_pointer = nullptr;
0475   tcdsInRange = false;
0476   uint16_t selectedTCDSFed = 0;
0477   int selectedTCDSFileIndex = -1;
0478   unsigned int fedsInEvent = 0;
0479   for (size_t index = 0; index < events_.size(); index++) {
0480     uint32_t eventSize = events_[index]->eventSize();
0481     unsigned char* event = (unsigned char*)events_[index]->payload();
0482     while (eventSize > 0) {
0483       assert(eventSize >= FEDTrailer::length);
0484       eventSize -= FEDTrailer::length;
0485       const FEDTrailer fedTrailer(event + eventSize);
0486       const uint32_t fedSize = fedTrailer.fragmentLength() << 3;  //trailer length counts in 8 bytes
0487       assert(eventSize >= fedSize - FEDHeader::length);
0488       eventSize -= (fedSize - FEDHeader::length);
0489       const FEDHeader fedHeader(event + eventSize);
0490       const uint16_t fedId = fedHeader.sourceID();
0491       if (fedId > FEDNumbering::MAXFEDID) {
0492         throw cms::Exception("DataModeFRDStriped:::fillFRDCollection") << "Out of range FED ID : " << fedId;
0493       } else if (fedId >= MINTCDSuTCAFEDID_ && fedId <= MAXTCDSuTCAFEDID_) {
0494         if (!selectedTCDSFed) {
0495           selectedTCDSFed = fedId;
0496           selectedTCDSFileIndex = index;
0497           tcds_pointer = event + eventSize;
0498           if (fedId >= FEDNumbering::MINTCDSuTCAFEDID && fedId <= FEDNumbering::MAXTCDSuTCAFEDID) {
0499             tcdsInRange = true;
0500           }
0501         } else if (!testing_)
0502           throw cms::Exception("DataModeFRDStriped:::fillFRDCollection")
0503               << "Second TCDS FED ID " << fedId << " found in file at index" << selectedTCDSFileIndex
0504               << ". First ID: " << selectedTCDSFed << " found in file at index " << (uint64_t)index;
0505       }
0506       FEDRawData& fedData = rawData.FEDData(fedId);
0507       fedData.resize(fedSize);
0508       memcpy(fedData.data(), event + eventSize, fedSize);
0509 
0510       fedsInEvent++;
0511       if (verifyFEDs_ || !expectedFedsInEvent_) {
0512         if (fedIdSet_.find(fedId) == fedIdSet_.end()) {
0513           if (expectedFedsInEvent_)
0514             throw cms::Exception("DataModeFRDStriped:::fillFRDCollection")
0515                 << "FEDID " << fedId << " from the file at index " << (uint64_t)index
0516                 << " was not found in previous events";
0517           else
0518             fedIdSet_.insert(fedId);
0519         }
0520       }
0521     }
0522     assert(eventSize == 0);
0523   }
0524 
0525   if (!fedsInEvent)
0526     throw cms::Exception("DataModeFRDStriped:::fillFRDCollection")
0527         << "Event " << events_.at(0)->event() << " does not contain any FEDs";
0528   else if (!expectedFedsInEvent_) {
0529     expectedFedsInEvent_ = fedsInEvent;
0530     if (fedIdSet_.size() != fedsInEvent) {
0531       throw cms::Exception("DataModeFRDStriped:::fillFRDCollection")
0532           << "First received event: " << events_.at(0)->event() << " contains duplicate FEDs";
0533     }
0534   } else if (fedsInEvent != expectedFedsInEvent_)
0535     throw cms::Exception("DataModeFRDStriped:::fillFRDCollection")
0536         << "Event " << events_.at(0)->event() << " does not contain same number of FEDs as previous: " << fedsInEvent
0537         << "/" << expectedFedsInEvent_;
0538 
0539   return tstamp;
0540 }
0541 
0542 std::vector<std::shared_ptr<const edm::DaqProvenanceHelper>>& DataModeFRDStriped::makeDaqProvenanceHelpers() {
0543   //set FRD data collection
0544   daqProvenanceHelpers_.clear();
0545   daqProvenanceHelpers_.emplace_back(std::make_shared<const edm::DaqProvenanceHelper>(
0546       edm::TypeID(typeid(FEDRawDataCollection)), "FEDRawDataCollection", "FEDRawDataCollection", "DAQSource"));
0547   return daqProvenanceHelpers_;
0548 }
0549 
0550 bool DataModeFRDStriped::checksumValid() {
0551   bool status = true;
0552   for (size_t i = 0; i < events_.size(); i++) {
0553     uint32_t crc = 0;
0554     auto const& event = events_[i];
0555     crc = crc32c(crc, (const unsigned char*)event->payload(), event->eventSize());
0556     if (crc != event->crc32c()) {
0557       std::ostringstream ss;
0558       ss << "Found a wrong crc32c checksum at readout index " << i << ": expected 0x" << std::hex << event->crc32c()
0559          << " but calculated 0x" << crc << ". ";
0560       crcMsg_ += ss.str();
0561       status = false;
0562     }
0563   }
0564   return status;
0565 }
0566 
0567 std::string DataModeFRDStriped::getChecksumError() const { return crcMsg_; }
0568 
0569 /*
0570   read multiple input files for this model
0571 */
0572 
0573 std::pair<bool, std::vector<std::string>> DataModeFRDStriped::defineAdditionalFiles(std::string const& primaryName,
0574                                                                                     bool fileListMode) const {
0575   std::vector<std::string> additionalFiles;
0576 
0577   if (fileListMode) {
0578     //for the unit test
0579     additionalFiles.push_back(primaryName + "_1");
0580     return std::make_pair(true, additionalFiles);
0581   }
0582 
0583   auto fullpath = std::filesystem::path(primaryName);
0584   auto fullname = fullpath.filename();
0585 
0586   for (size_t i = 1; i < buPaths_.size(); i++) {
0587     std::filesystem::path newPath = buPaths_[i] / fullname;
0588     additionalFiles.push_back(newPath.generic_string());
0589   }
0590   return std::make_pair(true, additionalFiles);
0591 }
0592 
0593 void DataModeFRDStriped::makeDataBlockView(unsigned char* addr, RawInputFile* rawFile) {
0594   fileHeaderSize_ = rawFile->rawHeaderSize_;
0595   std::vector<uint64_t> const& fileSizes = rawFile->fileSizes_;
0596   numFiles_ = fileSizes.size();
0597   //add offset address for each file payload
0598   dataBlockAddrs_.clear();
0599   dataBlockAddrs_.push_back(addr);
0600   dataBlockMaxAddrs_.clear();
0601   dataBlockMaxAddrs_.push_back(addr + fileSizes[0] - fileHeaderSize_);
0602   auto fileAddr = addr;
0603   for (unsigned int i = 1; i < fileSizes.size(); i++) {
0604     fileAddr += fileSizes[i - 1];
0605     dataBlockAddrs_.push_back(fileAddr);
0606     dataBlockMaxAddrs_.push_back(fileAddr + fileSizes[i] - fileHeaderSize_);
0607   }
0608 
0609   dataBlockMax_ = rawFile->currentChunkSize();
0610   blockCompleted_ = false;
0611   //set event cached as we set initial address here
0612   bool result = makeEvents();
0613   assert(result);
0614   eventCached_ = true;
0615   setDataBlockInitialized(true);
0616 }
0617 
0618 bool DataModeFRDStriped::nextEventView(RawInputFile*) {
0619   blockCompleted_ = false;
0620   if (eventCached_)
0621     return true;
0622   for (unsigned int i = 0; i < events_.size(); i++) {
0623     //add last event length to each stripe
0624     dataBlockAddrs_[i] += events_[i]->size();
0625   }
0626   return makeEvents();
0627 }
0628 
0629 bool DataModeFRDStriped::makeEvents() {
0630   events_.clear();
0631   assert(!blockCompleted_);
0632   int completed = 0;
0633   uint64_t testEvtId = 0;
0634 
0635   for (int i = 0; i < numFiles_; i++) {
0636     if (dataBlockAddrs_[i] >= dataBlockMaxAddrs_[i]) {
0637       //must be exact
0638       assert(dataBlockAddrs_[i] == dataBlockMaxAddrs_[i]);
0639       blockCompleted_ = true;
0640       completed++;
0641       continue;
0642     }
0643     if (blockCompleted_)
0644       continue;
0645     events_.emplace_back(std::make_unique<FRDEventMsgView>(dataBlockAddrs_[i]));
0646 
0647     if (testEvtId == 0)
0648       testEvtId = events_[i]->event();
0649     else if (testEvtId != events_[i]->event())
0650       throw cms::Exception("DAQSource::getNextEvent")
0651           << " event id mismatch:" << events_[i]->event()
0652           << " while in previously parsed RDEventMsgView (other file):" << testEvtId;
0653 
0654     if (dataBlockAddrs_[i] + events_[i]->size() > dataBlockMaxAddrs_[i])
0655       throw cms::Exception("DAQSource::getNextEvent")
0656           << " event id:" << events_[i]->event() << " lumi:" << events_[i]->lumi() << " run:" << events_[i]->run()
0657           << " of size:" << events_[i]->size() << " bytes does not fit into the buffer or has corrupted header";
0658 
0659     if (events_[i]->version() < 5)
0660       throw cms::Exception("DAQSource::getNextEvent")
0661           << "Unsupported FRD version " << events_[i]->version() << ". Minimum supported is v5.";
0662   }
0663   if (completed < numFiles_) {
0664     for (int i = 0; i < numFiles_; i++) {
0665       if (dataBlockAddrs_[i] == dataBlockMaxAddrs_[i]) {
0666         edm::LogError("dataModeFRDStriped::makeEvents") << "incomplete file block read from directory " << buPaths_[i];
0667         errorDetected_ = true;
0668       }
0669     }
0670   }
0671   return !blockCompleted_;
0672 }