Project
Loading...
Searching...
No Matches
RawReaderCRU.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15
16#include <fmt/format.h>
17#include <filesystem>
18#include "TSystem.h"
19#include "TObjArray.h"
20
21#include "Headers/DataHeader.h"
23#include "TPCBase/Mapper.h"
24#include "Framework/Logger.h"
27
28#define CHECK_BIT(var, pos) ((var) & (1 << (pos)))
29
30using namespace o2::tpc::rawreader;
31
32std::ostream& operator<<(std::ostream& output, const RDH& rdh);
33
34template <typename DataType>
35std::istream& operator>>(std::istream& input, DataType& data);
36
37void printHeader();
38void printHorizontal(const RDH& rdh);
39
40/*
41// putting this here instead of inside the header trigger unreasonably large compliation times for some reason
42RawReaderCRUEventSync::LinkInfo& RawReaderCRUEventSync::getLinkInfo(uint32_t heartbeatOrbit, int cru, uint8_t globalLinkID)
43{
44 // check if event is already registered. If not create a new one.
45 auto& event = createEvent(heartbeatOrbit);
46 return event.CRUInfoArray[cru].LinkInformation[globalLinkID];
47}
48*/
49
51{
52 // TODO: might be that reversing the loop below has the same effect as using mLastEvent
53 if (mLastEvent && mLastEvent->hasHearbeatOrbit(heartbeatOrbit)) {
54 return *mLastEvent;
55 }
56
57 for (auto& ev : mEventInformation) {
58 const auto hbMatch = ev.hasHearbeatOrbit(heartbeatOrbit);
59 const long hbDiff = long(heartbeatOrbit) - long(ev.HeartbeatOrbits.front());
60 if (hbMatch) {
61 mLastEvent = &ev;
62 return ev;
63 } else if ((hbDiff >= 0) && (hbDiff < 128)) {
64 ev.HeartbeatOrbits.emplace_back(heartbeatOrbit);
65 std::sort(ev.HeartbeatOrbits.begin(), ev.HeartbeatOrbits.end());
66 mLastEvent = &ev;
67 return ev;
68 }
69 }
70 auto& ev = mEventInformation.emplace_back(heartbeatOrbit);
71 mLastEvent = &ev;
72 return ev;
73}
74
76{
77 // expected number of packets in one HBorbit
78 const size_t numberOfPackets = ExpectedNumberOfPacketsPerHBFrame;
79
80 for (int iEvent = mEventInformation.size() - 1; iEvent >= 0; --iEvent) {
81 auto& event = mEventInformation[iEvent];
82 event.IsComplete = true;
83 size_t totalPayloadSize = 0;
84 for (size_t iCRU = 0; iCRU < event.CRUInfoArray.size(); ++iCRU) {
85 const auto& cruInfo = event.CRUInfoArray[iCRU];
86 if (!cruInfo.isPresent()) {
87 if (mCRUSeen[iCRU] >= 0) {
88 LOGP(info, "CRU {} missing in event {}", iCRU, iEvent);
89 event.IsComplete = false;
90 break;
91 }
92
93 continue;
94 } else {
95 totalPayloadSize += cruInfo.totalPayloadSize();
96 }
97
98 if (!cruInfo.isComplete(rawDataType)) {
99 LOGP(info, "CRU info is incomplete");
100 event.IsComplete = false;
101 break;
102 }
103 }
104
105 // remove empty events
106 // can be problems in the filtering in readout
107 // typically these are empty HB frame with HB start and HB stop packets only
108 if (totalPayloadSize == 0) {
109 O2INFO("Removing empty event with HB Orbit %u", event.HeartbeatOrbits[0]);
110 mEventInformation.erase(mEventInformation.begin() + iEvent);
111 }
112 }
113}
114
115void RawReaderCRUEventSync::setLinksSeen(const CRU cru, const std::bitset<MaxNumberOfLinks>& links)
116{
117 for (auto& ev : mEventInformation) {
118 auto& cruInfo = ev.CRUInfoArray[cru];
119 for (int ilink = 0; ilink < cruInfo.LinkInformation.size(); ++ilink) {
120 auto& linkInfo = cruInfo.LinkInformation[ilink];
121 linkInfo.WasSeen = links[ilink];
122 }
123 }
124}
125
127{
128 const std::string redBG("\033[41m");
129 const std::string red("\033[31m");
130 const std::string green("\033[32m");
131 const std::string bold("\033[1m");
132 const std::string clear("\033[0m");
133
134 std::cout << "CRU information";
135 for (size_t iCRU = 0; iCRU < mCRUSeen.size(); ++iCRU) {
136 const auto readerNumber = mCRUSeen[iCRU];
137 if (readerNumber >= 0) {
138 std::cout << fmt::format("CRU {:2} found in reader {}\n", iCRU, readerNumber);
139 }
140 }
141
142 std::cout << "Detailed event information\n";
143 // event loop
144 for (int i = 0; i < mEventInformation.size(); ++i) {
145 const auto& event = mEventInformation[i];
146 const bool isComplete = event.IsComplete;
147 if (!isComplete) {
148 std::cout << redBG;
149 } else {
150 std::cout << green;
151 }
152 std::cout << "Event " << i << " \n"
153 << clear << " heartbeatOrbits: ";
154 for (const auto& orbit : event.HeartbeatOrbits) {
155 std::cout << orbit << " ";
156 }
157 std::cout << "\n"
158 << " firstOrbit: " << event.getFirstOrbit() << "\n"
159 << " Is complete: " << isComplete << "\n";
160
161 // cru loop
162 for (size_t iCRU = 0; iCRU < event.CRUInfoArray.size(); ++iCRU) {
163 const auto& cruInfo = event.CRUInfoArray[iCRU];
164 if (!cruInfo.isPresent()) {
165 continue;
166 }
167 std::cout << " ";
168 if (!cruInfo.isComplete()) {
169 std::cout << bold + red;
170 }
171 std::cout << "CRU " << iCRU << clear << "\n";
172 const auto& cruLinks = cruInfo.LinkInformation;
173
174 // link loop
175 for (size_t iLink = 0; iLink < cruLinks.size(); ++iLink) {
176 const auto& linkInfo = event.CRUInfoArray[iCRU].LinkInformation[iLink];
177 if (!linkInfo.IsPresent) {
178 continue;
179 }
180 std::cout << " ";
181 if (!linkInfo.isComplete()) {
182 std::cout << red;
183 }
184 std::cout << "Link " << iLink << clear << "\n";
185 if (!linkInfo.HBEndSeen) {
186 std::cout << red;
187 }
188 std::cout << " HBEndSeen: " << linkInfo.HBEndSeen << clear << "\n";
189 if (linkInfo.PacketPositions.size() != ExpectedNumberOfPacketsPerHBFrame) {
190 std::cout << red;
191 }
192 std::cout << " Number of Packets: " << linkInfo.PacketPositions.size() << " (" << ExpectedNumberOfPacketsPerHBFrame << ")" << clear << "\n";
193 std::cout << " Payload size : " << linkInfo.PayloadSize << " (" << linkInfo.PayloadSize / 16 << " GBT frames)"
194 << "\n";
195 std::cout << " Packets: ";
196 for (const auto& packet : linkInfo.PacketPositions) {
197 std::cout << packet << " ";
198 }
199 std::cout << "\n";
200 }
201 std::cout << "\n\n";
202 }
203 }
204}
205
206//==============================================================================
208{
209 if (mFileIsScanned) {
210 return 0;
211 }
212
213 // std::vector<PacketDescriptor> mPacketDescriptorMap;
214 // const uint64_t RDH_HEADERWORD0 = 0x1ea04003;
215 // const uint64_t RDH_HEADERWORD0 = 0x00004003;
216 const uint64_t RDH_HEADERWORD0 = 0x00004000; // + RDHUtils::getVersion<o2::header::RAWDataHeader>();
217
218 auto& file = getFileHandle();
219
220 LOGP(info, "scanning file {}", mInputFileName);
221 // get length of file in bytes
222 file.seekg(0, file.end);
223 mFileSize = file.tellg();
224 file.seekg(0, file.beg);
225
226 const bool isTFfile = (mInputFileName.rfind(".tf") == mInputFileName.size() - 3);
227
228 // the file is supposed to contain N x 8kB packets. So the number of packets
229 // can be determined by the file-size. Ideally, this is not required but the
230 // information is derived directly from the header size and payload size.
231 // *** to be adapted to header info ***
232 const uint32_t numPackets = mFileSize / (8 * 1024);
233
234 // read in the RDH, then jump to the next RDH position
235 RDH rdh;
237 uint32_t currentPacket = 0;
238 uint32_t lastHeartbeatOrbit = 0;
239
240 if (isTFfile) {
241 // skip the StfBuilder meta data information
242 file >> dh;
243 file.seekg(dh.payloadSize, std::ios::cur);
244 file >> dh;
245 file.seekg(dh.payloadSize, std::ios::cur);
246 }
247
248 size_t currentPos = file.tellg();
249 size_t dhPayloadSize{};
250 size_t dhPayloadSizeSeen{};
251
252 while ((currentPos < mFileSize) && !file.eof()) {
253 // ===| in case of TF data file read data header |===
254 if (isTFfile && (!dhPayloadSize || (dhPayloadSizeSeen == dhPayloadSize))) {
255 file >> dh;
256 dhPayloadSize = dh.payloadSize;
257 if (dh.dataOrigin != o2::header::gDataOriginTPC) {
258 file.seekg(dhPayloadSize, file.cur);
259 currentPos = file.tellg();
260 dhPayloadSize = 0;
261 continue;
262 }
263 dhPayloadSizeSeen = 0;
264 currentPos = file.tellg();
265 }
266
267 // ===| read in the RawDataHeader at the current position |=================
268 file >> rdh;
269
270 const size_t packetSize = RDHUtils::getOffsetToNext(rdh);
271 const size_t offset = packetSize - RDHUtils::getHeaderSize(rdh);
272 const auto memorySize = RDHUtils::getMemorySize(rdh);
273 const auto payloadSize = memorySize - RDHUtils::getHeaderSize(rdh);
274 dhPayloadSizeSeen += packetSize;
275
276 // ===| check for truncated file |==========================================
277 const size_t curPos = file.tellg();
278 if ((curPos + offset) > mFileSize) {
279 LOGP(error, "File truncated at {}, offset {} would exceed file size of {}", curPos, offset, mFileSize);
280 break;
281 }
282
283 // ===| skip IDC data |=====================================================
284 const auto detField = o2::raw::RDHUtils::getDetectorField(rdh);
285 if (((detField != 0xdeadbeef) && (detField > 2)) || (payloadSize == 0)) {
286 file.seekg(offset, file.cur);
287 ++currentPacket;
288 currentPos = file.tellg();
289 continue;
290 }
291
292 // ===| try to detect data type if not already set |========================
293 //
294 // for now we assume only HB scaling and triggered mode
295 //
296 // in case of triggered data we assume that that the for pageCnt == 1 we have
297 // triggerType == 0x10 in the firt packet
298 //
299 if (mManager) {
300 if (mManager->mDetectDataType) {
301 const uint64_t triggerTypeForTriggeredData = 0x10;
302 const uint64_t triggerType = RDHUtils::getTriggerType(rdh);
303 const uint64_t pageCnt = RDHUtils::getPageCounter(rdh);
304 const uint64_t linkID = RDHUtils::getLinkID(rdh);
305
306 // if (pageCnt == 0) {
307 if ((linkID == 15) || (detField == 0x1) || (detField == 0x2)) {
308 mManager->mRawDataType = RAWDataType::LinkZS;
309 LOGP(info, "Detected LinkZS data");
310 mManager->mDetectDataType = false;
311 }
312 //}
313
314 if (pageCnt == 1) {
315 if (triggerType == triggerTypeForTriggeredData) {
316 mManager->mDataType = DataType::Triggered;
317 O2INFO("Detected triggered data");
318 } else {
319 mManager->mDataType = DataType::HBScaling;
320 O2INFO("Detected HB scaling");
321 }
322 mManager->mDetectDataType = false;
323 }
324 }
325 }
326
327 // ===| get relavant data information |=====================================
328 auto feeId = RDHUtils::getFEEID(rdh);
329 // treat old RDH where feeId was not set properly
330 if (feeId == 4844) {
331 const rdh_utils::FEEIDType cru = RDHUtils::getCRUID(rdh);
332 const rdh_utils::FEEIDType link = RDHUtils::getLinkID(rdh);
333 const rdh_utils::FEEIDType endPoint = RDHUtils::getEndPointID(rdh);
334 feeId = rdh_utils::getFEEID(cru, endPoint, link);
335
336 RDHUtils::setFEEID(rdh, feeId);
337 }
338 const auto heartbeatOrbit = RDHUtils::getHeartBeatOrbit(rdh);
339 const auto heartbeatOrbitEvent = isTFfile ? dh.firstTForbit : RDHUtils::getHeartBeatOrbit(rdh);
340 const auto endPoint = rdh_utils::getEndPoint(feeId);
341 auto linkID = rdh_utils::getLink(feeId);
342 if (linkID == 21 && detField == 0x02) {
343 linkID = 0;
344 }
345 const auto globalLinkID = linkID + endPoint * 12;
346
347 // ===| check if cru should be forced |=====================================
348 if (!mForceCRU) {
349 mCRU = rdh_utils::getCRU(feeId);
350 // mCRU = RDHUtils::getCRUID(rdh); // work-around for MW2 data
351 } else {
352 // overwrite cru id in rdh for further processing
353 RDHUtils::setCRUID(rdh, mCRU);
354 }
355
356 // ===| find evnet info or create a new one |===============================
357 RawReaderCRUEventSync::LinkInfo* linkInfo = nullptr;
358 if (mManager) {
359 // in case of triggered mode, we use the first heartbeat orbit as event identifier
360 if ((lastHeartbeatOrbit == 0) || (heartbeatOrbitEvent != lastHeartbeatOrbit)) {
361 mManager->mEventSync.createEvent(heartbeatOrbitEvent, mManager->getDataType());
362 lastHeartbeatOrbit = heartbeatOrbitEvent;
363 }
364 linkInfo = &mManager->mEventSync.getLinkInfo(rdh, mManager->getDataType());
365 mManager->mEventSync.setCRUSeen(mCRU, mReaderNumber);
366 }
367
368 // ===| set up packet descriptor map for GBT frames |=======================
369 //
370 // * check Header for Header ID
371 // * create the packet descriptor
372 // * set the mLinkPresent flag
373 //
374 if ((rdh.word0 & 0x0000FFF0) == RDH_HEADERWORD0) {
375 // non 0 stop bit means data with payload
376 if (payloadSize) {
377 mPacketDescriptorMaps[globalLinkID].emplace_back(currentPos, mCRU, linkID, endPoint, memorySize, packetSize, heartbeatOrbit);
378 mLinkPresent[globalLinkID] = true;
379 mPacketsPerLink[globalLinkID]++;
380 if (linkInfo) {
381 linkInfo->PacketPositions.emplace_back(mPacketsPerLink[globalLinkID] - 1);
382 linkInfo->IsPresent = true;
383 linkInfo->PayloadSize += payloadSize;
384 }
385 }
386 if (RDHUtils::getStop(rdh) == 1) {
387 // stop bit 1 means we hit the HB end frame without payload.
388 // This marks the end of an "event" in HB scaling mode.
389 if (linkInfo) {
390 linkInfo->HBEndSeen = true;
391 }
392 }
393 } else {
394 O2ERROR("Found header word %x and required header word %x don't match, at %zu, stopping file scan", rdh.word0, RDH_HEADERWORD0, currentPos);
395 break;
396 }
397
398 // debug output
399 if (mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::RDHDump)) {
400 printHorizontal(rdh);
401 if (RDHUtils::getStop(rdh)) {
402 std::cout << "\n";
403 printHeader();
404 }
405 }
406
407 file.seekg(offset, file.cur);
408 ++currentPacket;
409 currentPos = file.tellg();
410 }
411
412 // close the File
413 file.close();
414
415 // go through events and set the status if links were seen
416 if (mManager) {
417 // in case of triggered mode, we use the first heartbeat orbit as event identifier
418 mManager->mEventSync.setLinksSeen(mCRU, mLinkPresent);
419 }
420
421 if (mVerbosity) {
422 // show the mLinkPresent map
423 std::cout << "Links present" << std::endl;
424 for (int i = 0; i < MaxNumberOfLinks; i++) {
425 mLinkPresent[i] == true ? std::cout << "1 " : std::cout << "0 ";
426 };
427 std::cout << '\n';
428
429 std::cout << std::dec
430 << "File Name : " << mInputFileName << "\n"
431 << "File size [bytes] : " << mFileSize << "\n"
432 << "Packets : " << numPackets << "\n"
433 << "\n";
434
435 if (mVerbosity > 1) {
436 // ===| display packet statistics |===
437 for (int i = 0; i < MaxNumberOfLinks; i++) {
438 if (mLinkPresent[i]) {
439 std::cout << "Packets for link " << i << ": " << mPacketsPerLink[i] << "\n";
440 //
441 // ===| display the packet descriptor map |===
442 for (const auto& pd : mPacketDescriptorMaps[i]) {
443 std::cout << pd << "\n";
444 }
445 }
446 }
447 std::cout << "\n";
448 }
449 }
450
451 mFileIsScanned = true;
452
453 return 0;
454}
455
457{
458 auto& file = getFileHandle();
459
460 // loop over the MaxNumberOfLinks potential links in the data
461 // only if data from the link is present and selected
462 // for decoding it will be decoded.
463 for (int link = 0; link < MaxNumberOfLinks; link++) {
464 // all links have been selected
465 if (!checkLinkPresent(link)) {
466 continue;
467 }
468
469 if (mVerbosity) {
470 std::cout << "Finding sync pattern for link " << link << std::endl;
471 std::cout << "Num packets : " << mPacketsPerLink[link] << std::endl;
472 }
473
474 GBTFrame gFrame;
475 uint32_t packetID{0};
476
477 // loop over the packets for each link and process them
478 for (auto packet : mPacketDescriptorMaps[link]) {
479 gFrame.setPacketNumber(packetID);
480
481 file.seekg(packet.getPayloadOffset(), file.beg);
482
483 // read in the data frame by frame, extract the 5-bit halfwords for
484 // the two data streams and store them in the corresponding half-word
485 // vectors
486 for (int frames = 0; frames < packet.getPayloadSize() / 16; frames++) {
487 file >> gFrame;
488 // extract the half words from the 4 32-bit words
489 gFrame.getFrameHalfWords();
490 gFrame.updateSyncCheck(mSyncPositions[link]);
491 };
492
493 packetID++;
494
495 // TODO: In future there might be more then one sync in the stream
496 // this should be takein into account
497 if (syncFoundForLink(link)) {
498 if (mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::SyncPositions)) {
499 std::cout << "Sync positions for link " << link << '\n';
500 const auto& syncs = mSyncPositions[link];
501 for (int i = 0; i < syncs.size(); ++i) {
502 std::cout << i << " : " << syncs[i];
503 }
504 std::cout << '\n';
505 }
506 break;
507 }
508 }
509 }
510}
511
512int RawReaderCRU::processPacket(GBTFrame& gFrame, uint32_t startPos, uint32_t size, ADCRawData& rawData)
513{
514 // open the data file
515 auto& file = getFileHandle();
516
517 // jump to the start position of the packet
518 file.seekg(startPos, file.beg);
519
520 // read in the data frame by frame, extract the 5-bit halfwords for
521 // the two data streams and store them in the corresponding half-word
522 // vectors
523 for (int frames = 0; frames < size / 16; frames++) {
524 file >> gFrame;
525
526 // extract the half words from the 4 32-bit words
527 gFrame.getFrameHalfWords();
528
529 // debug output
530 if (mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::GBTFrames)) {
531 std::cout << gFrame;
532 }
533
534 gFrame.getAdcValues(rawData);
535 gFrame.updateSyncCheck(mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::SyncPositions));
536 if (!(rawData.getNumTimebins() % 16) && (rawData.getNumTimebins() >= mNumTimeBins * 16)) {
537 return 1;
538 }
539 };
540 return 0;
541}
542
543int RawReaderCRU::processMemory(const std::vector<std::byte>& data, ADCRawData& rawData)
544{
545 GBTFrame gFrame;
546
547 const bool dumpSyncPositoins = CHECK_BIT(mDebugLevel, DebugLevel::SyncPositions);
548
549 // 16 bytes is the size of a GBT frame
550 for (int iFrame = 0; iFrame < data.size() / 16; ++iFrame) {
551 gFrame.setFrameNumber(iFrame);
552 gFrame.setPacketNumber(iFrame / 508);
553
554 // in readFromMemory a simple memcopy to the internal data structure is done
555 // I tried using the memory block directly, storing in an internal data member
556 // reinterpret_cast<const uint32_t*>(data.data() + iFrame * 16), so it could be accessed the
557 // same way as the mData array.
558 // however, this was ~5% slower in execution time. I suspect due to cache misses
559 gFrame.readFromMemory(gsl::span<const std::byte>(data.data() + iFrame * 16, 16));
560
561 // extract the half words from the 4 32-bit words
562 gFrame.getFrameHalfWords();
563
564 // debug output
565 if (mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::GBTFrames)) {
566 std::cout << gFrame;
567 }
568
569 gFrame.getAdcValues(rawData);
570 gFrame.updateSyncCheck(mVerbosity && dumpSyncPositoins);
571 if (!(rawData.getNumTimebins() % 16) && (rawData.getNumTimebins() >= mNumTimeBins * 16)) {
572 break;
573 }
574 };
575
576 if (mDumpTextFiles && dumpSyncPositoins) {
577 const auto fileName = mOutputFilePrefix + "/LinkPositions.txt";
578 std::ofstream file(fileName, std::ofstream::app);
579 auto& syncPositions = gFrame.getSyncArray();
580
581 for (int s = 0; s < 5; ++s) {
582 auto& syncPos = syncPositions[s];
583 if (syncPos.synched()) {
584 file << mEventNumber << "\t"
585 << mCRU.number() << "\t"
586 << mLink << "\t"
587 << s << "\t"
588 << syncPos.getPacketNumber() << "\t"
589 << syncPos.getFrameNumber() << "\t"
590 << syncPos.getHalfWordPosition() << "\n";
591 }
592 }
593 }
594 return 0;
595}
596
598{
599 return mManager ? mManager->mEventSync.getNumberOfEvents(mCRU) : 0;
600}
601
602void RawReaderCRU::fillADCdataMap(const ADCRawData& rawData)
603{
604 // TODO: Ugly copy below in runADCDataCallback. Modification in here should be also refected there
605 const auto& mapper = Mapper::instance();
606
607 // cru and link must be set correctly before
608 const CRU cru(mCRU);
609 const int fecLinkOffsetCRU = (mapper.getPartitionInfo(cru.partition()).getNumberOfFECs() + 1) / 2;
610 const int fecInPartition = (mLink % 12) + (mLink > 11) * fecLinkOffsetCRU;
611 const int regionIter = mCRU % 2;
612
613 const int sampaMapping[10] = {0, 0, 1, 1, 2, 3, 3, 4, 4, 2};
614 const int channelOffset[10] = {0, 16, 0, 16, 0, 0, 16, 0, 16, 16};
615
616 for (int istreamm = 0; istreamm < 5; ++istreamm) {
617 const int partitionStream = istreamm + regionIter * 5;
618 const int sampa = sampaMapping[partitionStream];
619
620 const auto& dataVector = rawData.getDataVector(istreamm);
621
622 // loop over all data. Each stream has 16 ADC values for each sampa channel times nTimeBins
623 for (int idata = 0; idata < dataVector.size(); ++idata) {
624 const int ichannel = idata % 16;
625 const int sampaChannel = ichannel + channelOffset[partitionStream];
626 const auto& padPos = mapper.padPosRegion(cru.region(), fecInPartition, sampa, sampaChannel);
627 mADCdata[padPos].emplace_back(dataVector[idata]);
628 }
629 }
630}
631
633{
634 // TODO: Ugly copy below in runADCDataCallback. Modification in here should be also refected there
635 const auto& mapper = Mapper::instance();
636
637 // cru and link must be set correctly before
638 const CRU cru(mCRU);
639 const int fecLinkOffsetCRU = (mapper.getPartitionInfo(cru.partition()).getNumberOfFECs() + 1) / 2;
640 const int fecInPartition = (mLink % 12) + (mLink > 11) * fecLinkOffsetCRU;
641 const int regionIter = mCRU % 2;
642
643 const int sampaMapping[10] = {0, 0, 1, 1, 2, 3, 3, 4, 4, 2};
644 const int channelOffset[10] = {0, 16, 0, 16, 0, 0, 16, 0, 16, 16};
645
646 for (int istreamm = 0; istreamm < 5; ++istreamm) {
647 const int partitionStream = istreamm + regionIter * 5;
648 const int sampa = sampaMapping[partitionStream];
649
650 const auto& dataVector = rawData.getDataVector(istreamm);
651 if (dataVector.size() < 16) {
652 continue;
653 }
654
655 // loop over all data. Each stream has 16 ADC values for each sampa channel times nTimeBins
656 for (int ichannel = 0; ichannel < 16; ++ichannel) {
657 const int sampaChannel = ichannel + channelOffset[partitionStream];
658 const auto& padPos = mapper.padROCPos(cru, fecInPartition, sampa, sampaChannel);
659 // printf("Fill: %d %d %d %d / %d %d %d\n", int(mCRU), int(cru.roc()), ichannel, sampaChannel, int(padPos.getROC()), int(padPos.getRow()), int(padPos.getPad()));
660 mManager->mADCDataCallback(padPos, cru, gsl::span<const uint32_t>(dataVector.data() + ichannel, dataVector.size() - ichannel));
661 }
662 }
663}
664
666{
667 GBTFrame gFrame;
668 // gFrame.setSyncPositions(mSyncPositions[mLink]);
669
670 ADCRawData rawData;
671
672 if (mVerbosity) {
673 std::cout << "Processing data for link " << mLink << std::endl;
674 std::cout << "Num packets : " << mPacketsPerLink[mLink] << std::endl;
675 }
676
677 // ===| mapping to be updated |===============================================
678 // CRU cru; // assuming each decoder only hast once CRU
679 const int link = mLink;
680
681 const auto& linkInfoArray = mManager->mEventSync.getLinkInfoArrayForEvent(mEventNumber, mCRU);
682
683 // loop over the packets for each link and process them
684 // for (const auto& packet : mPacketDescriptorMaps[link]) {
685 for (auto packetNumber : linkInfoArray[link].PacketPositions) {
686 const auto& packet = mPacketDescriptorMaps[link][packetNumber];
687
688 // cru = packet.getCRUID();
689 // std::cout << "Packet : " << packetID << std::endl;
690 gFrame.setPacketNumber(packetNumber);
691 int retCode = processPacket(gFrame, packet.getPayloadOffset(), packet.getPayloadSize(), rawData);
692 if (retCode) {
693 break;
694 }
695
696 // if (rawData.getNumTimebins() >= mNumTimeBins * 16)
697 // break;
698 };
699
700 // ===| fill ADC data to the output structure |===
701 if (mFillADCdataMap) {
702 fillADCdataMap(rawData);
703 }
704 if (mManager && mManager->mADCDataCallback) {
705 runADCDataCallback(rawData);
706 }
707
708 // std::cout << "Output Data" << std::endl;
709
711 if (mDumpTextFiles) {
712 std::ofstream file;
713 std::string fileName;
714 for (int s = 0; s < 5; s++) {
715 if (mStream == 0x0 or ((mStream >> s) & 0x1) == 0x1) {
716 if (gFrame.syncFound(s) == false) {
717 std::cout << "No sync found" << std::endl;
718 }
719 // debug output
720 rawData.setOutputStream(s);
721 rawData.setNumTimebins(mNumTimeBins);
722 if (mVerbosity && CHECK_BIT(mDebugLevel, DebugLevel::ADCValues)) {
723 std::cout << rawData << std::endl;
724 };
725 // write the data to file
726 // fileName = "ADC_" + std::to_string(mLink) + "_" + std::to_string(s);
727 fileName = mOutputFilePrefix + "_ADC_" + std::to_string(mLink) + "_" + std::to_string(s) + ".txt";
728 file.open(fileName, std::ofstream::out);
729 file << rawData;
730 file.close();
731 }
732 }
733 }
734 return 0;
735}
736
738{
739
740 if (mVerbosity) {
741 std::cout << "Processing data for link " << mLink << std::endl;
742 std::cout << "Num packets : " << mPacketsPerLink[mLink] << std::endl;
743 }
744
745 size_t dataSize = 4000 * 16;
746 // if (mDataType == DataType::HBScaling) {
747 // dataSize =
748 // } else if (mDataType == DataType::Triggered) {
751 // dataSize = 4000 * 16;
752 // }
753
754 std::vector<std::byte> data;
755 data.reserve(dataSize);
756 collectGBTData(data);
757
758 ADCRawData rawData;
759 processMemory(data, rawData);
760
761 // ===| fill ADC data to the output structure |===
762 if (mFillADCdataMap) {
763 fillADCdataMap(rawData);
764 }
765 if (mManager->mADCDataCallback) {
766 runADCDataCallback(rawData);
767 }
768}
769
770void RawReaderCRU::collectGBTData(std::vector<std::byte>& data)
771{
772 const auto& linkInfoArray = mManager->mEventSync.getLinkInfoArrayForEvent(mEventNumber, mCRU);
773 auto& file = getFileHandle();
774
775 size_t presentDataPosition = 0;
776
777 // loop over the packets for each link and process them
778 // for (const auto& packet : mPacketDescriptorMaps[link]) {
779 for (auto packetNumber : linkInfoArray[mLink].PacketPositions) {
780 const auto& packet = mPacketDescriptorMaps[mLink][packetNumber];
781
782 const auto payloadStart = packet.getPayloadOffset();
783 const auto payloadSize = size_t(packet.getPayloadSize());
784 data.insert(data.end(), payloadSize, (std::byte)0);
785 // jump to the start position of the packet
786 file.seekg(payloadStart, std::ios::beg);
787
788 // read data
789 file.read(((char*)data.data()) + presentDataPosition, payloadSize);
790
791 presentDataPosition += payloadSize;
792 }
793}
794
796{
797 const auto& eventInfo = mManager->mEventSync.getEventInfo(mEventNumber);
798 const auto& linkInfoArray = eventInfo.CRUInfoArray[mCRU].LinkInformation;
799 const auto firstOrbitInEvent = eventInfo.getFirstOrbit();
800
801 auto& file = getFileHandle();
802
803 char buffer[8192];
804
805 // loop over the packets for each link and process them
806 for (const auto packetNumber : linkInfoArray[mLink].PacketPositions) {
807 const auto& packet = mPacketDescriptorMaps[mLink][packetNumber];
808 const size_t payloadOffset = packet.getPayloadOffset();
809 const size_t payloadSize = packet.getPayloadSize();
810 if ((payloadOffset + payloadSize) > mFileSize) {
811 LOGP(error, "File truncated at {}, size {} would exceed file size of {}", payloadOffset, payloadSize, mFileSize);
812 break;
813 }
814 file.seekg(payloadOffset, file.beg);
815 file.read(buffer, payloadSize);
816 const uint32_t syncOffsetReference = 144; // <<< TODO: fix value as max offset over all links
817 o2::tpc::raw_processing_helpers::processZSdata(buffer, payloadSize, packet.getFEEID(), packet.getHeartBeatOrbit(), firstOrbitInEvent, syncOffsetReference, mManager->mLinkZSCallback); // last parameter should be true for MW2 data
818 }
819}
820
821void RawReaderCRU::processLinks(const uint32_t linkMask)
822{
823 if (!mManager) {
824 LOGP(error, "cannont run without manager");
825 return;
826 }
827
828 try {
829 // read the input data from file, create the packet descriptor map
830 // and the mLinkPresent map.
831 scanFile();
832
833 // check if selected event is valid
834 if (mEventNumber >= mManager->mEventSync.getNumberOfEvents()) {
835 O2ERROR("Selected event number %u is larger than the events in file %lu", mEventNumber, mManager->mEventSync.getNumberOfEvents());
836 return;
837 }
838
839 // loop over the MaxNumberOfLinks potential links in the data
840 // only if data from the link is present and selected
841 // for decoding it will be decoded.
842 for (int lnk = 0; lnk < MaxNumberOfLinks; lnk++) {
843 // all links have been selected
844 if (((linkMask == 0) || ((linkMask >> lnk) & 1)) && checkLinkPresent(lnk) == true) {
845 // set the active link variable and process the data
846 if (mDebugLevel) {
847 fmt::print("Processing link {}\n", lnk);
848 }
849 setLink(lnk);
850 if (mManager->mRawDataType == RAWDataType::GBT) {
851 // processDataFile();
853 } else {
855 }
856 }
857 }
858
859 } catch (const RawReaderCRU::Error& e) {
860 std::cout << e.what() << std::endl;
861 exit(10);
862 } catch (const std::exception& e) {
863 std::cerr << "ERROR: " << e.what() << std::endl;
864 exit(100);
865 } catch (...) {
866 std::cerr << "ERROR: Unknown error" << std::endl;
867 exit(1000);
868 }
869}
870
871void RawReaderCRU::processFile(const std::string_view inputFile, uint32_t timeBins, uint32_t linkMask, uint32_t stream, uint32_t debugLevel, uint32_t verbosity, const std::string_view outputFilePrefix)
872{
873 // Instantiate the RawReaderCRU
874 RawReaderCRUManager cruManager;
875 cruManager.setDebugLevel(debugLevel);
876 // RawReaderCRU& rawReaderCRU = cruManager.createReader(inputFile, timeBins, 0, stream, debugLevel, verbosity, outputFilePrefix);
877 cruManager.setupReaders(inputFile, timeBins, debugLevel, verbosity, outputFilePrefix);
878 for (auto& reader : cruManager.getReaders()) {
879 reader->mDumpTextFiles = true;
880 reader->mFillADCdataMap = false;
881 }
882 cruManager.init();
883
884 if (CHECK_BIT(debugLevel, DebugLevel::SyncPositions)) {
885 std::string outPrefix = outputFilePrefix.length() ? outputFilePrefix.data() : "./";
886 const auto fileName = outPrefix + "/LinkPositions.txt";
887 std::ofstream file(fileName, std::ofstream::out);
888 file << "EventNumber/I"
889 << ":"
890 << "CRU"
891 << ":"
892 << "Link"
893 << ":"
894 << "Stream"
895 << ":"
896 << "PacketNumber"
897 << ":"
898 << "FrameNumber"
899 << ":"
900 << "HalfWordPosition"
901 << "\n";
902 file.close();
903 }
904
905 for (int ievent = 0; ievent < cruManager.getNumberOfEvents(); ++ievent) {
906 fmt::print("=============| event {: 5d} |===============\n", ievent);
907
908 cruManager.processEvent(ievent);
909 }
910}
911
912void RawReaderCRU::copyEvents(const std::vector<uint32_t>& eventNumbers, std::string outputDirectory, std::ios_base::openmode mode)
913{
914 // assemble output file name
915 std::string outputFileName(gSystem->BaseName(mInputFileName.data()));
916 if (outputDirectory.empty()) {
917 outputFileName.insert(0, "filtered.");
918 outputDirectory = gSystem->DirName(mInputFileName.data());
919 }
920 outputFileName.insert(0, "/");
921 outputFileName.insert(0, outputDirectory);
922
923 std::ofstream outputFile(outputFileName, std::ios_base::binary | mode);
924
925 // open the input file
926 auto& file = getFileHandle();
927
928 // data buffer. Maximum size is 8k
929 char buffer[8192];
930
931 // loop over events
932 for (const auto eventNumber : eventNumbers) {
933
934 const auto& linkInfoArray = mManager->mEventSync.getLinkInfoArrayForEvent(eventNumber, mCRU);
935
936 for (int iLink = 0; iLink < MaxNumberOfLinks; ++iLink) {
937 const auto& linkInfo = linkInfoArray[iLink];
938 if (!linkInfo.IsPresent) {
939 continue;
940 }
941 for (auto packetNumber : linkInfo.PacketPositions) {
942 const auto& packet = mPacketDescriptorMaps[iLink][packetNumber];
943 file.seekg(packet.getHeaderOffset(), file.beg);
944 file.read(buffer, packet.getPacketSize());
945 outputFile.write(buffer, packet.getPacketSize());
946 }
947 }
948 }
949}
950
951void RawReaderCRU::writeGBTDataPerLink(std::string_view outputDirectory, int maxEvents)
952{
953 // open the input file
954 auto& file = getFileHandle();
955
956 // data buffer. Maximum size is 8k
957 char buffer[8192];
958
959 // loop over events
960 for (int eventNumber = 0; eventNumber < getNumberOfEvents(); ++eventNumber) {
961 if ((maxEvents > -1) && (eventNumber >= maxEvents)) {
962 break;
963 }
964
965 const auto& linkInfoArray = mManager->mEventSync.getLinkInfoArrayForEvent(eventNumber, mCRU);
966
967 for (int iLink = 0; iLink < MaxNumberOfLinks; ++iLink) {
968 const auto& linkInfo = linkInfoArray[iLink];
969 if (!linkInfo.IsPresent) {
970 continue;
971 }
972 printf("Event %4d, Link %2d\n", eventNumber, iLink);
973
974 const int ep = iLink >= 12;
975 const int link = iLink - (ep)*12;
976 auto outputFileName = fmt::format("{}/CRU_{:02}_EP_{}_Link_{:02}", outputDirectory.data(), (int)mCRU, ep, link);
977 std::ofstream outputFile(outputFileName, std::ios_base::binary | std::ios_base::app);
978
979 for (auto packetNumber : linkInfo.PacketPositions) {
980 const auto& packet = mPacketDescriptorMaps[iLink][packetNumber];
981 file.seekg(packet.getPayloadOffset(), file.beg);
982 file.read(buffer, packet.getPayloadSize());
983 outputFile.write(buffer, packet.getPayloadSize());
984 }
985 }
986 }
987}
988//==============================================================================
989//===| stream overloads for helper classes |====================================
990//
991
992void ADCRawData::streamTo(std::ostream& output) const
993{
994 const auto numTimeBins = std::min(getNumTimebins(), mNumTimeBins);
995 for (int i = 0; i < numTimeBins * 16; i++) {
996 if (i % 16 == 0) {
997 output << std::setw(4) << std::to_string(i / 16) << " : ";
998 }
999 output << std::setw(4) << mADCRaw[(mOutputStream)][i];
1000 output << (((i + 1) % 16 == 0) ? "\n" : " ");
1001 };
1002};
1003
1004void GBTFrame::streamFrom(std::istream& input)
1005{
1006 mFilePos = input.tellg();
1007 mFrameNum++;
1008 for (int i = 0; i < 4; i++) {
1009 input.read(reinterpret_cast<char*>(&(mData[i])), sizeof(mData[i]));
1010 };
1011}
1012
1013// std::ostream& operator<<(std::ostream& output, const RawReaderCRU::GBTFrame& frame)
1014void GBTFrame::streamTo(std::ostream& output) const
1015{
1016 const auto offset = mPrevHWpos ^ 4;
1017 output << std::dec << "\033[94m"
1018 << std::setfill('0') << std::setw(8) << mPacketNum << " "
1019 << std::setfill('0') << std::setw(8) << mFilePos << " "
1020 << std::setfill('0') << std::setw(8) << mFrameNum << " : "
1021 << std::hex
1022 << std::setfill('0') << std::setw(8) << mData[3] << "."
1023 << std::setfill('0') << std::setw(8) << mData[2] << "."
1024 << std::setfill('0') << std::setw(8) << mData[1] << "."
1025 << std::setfill('0') << std::setw(8) << mData[0] << " : "
1026 << "\033[0m";
1027 for (int i = 0; i < 5; i++) {
1028 for (int j = 0; j < 4; j++) {
1029 output << std::hex << std::setw(4) << mFrameHalfWords[i][j + offset] << " ";
1030 }
1031 output << "| ";
1032 };
1033 output << std::endl;
1034}
1035
1036// friend std::ostream& operator<<(std::ostream& output, const PacketDescriptor& PD)
1038{
1039 output << "===| Packet Descriptor |================================="
1040 << "\n";
1041 output << "CRU ID : " << std::dec << getCRUID() << "\n";
1042 output << "Link ID : " << getLinkID() << "\n";
1043 output << "Global Link ID : " << getGlobalLinkID() << "\n";
1044 output << "Header Offset : " << getHeaderOffset() << "\n";
1045 output << "Payload Offset : " << getPayloadOffset() << "\n";
1046 output << "Payload Size : " << getPayloadSize() << "\n";
1047 output << "\n";
1048}
1049
1050std::ostream& operator<<(std::ostream& output, const RDH& rdh)
1051{
1052 output << "word0 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word0 << "\n"
1053 << std::dec;
1054 output << " version : " << RDHUtils::getVersion(rdh) << "\n";
1055 output << " headerSize : " << RDHUtils::getHeaderSize(rdh) << "\n";
1056 if (RDHUtils::getVersion(rdh) == 4) {
1057 output << " blockLength : " << RDHUtils::getBlockLength(rdh) << "\n";
1058 }
1059 output << " feeId : " << RDHUtils::getFEEID(rdh) << "\n";
1060 output << " priority : " << RDHUtils::getPriorityBit(rdh) << "\n";
1061 output << "\n";
1062
1063 output << "word1 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word1 << "\n"
1064 << std::dec;
1065 output << " Offset to next : " << int(RDHUtils::getOffsetToNext(rdh)) << "\n";
1066 output << " Memory size : " << int(RDHUtils::getMemorySize(rdh)) << "\n";
1067 output << " LinkID : " << int(RDHUtils::getLinkID(rdh)) << "\n";
1068 output << " Global LinkID : " << int(RDHUtils::getLinkID(rdh)) + (((rdh.word1 >> 32) >> 28) * 12) << "\n";
1069 output << " CRUid : " << RDHUtils::getCRUID(rdh) << "\n";
1070 output << " Packet Counter : " << RDHUtils::getPacketCounter(rdh) << "\n";
1071 output << " DataWrapper-ID : " << (((rdh.word1 >> 32) >> 28) & 0x01) << "\n";
1072 output << "\n";
1073
1074 output << "word2 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word2 << "\n"
1075 << std::dec;
1076 output << " triggerOrbit : " << RDHUtils::getTriggerOrbit(rdh) << "\n";
1077 output << " heartbeatOrbit : " << RDHUtils::getHeartBeatOrbit(rdh) << "\n";
1078 output << "\n";
1079
1080 output << "word3 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word3 << "\n"
1081 << std::dec;
1082 output << "\n";
1083
1084 output << "word4 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word4 << "\n"
1085 << std::dec;
1086 output << " triggerBC : " << RDHUtils::getTriggerBC(rdh) << "\n";
1087 output << " heartbeatBC : " << RDHUtils::getHeartBeatBC(rdh) << "\n";
1088 output << " triggerType : " << RDHUtils::getTriggerType(rdh) << "\n";
1089 output << "\n";
1090
1091 output << "word5 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word5 << "\n"
1092 << std::dec;
1093 output << "\n";
1094
1095 output << "word6 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word6 << "\n"
1096 << std::dec;
1097 output << " detectorField : " << RDHUtils::getDetectorField(rdh) << "\n";
1098 output << " par : " << RDHUtils::getDetectorPAR(rdh) << "\n";
1099 output << " stop : " << RDHUtils::getStop(rdh) << "\n";
1100 output << " pageCnt : " << RDHUtils::getPageCounter(rdh) << "\n";
1101 output << "\n";
1102
1103 output << "word7 : 0x" << std::setfill('0') << std::setw(16) << std::hex << rdh.word7 << "\n"
1104 << std::dec;
1105 output << "\n";
1106
1107 return output;
1108}
1109
1111{
1112 fmt::print("{:>5} {:>4} {:>4} {:>4} {:>6} {:>4} {:>3} {:>4} {:>10} {:>5} {:>1}\n",
1113 "PkC", "pCnt", "trg", "fId", "OffN", "Mem", "CRU", "GLID", "HBOrbit", "HBBC", "s");
1114}
1115
1116void printHorizontal(const RDH& rdh)
1117{
1118 const int globalLinkID = int(RDHUtils::getLinkID(rdh)) + (((rdh.word1 >> 32) >> 28) * 12);
1119
1120 fmt::print("{:>5} {:>4} {:>4} {:>4} {:>6} {:>4} {:>3} {:>4} {:>10} {:>5} {:>1}\n",
1121 (uint64_t)RDHUtils::getPacketCounter(rdh),
1122 (uint64_t)RDHUtils::getPageCounter(rdh),
1123 (uint64_t)RDHUtils::getTriggerType(rdh),
1124 (uint64_t)RDHUtils::getFEEID(rdh),
1125 (uint64_t)RDHUtils::getOffsetToNext(rdh),
1126 (uint64_t)RDHUtils::getMemorySize(rdh),
1127 (uint64_t)RDHUtils::getCRUID(rdh),
1128 (uint64_t)globalLinkID,
1129 (uint64_t)RDHUtils::getHeartBeatOrbit(rdh),
1130 (uint64_t)RDHUtils::getHeartBeatBC(rdh),
1131 (uint64_t)RDHUtils::getStop(rdh));
1132}
1133
1134template <typename DataType>
1135std::istream& operator>>(std::istream& input, DataType& data)
1136{
1137 const int dataTypeSize = sizeof(data);
1138 auto charPtr = reinterpret_cast<char*>(&data);
1139 input.read(charPtr, dataTypeSize);
1140 return input;
1141}
1142
1143//==============================================================================
1145{
1146 if (mIsInitialized) {
1147 return;
1148 }
1149
1150 for (auto& reader : mRawReadersCRU) {
1151 reader->scanFile();
1152 }
1153
1154 mEventSync.sortEvents();
1155 mEventSync.analyse(mRawDataType);
1156
1157 O2INFO("Event information:");
1158 O2INFO(" Number of all events: %lu", getNumberOfEvents());
1159 O2INFO(" Number of complete events: %lu", getNumberOfCompleteEvents());
1160
1161 if (CHECK_BIT(mDebugLevel, DebugLevel::EventInfo)) {
1162 std::cout << mEventSync;
1163 }
1164
1165 mIsInitialized = true;
1166}
1167
1168void RawReaderCRUManager::setupReaders(const std::string_view inputFileNames,
1169 uint32_t numTimeBins,
1170 uint32_t debugLevel,
1171 uint32_t verbosity,
1172 const std::string_view outputFilePrefix)
1173{
1174 reset();
1175 const TString files = gSystem->GetFromPipe(TString::Format("ls %s", inputFileNames.data()));
1176 std::unique_ptr<TObjArray> arr(files.Tokenize("\n"));
1177 setDebugLevel(debugLevel);
1178
1179 for (auto file : *arr) {
1180 // fix the number of time bins
1181 auto& reader = createReader(file->GetName(), numTimeBins);
1182 reader.setReaderNumber(mRawReadersCRU.size() - 1);
1183 reader.setVerbosity(verbosity);
1184 reader.setDebugLevel(debugLevel);
1185 reader.setOutputFilePrefix(outputFilePrefix);
1186 O2INFO("Adding file: %s\n", file->GetName());
1187 }
1188}
1189
1190void RawReaderCRUManager::copyEvents(const std::vector<uint32_t> eventNumbers, std::string_view outputDirectory, std::ios_base::openmode mode)
1191{
1192 // make sure events have been built
1193 init();
1194 const auto& cruSeen = mEventSync.getCRUSeen();
1195
1196 for (size_t iCRU = 0; iCRU < cruSeen.size(); ++iCRU) {
1197 const auto readerNumber = cruSeen[iCRU];
1198 if (readerNumber >= 0) {
1199 auto& reader = mRawReadersCRU[readerNumber];
1200 reader->forceCRU(iCRU);
1201 reader->copyEvents(eventNumbers, outputDirectory.data(), mode);
1202 }
1203 }
1204}
1205
1206void RawReaderCRUManager::copyEvents(const std::string_view inputFileNames, const std::vector<uint32_t> eventNumbers, std::string_view outputDirectory, std::ios_base::openmode mode)
1207{
1208 RawReaderCRUManager manager;
1209 manager.setupReaders(inputFileNames);
1210 manager.copyEvents(eventNumbers, outputDirectory, mode);
1211}
1212
1213void RawReaderCRUManager::writeGBTDataPerLink(std::string_view outputDirectory, int maxEvents)
1214{
1215 init();
1216 const auto& cruSeen = mEventSync.getCRUSeen();
1217
1218 for (size_t iCRU = 0; iCRU < cruSeen.size(); ++iCRU) {
1219 const auto readerNumber = cruSeen[iCRU];
1220 if (readerNumber >= 0) {
1221 auto& reader = mRawReadersCRU[readerNumber];
1222 reader->forceCRU(iCRU);
1223 reader->writeGBTDataPerLink(outputDirectory, maxEvents);
1224 }
1225 }
1226}
1227
1228void RawReaderCRUManager::writeGBTDataPerLink(const std::string_view inputFileNames, std::string_view outputDirectory, int maxEvents)
1229{
1230 if (!std::filesystem::exists(outputDirectory)) {
1231 if (!std::filesystem::create_directories(outputDirectory)) {
1232 LOG(fatal) << "could not create output directory " << outputDirectory;
1233 } else {
1234 LOG(info) << "created output directory " << outputDirectory;
1235 }
1236 }
1237
1238 RawReaderCRUManager manager;
1239 manager.setupReaders(inputFileNames);
1240 manager.writeGBTDataPerLink(outputDirectory, maxEvents);
1241}
1242
1243void RawReaderCRUManager::processEvent(uint32_t eventNumber, EndReaderCallback endReader)
1244{
1245 const auto& cruSeen = mEventSync.getCRUSeen();
1246
1247 for (size_t iCRU = 0; iCRU < cruSeen.size(); ++iCRU) {
1248 const auto readerNumber = cruSeen[iCRU];
1249 if (readerNumber >= 0) {
1250 auto& reader = mRawReadersCRU[readerNumber];
1251 if (reader->getFillADCdataMap()) {
1252 LOGF(warning, "Filling of ADC data map not supported in RawReaderCRUManager::processEvent, it is disabled now. use ADCDataCallback");
1253 reader->setFillADCdataMap(false);
1254 }
1255 reader->setEventNumber(eventNumber);
1256 reader->forceCRU(iCRU);
1257 reader->processLinks();
1258 if (endReader) {
1259 endReader();
1260 }
1261 }
1262 }
1263}
#define verbosity
uint64_t orbit
Definition RawEventData.h:6
int32_t i
#define O2ERROR(...)
Definition Logger.h:18
#define O2INFO(...)
Definition Logger.h:17
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
uint32_t j
Definition RawData.h:0
void printHorizontal(const RDH &rdh)
std::istream & operator>>(std::istream &input, DataType &data)
void printHeader()
#define CHECK_BIT(var, pos)
std::ostream & operator<<(std::ostream &output, const RDH &rdh)
unsigned short number() const
Definition CRU.h:60
unsigned char partition() const
Definition CRU.h:63
static Mapper & instance(const std::string mappingDir="")
Definition Mapper.h:44
helper to store the ADC raw data
void streamTo(std::ostream &output) const
write data to ostream
const DataVector & getDataVector(int stream) const
get the data vector for a specific stream
uint32_t getNumTimebins() const
number of time bin for selected stream
void setNumTimebins(uint32_t numTB)
void setOutputStream(uint32_t outStream)
select the stream for which data should be processed
helper to encapsulate a GBTFrame
void setFrameNumber(uint32_t frameNumber)
set packet number
void readFromMemory(gsl::span< const std::byte > data)
read from memory
void updateSyncCheck(SyncArray &syncArray)
update the sync check
void getFrameHalfWords()
extract the 4 5b halfwords for the 5 data streams from one GBT frame
void getAdcValues(ADCRawData &rawData)
decode ADC values
void setPacketNumber(uint32_t packetNumber)
set packet number
void streamFrom(std::istream &input)
read from istream
void streamTo(std::ostream &output) const
write data to ostream
const SyncArray & getSyncArray() const
get the syncronisation array
const EventInfo & getEventInfo(size_t eventNumber) const
information of specific event
LinkInfo & getLinkInfo(const RDH &rdh, DataType dataType)
get link information for a specific event and cru
static constexpr size_t ExpectedNumberOfPacketsPerHBFrame
expected number of packets in one HB frame per link for HB scaling mode
size_t getNumberOfEvents() const
get number of all events
void setCRUSeen(const CRU cru, const uint16_t reader=0)
set a cru as seen
void analyse(RAWDataType rawDataType=RAWDataType::GBT)
analyse events and mark complete events
void setLinksSeen(const CRU cru, const std::bitset< MaxNumberOfLinks > &links)
set links that were seen for a CRU
const LinkInfoArray_t & getLinkInfoArrayForEvent(size_t eventNumber, int cru) const
get array with all link informaiton for a specific event number and cru
void streamTo(std::ostream &output) const
write data to ostream
EventInfo & createEvent(const uint32_t heartbeatOrbit, DataType dataType)
create a new event or return the one with the given HB orbit
void writeGBTDataPerLink(std::string_view outputDirectory, int maxEvents=-1)
copy single events from raw input files to another file
DataType getDataType() const
get data type
void init()
initialize all readers
std::function< void()> EndReaderCallback
void processEvent(uint32_t eventNumber, EndReaderCallback endReader=nullptr)
process event calling mADCDataCallback to process values
void setupReaders(const std::string_view inputFileNames, uint32_t numTimeBins=1000, uint32_t debugLevel=0, uint32_t verbosity=0, const std::string_view outputFilePrefix="")
size_t getNumberOfEvents() const
get number of all events
void copyEvents(const std::vector< uint32_t > eventNumbers, std::string_view outputDirectory, std::ios_base::openmode mode=std::ios_base::openmode(0))
copy single events to another file
void setDebugLevel(uint32_t debugLevel)
set debug level
auto & getReaders()
return vector of readers
void streamTo(std::ostream &output) const
write data to ostream
int processDataFile()
process all data for the selected link reading single 8k packet from file
void findSyncPositions()
find sync positions for all links
static void processFile(const std::string_view inputFile, uint32_t timeBins=0, uint32_t linkMask=0, uint32_t stream=0, uint32_t debugLevel=0, uint32_t verbosity=0, const std::string_view outputFilePrefix="")
bool checkLinkPresent(uint32_t link)
status bits of present links
void processLinks(const uint32_t linkMask=0)
process links
void processDataMemory()
Collect data to memory and process data.
void runADCDataCallback(const ADCRawData &rawData)
run a data filling callback function
std::ifstream & getFileHandle()
file handling
void setDebugLevel(uint32_t debugLevel=1)
set debug level
void copyEvents(const std::vector< uint32_t > &eventNumbers, std::string outputDirectory, std::ios_base::openmode mode=std::ios_base::openmode(0))
copy single events to another file
int processMemory(const std::vector< std::byte > &data, ADCRawData &rawData)
size_t getNumberOfEvents() const
get number of events
bool syncFoundForLink(const int link) const
void processLinkZS()
process data in case of Link based zero suppression
int processPacket(GBTFrame &gFrame, uint32_t startPos, uint32_t size, ADCRawData &rawData)
process single packet
void writeGBTDataPerLink(std::string_view outputDirectory, int maxEvents=-1)
write GBT data into separate files per link
struct _cl_event * event
Definition glcorearb.h:2982
GLenum mode
Definition glcorearb.h:266
GLuint buffer
Definition glcorearb.h:655
GLsizeiptr size
Definition glcorearb.h:659
GLenum GLsizei dataSize
Definition glcorearb.h:3994
GLboolean * data
Definition glcorearb.h:298
GLintptr offset
Definition glcorearb.h:660
GLuint GLuint stream
Definition glcorearb.h:1806
GLfloat green
Definition glcorearb.h:279
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
std::istream & operator>>(std::istream &stream, o2::header::RAWDataHeaderV4 &header)
bool processZSdata(const char *data, size_t size, rdh_utils::FEEIDType feeId, uint32_t orbit, uint32_t referenceOrbit, uint32_t syncOffsetReference, ADCCallback fillADC)
@ HBScaling
heart beat sclaing mode
@ GBT
GBT encoded raw data.
@ LinkZS
Link based zero suppression.
@ RDHDump
dump full RDH
@ EventInfo
dump event synchronisation information
@ GBTFrames
dump all GBT frames
@ SyncPositions
dump sync positions
@ ADCValues
dump extracted ADC values
uint16_t FEEIDType
Definition RDHUtils.h:26
std::ostream & operator<<(std::ostream &out, const tpc::FECInfo &fec)
Definition FECInfo.cxx:27
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
the main header struct
Definition DataHeader.h:618
TForbitType firstTForbit
Definition DataHeader.h:674
PayloadSizeType payloadSize
Definition DataHeader.h:666
static void setFEEID(RDHv4 &rdh, uint16_t v)
Definition RDHUtils.h:146
static constexpr int getVersion()
get numeric version of the RDH
Definition RDHUtils.h:58
static void setCRUID(H &rdh, uint16_t v, NOTPTR(H))
Definition RDHUtils.h:288
CRUInfoArray_t CRUInfoArray
Link information for each cru.
bool hasHearbeatOrbit(uint32_t heartbeatOrbit) const
check if heartbeatOrbit contributes to the event
std::vector< uint32_t > HeartbeatOrbits
vector of heartbeat orbits contributing to the event
helper to store link information in an event
bool IsPresent
if the link is present in the current event
size_t PayloadSize
total payload size of the link in the present event
std::vector< size_t > PacketPositions
all packet positions of this link in an event
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
vec clear()