CbmRoot
Loading...
Searching...
No Matches
CbmRichUnpackAlgo.cxx
Go to the documentation of this file.
1/* Copyright (C) 2021 Goethe-University Frankfurt, Frankfurt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pascal Raisig [committer] */
4
5#include "CbmRichUnpackAlgo.h"
6
7#include <FairParGenericSet.h>
8#include <FairTask.h>
9#include <Logger.h>
10
11#include <Rtypes.h>
12#include <RtypesCore.h>
13
14#include <cstdint>
15
16
18
20
21// ---- unpack
22bool CbmRichUnpackAlgo::unpack(const fles::Timeslice* ts, std::uint16_t icomp, UInt_t imslice)
23{
24
25 const fles::MicrosliceView mv = ts->get_microslice(icomp, imslice);
26 const fles::MicrosliceDescriptor& msDesc = mv.desc();
27
29 reader.SetData(mv.content(), msDesc.size);
30
31 auto mstime = msDesc.idx;
32 fMsRefTime = mstime - fTsStartTime;
33
34 // There are a lot of MS with 8 bytes size
35 // Does one need these MS?
36 if (reader.GetSize() <= 8) return true;
37
38 while (true) {
39 processTrbPacket(reader);
40 // -4*2 for 2 last words which contain microslice index
41 if (reader.GetOffset() >= reader.GetSize() - 8) break;
42 // -4*3 for 0xffffffff padding and 2 last words which contain microslice index
43 if (reader.IsNextPadding() && reader.GetOffset() >= reader.GetSize() - 12) break;
44 }
45
46 uint32_t msIndexWord1 = reader.NextWord();
47 if (isLog()) LOG(debug4) << getLogHeader(reader) << "Microslice Index 1:" << reader.GetWordAsHexString(msIndexWord1);
48
49 uint32_t msIndexWord2 = reader.NextWord();
50 if (isLog()) LOG(debug4) << getLogHeader(reader) << "Microslice Index 2:" << reader.GetWordAsHexString(msIndexWord2);
51
52
53 return true;
54}
55
57{
58 uint32_t word = reader.NextWord();
59 uint32_t hubId = word & 0xffff; // 16 bits
60 uint32_t hubSize = (word >> 16) & 0xffff; // 16 bits
61 if (isLog())
62 LOG(debug4) << getLogHeader(reader) << "hubId:0x" << std::hex << hubId << std::dec << " hubSize:" << hubSize;
63
64 //if ((HubId == 0xc001) || (HubId == 0xc000)) //CTS subevent?
65 //if (HubId == 0x5555)
66 //if (((HubId >> 8) & 0xff) == 0x82) // TRB subevent? // TODO: check if it starts from 0x82
67
68 // if true then it is CTS sub-sub-event
69 bool isLast = false;
70 size_t counter = 0;
71 size_t totalSize = 0;
72 while (!isLast) {
73 word = reader.NextWord();
74 uint32_t subSubEventId = word & 0xffff; // 16 bits
75 uint32_t subSubEventSize = (word >> 16) & 0xffff; // 16 bits
76 isLast = reader.IsLastSubSubEvent(subSubEventSize); // if true then it is CTS sub-sub-event
77 counter++;
78 totalSize += (1 + subSubEventSize);
79
80 if (isLog())
81 LOG(debug4) << getLogHeader(reader) << counter << ((isLast) ? " CTS" : " DiRICH") << " subSubEventId:0x"
82 << std::hex << subSubEventId << std::dec << " subSubEventSize:" << subSubEventSize;
83
84 if (!isLast) { // DiRICH event
85 // check correctness of subsub event, for safety reasons
86 if (((subSubEventId >> 12) & 0xF) != 0x7) {
87 LOG(error) << getLogHeader(reader) << "ERROR: subSubEventId has strange value:0x" << std::hex << subSubEventId
88 << std::dec;
89 }
90 processSubSubEvent(reader, subSubEventSize, subSubEventId);
91 }
92 else { // CTS event
93 processCtsSubSubEvent(reader, subSubEventSize, subSubEventId);
94 }
95
96 // if (fbDebugMonitorMode) {
97 // //This address calculation is just for mCBM; will be a problem when using full CBM RICH acceptance
98 // uint16_t histAddr = ((subSubEventId >> 8) & 0xF) * 18 + ((subSubEventId >> 4) & 0xF) * 2 + (subSubEventId & 0xF);
99 // fhSubSubEventSize->Fill(histAddr, subSubEventSize); // Words in a DiRICH
100 // }
101
102 if ((totalSize == hubSize && !isLast) || (totalSize != hubSize && isLast)) {
103 if (isLog()) LOG(error) << "ERROR: totalSize OR isLast is wrong";
104 }
105
106 if (totalSize >= hubSize || isLast) break;
107 }
108
109 // read last words
110 int lastWordsCounter = 0;
111 while (true) {
112 lastWordsCounter++;
113 word = reader.NextWord();
114 if (isLog()) LOG(debug4) << getLogHeader(reader);
115 if (word == 0x600dda7a) break;
116 if (lastWordsCounter >= 7) {
117 LOG(error) << getLogHeader(reader)
118 << "CbmMcbm2018UnpackerAlgoRich::ProcessHubBlock() ERROR: No word == 0x600dda7a";
119 }
120 }
121}
122
124 uint32_t subSubEventId)
125{
126 uint32_t word = reader.NextWord();
127 uint32_t ctsState = word & 0xffff; // 16 bits
128 uint32_t nofInputs = (word >> 16) & 0xf; // 4 bits
129 uint32_t nofTrigCh = (word >> 20) & 0x1f; // 5 bits
130 uint32_t inclLastIdle = (word >> 25) & 0x1; // 1 bit
131 uint32_t inclTrigInfo = (word >> 26) & 0x1; // 1 bit
132 uint32_t inclTime = (word >> 27) & 0x1; // 1 bit
133 uint32_t ETM = (word >> 28) & 0x3; // 2 bits
134 uint32_t ctsInfoSize = 2 * nofInputs + 2 * nofTrigCh + 2 * inclLastIdle + 3 * inclTrigInfo + inclTime; // in words
135 switch (ETM) {
136 case 0: break;
137 case 1: ctsInfoSize += 1; break;
138 case 2: ctsInfoSize += 4; break;
139 case 3: break;
140 }
141 if (isLog()) LOG(debug4) << getLogHeader(reader) << "CTS ctsState:" << ctsState << " ctsInfoSize:" << ctsInfoSize;
142 for (uint32_t i = 0; i < ctsInfoSize; i++) {
143 word = reader.NextWord(); // do nothing?
144 if (isLog()) LOG(debug4) << getLogHeader(reader) << "CTS info words";
145 }
146 int nofTimeWords = subSubEventSize - ctsInfoSize - 1;
147 processSubSubEvent(reader, nofTimeWords, subSubEventId);
148}
149
151 uint32_t subSubEventId)
152{
153 // Store if a certain TDC word type was analysed,
154 // later one can check if the order is correct
155 bool wasHeader = false;
156 bool wasEpoch = false;
157 bool wasTime = false;
158 bool wasTrailer = false;
159 uint32_t epoch = 0; // store last epoch obtained in sub-sub-event
160 bool errorInData = false;
161
162 // Store last raising edge time for each channel or -1. if no time
163 // this array is used to match raising and falling edges
164 std::vector<double> raisingTime(33, -1.);
165
166 // check if DiRICH (SubSubEvId) is masked
167 bool DiRICH_masked = false;
168 if (checkMaskedDiRICH(subSubEventId)) {
169 DiRICH_masked = true;
170 }
171
172 for (int i = 0; i < nofTimeWords; i++) {
173 uint32_t word = reader.NextWord();
174 if (DiRICH_masked) continue;
176
178 if (!wasHeader || !wasEpoch || wasTrailer) {
179 LOG(error) << getLogHeader(reader) << "DiRICH 0x" << std::hex << subSubEventId << std::dec
180 << ": illegal position of TDC Time (before header/epoch or after trailer)";
181 errorInData = true;
182 continue;
183 }
184 wasTime = true;
185 processTimeDataWord(reader, i, epoch, word, subSubEventId, raisingTime);
186 }
187 else if (type == CbmRichUnpackAlgoTdcWordType::Epoch) {
188 if (!wasHeader || wasTrailer) {
189 LOG(error) << getLogHeader(reader) << "DiRICH 0x" << std::hex << subSubEventId << std::dec
190 << ": illegal position of TDC Epoch (before header or after trailer)";
191 errorInData = true;
192 continue;
193 }
194 wasEpoch = true;
196 if (isLog()) LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << i << "] epoch:" << epoch;
197 }
198 else if (type == CbmRichUnpackAlgoTdcWordType::Header) {
199 if (wasEpoch || wasTime || wasTrailer) {
200 LOG(error) << getLogHeader(reader) << "DiRICH 0x" << std::hex << subSubEventId << std::dec
201 << ": illegal position of TDC Header (after time/epoch/trailer)";
202 errorInData = true;
203 continue;
204 }
205 wasHeader = true;
206 // uint16_t errorBits = CbmRichUnpackAlgoTdcWordReader::ProcessHeader(word);
207 // ErrorMsg(errorBits, CbmRichUnpackAlgoErrorType::tdcHeader, subSubEventId);
208 if (isLog()) LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << i << "] header";
209 }
210 else if (type == CbmRichUnpackAlgoTdcWordType::Trailer) {
211 if (!wasEpoch || !wasTime || !wasHeader) {
212 LOG(error) << getLogHeader(reader) << "DiRICH 0x" << std::hex << subSubEventId << std::dec
213 << ": illegal position of TDC Trailer (before time/epoch/header)";
214 errorInData = true;
215 continue;
216 }
217 wasTrailer = true;
218 // uint16_t errorBits = CbmRichUnpackAlgoTdcWordReader::ProcessTrailer(word);
219 // ErrorMsg(errorBits, CbmRichUnpackAlgoErrorType::tdcTrailer, subSubEventId);
220 if (isLog()) LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << i << "] trailer";
221 }
222 else if (type == CbmRichUnpackAlgoTdcWordType::Debug) {
223 // for the moment do nothing
224 }
225 else if (type == CbmRichUnpackAlgoTdcWordType::Error) {
226 LOG(error) << getLogHeader(reader) << "Wrong TDC word!!! marker:" << ((word >> 29) & 0x7);
227 errorInData = true;
228 }
229 }
230
231 if (errorInData) {
232 //TODO:
233 }
234}
235
237{
238 processMbs(reader, false); // Current MBS
239 processMbs(reader, true); // Previous MBS
240
241 uint32_t trbNum = reader.NextWord(); // TRB trigger number
242 if (isLog()) LOG(debug4) << getLogHeader(reader) << "TRB Num:" << reader.GetWordAsHexString(trbNum);
243
244 processHubBlock(reader);
245}
246
248{
249 uint32_t word = reader.NextWord();
250 uint32_t mbsNum = word & 0xffffff; //24 bits
251 uint32_t nofCtsCh = (word >> 24) & 0xff; // 8 bits
252 if (isLog())
253 LOG(debug4) << getLogHeader(reader) << "MBS mbsNum:0x" << std::hex << mbsNum << std::dec
254 << " nofCtsCh:" << nofCtsCh;
255
256 for (uint32_t i = 0; i < nofCtsCh; i++) {
257 uint32_t wordEpoch = reader.NextWord();
258 uint32_t epoch = CbmRichUnpackAlgoTdcWordReader::ProcessEpoch(wordEpoch);
259 if (isLog()) LOG(debug4) << getLogHeader(reader) << "MBS ch:" << i << " epoch:" << epoch;
260
261 uint32_t wordTime = reader.NextWord();
264 if (isLog()) LOG(debug4) << getLogHeader(reader) << "MBS ch:" << i << " " << td.ToString();
265
266 double fullTime = calculateTime(epoch, td.fCoarse, td.fFine);
267
268 if (isPrev && td.fChannel == 0) fMbsPrevTimeCh0 = fullTime;
269 if (isPrev && td.fChannel == 1) fMbsPrevTimeCh1 = fullTime;
270 }
271
272 double mbsCorr = fMbsPrevTimeCh1 - fMbsPrevTimeCh0;
273 if (isLog())
274 LOG(debug4) << getLogHeader(reader) << "MBS Prev ch1:" << std::setprecision(15) << fMbsPrevTimeCh1
275 << " ch0:" << fMbsPrevTimeCh0 << " corr:" << mbsCorr;
276}
277
278// ---- processTimeDataWord ----
280 uint32_t tdcWord, uint32_t subSubEventId, std::vector<double>& raisingTime)
281{
284 double fullTime = calculateTime(epoch, td.fCoarse, td.fFine);
285
286
287 if (td.fChannel == 0) {
288 if (td.IsRisingEdge()) {
289 fPrevLastCh0ReTime[subSubEventId] = fLastCh0ReTime[subSubEventId];
290 fLastCh0ReTime[subSubEventId] = fullTime;
291 if (isLog())
292 LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << iTdc << "] " << td.ToString()
293 << " CH0 Last:" << std::setprecision(15) << fLastCh0ReTime[subSubEventId]
294 << " PrevLast:" << fPrevLastCh0ReTime[subSubEventId]
295 << " diff:" << fLastCh0ReTime[subSubEventId] - fPrevLastCh0ReTime[subSubEventId];
296 }
297 }
298 else {
299 double dT = fullTime - fPrevLastCh0ReTime[subSubEventId];
300 double mbsCorr = fMbsPrevTimeCh1 - fMbsPrevTimeCh0;
301 double fullTimeCorr = dT - mbsCorr;
302 if (isLog())
303 LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << iTdc << "] " << td.ToString()
304 << " time:" << std::setprecision(15) << fullTime << " fullTimeCorr:" << fullTimeCorr;
305
306 if (td.fChannel < 1 || td.fChannel >= raisingTime.size()) {
307 LOG(error) << "ERROR: channel number is out of limit. Channel:" << td.fChannel;
308 }
309
310 if (td.IsRisingEdge()) {
311 // always store the latest raising edge. It means that in case RRFF situation only middle RF will be matched.
312 raisingTime[td.fChannel] = fullTimeCorr;
313 }
314 else {
315 if (raisingTime[td.fChannel] == -1.) {
316 //No raising channel was found before. Skip this falling edge time.
317 if (isLog())
318 LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << iTdc << "] "
319 << "No raising channel was found before. Skip this falling edge time.";
320 }
321 else {
322 // Matching was found, calculate ToT, if tot is in a good range -> create digi
323 double ToT = fullTimeCorr - raisingTime[td.fChannel];
324 if (isLog())
325 LOG(debug4) << getLogHeader(reader) << "SubSubEv[" << iTdc << "] "
326 << "ToT:" << ToT;
327 if (ToT >= fToTMin && ToT <= fToTMax) {
328 // if (fbMonitorMode) {
329 // TH1D* h = GetTotH1(subSubEventId, td.fChannel);
330 // if (h != nullptr) h->Fill(ToT);
331
332 // TH2D* h2 = GetTotH2(subSubEventId);
333 // if (h2 != nullptr) h2->Fill(td.fChannel, ToT);
334 // }
335 if (fullTimeCorr >= 0.0) {
336 writeOutputDigi(subSubEventId, td.fChannel, raisingTime[td.fChannel], ToT);
337 }
338 }
339 // pair was created, set raising edge to -1.
340 raisingTime[td.fChannel] = -1.;
341 }
342 }
343 }
344}
345
346void CbmRichUnpackAlgo::writeOutputDigi(Int_t fpgaID, Int_t channel, Double_t time, Double_t tot)
347{
348 Double_t ToTcorr = fbDoToTCorr ? fUnpackPar.GetToTshift(fpgaID, channel) : 0.;
349 Int_t pixelUID = this->getPixelUID(fpgaID, channel);
350 //check ordering
351 Double_t finalTime = time + (Double_t) fMsRefTime - fSystemTimeOffset;
352 // Double_t finalTime = time - fSystemTimeOffset - fTsStartTime;
353
354 fOutputVec.emplace_back(pixelUID, finalTime, tot - ToTcorr);
355
356 // Sorting of Digis (Now done in general part)
357 /*
358 Double_t lastTime = 0.;
359
360 if (fOutputVec.size() < 1) { fOutputVec.emplace_back(pixelUID, finalTime, tot - ToTcorr); }
361 else {
362 lastTime = fOutputVec[fOutputVec.size() - 1].GetTime();
363 if (fOutputVec[0].GetTime() > finalTime) {
364 fOutputVec.emplace(fOutputVec.begin(), pixelUID, finalTime, tot - ToTcorr);
365 }
366 else if (lastTime > finalTime) {
367 for (int i = fOutputVec.size() - 1; i >= 0; i--) {
368 lastTime = fOutputVec[i].GetTime();
369 if (lastTime <= finalTime) {
370 fOutputVec.emplace(fOutputVec.begin() + i + 1, pixelUID, finalTime, tot - ToTcorr);
371 break;
372 }
373 }
374 }
375 else {
376 fOutputVec.emplace_back(pixelUID, finalTime, tot - ToTcorr);
377 }
378 }
379 */
380}
381
382
ClassImp(CbmConverterManager)
CbmRichUnpackAlgoTdcWordType
Baseclass for the TrdR unpacker algorithms.
Double_t GetToTshift(Int_t tdc, Int_t ch) const
bool checkMaskedDiRICH(Int_t subSubEventId)
Unpack a given microslice. To be implemented in the derived unpacker algos.
CbmMcbm2018RichPar fUnpackPar
Parameters for the unpacking.
double calculateTime(uint32_t epoch, uint32_t coarse, uint32_t fine)
Int_t getPixelUID(Int_t fpgaID, Int_t ch) const
std::string getLogHeader(CbmRichUnpackAlgoMicrosliceReader &reader)
std::string GetWordAsHexString(uint32_t word)
bool IsLastSubSubEvent(uint32_t subSubEventSize)
void SetData(const uint8_t *data, size_t size)
static CbmRichUnpackAlgoTdcWordType GetTdcWordType(uint32_t tdcWord)
static void ProcessTimeData(uint32_t tdcWord, CbmRichUnpackAlgoTdcTimeData &outData)
static uint32_t ProcessEpoch(uint32_t tdcWord)
void processMbs(CbmRichUnpackAlgoMicrosliceReader &reader, bool isPrev)
CbmRichUnpackAlgo()
Create the Cbm Trd Unpack AlgoBase object.
virtual ~CbmRichUnpackAlgo()
Destroy the Cbm Trd Unpack Task object.
void processSubSubEvent(CbmRichUnpackAlgoMicrosliceReader &reader, int nofTimeWords, uint32_t subSubEventId)
std::map< uint32_t, double > fPrevLastCh0ReTime
void processTrbPacket(CbmRichUnpackAlgoMicrosliceReader &reader)
void processTimeDataWord(CbmRichUnpackAlgoMicrosliceReader &reader, int iTdc, uint32_t epoch, uint32_t tdcWord, uint32_t subSubEventId, std::vector< double > &raisingTime)
void processCtsSubSubEvent(CbmRichUnpackAlgoMicrosliceReader &reader, uint32_t subSubEventSize, uint32_t subSubEventId)
bool unpack(const fles::Timeslice *ts, std::uint16_t icomp, UInt_t imslice)
Unpack a given microslice. To be implemented in the derived unpacker algos.
void processHubBlock(CbmRichUnpackAlgoMicrosliceReader &reader)
void writeOutputDigi(Int_t fpgaID, Int_t channel, Double_t time, Double_t tot)
std::map< uint32_t, double > fLastCh0ReTime