REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawMultiCoBoAsAdToSignalProcess.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
57
58#include "TRestRawMultiCoBoAsAdToSignalProcess.h"
59
60#include "TRestDataBase.h"
61
62using namespace std;
63
64#include <bitset>
65
66#include "TTimeStamp.h"
67
69
70TRestRawMultiCoBoAsAdToSignalProcess::TRestRawMultiCoBoAsAdToSignalProcess() { Initialize(); }
71
72TRestRawMultiCoBoAsAdToSignalProcess::TRestRawMultiCoBoAsAdToSignalProcess(const char* configFilename) {
73 Initialize();
74}
75
76TRestRawMultiCoBoAsAdToSignalProcess::~TRestRawMultiCoBoAsAdToSignalProcess() {
77 // TRestRawMultiCoBoAsAdToSignalProcess destructor
78}
79
82
83 SetSectionName(this->ClassName());
84}
85
86Bool_t TRestRawMultiCoBoAsAdToSignalProcess::InitializeStartTimeStampFromFilename(TString fName) {
87 // these parameters have to be extracted from the file name. So do not change
88 // the origin binary file name.
89 int year, month, day, hour, minute, second, millisecond;
90
91 const Ssiz_t fnOffset = fName.Index(".graw");
92
93 if (fName.Length() != fnOffset + 5 || fnOffset < 28) {
94 cout << "Input binary file name type unknown!" << endl;
95 return kFALSE;
96 }
97
98 if (fName[fnOffset - 24] != '-' || fName[fnOffset - 21] != '-' || fName[fnOffset - 18] != 'T' ||
99 fName[fnOffset - 15] != ':' || fName[fnOffset - 12] != ':' || fName[fnOffset - 9] != '.' ||
100 fName[fnOffset - 5] != '_') {
101 cout << "Input binary file name unknown!" << endl;
102 return kFALSE;
103 }
104
105 year = (int)(fName[fnOffset - 28] - 48) * 1000 + (int)(fName[fnOffset - 27] - 48) * 100 +
106 (int)(fName[fnOffset - 26] - 48) * 10 + (int)(fName[fnOffset - 25] - 48) * 1;
107 month = (int)(fName[fnOffset - 23] - 48) * 10 + (int)(fName[fnOffset - 22] - 48) * 1;
108 day = (int)(fName[fnOffset - 20] - 48) * 10 + (int)(fName[fnOffset - 19] - 48) * 1;
109 hour = (int)(fName[fnOffset - 17] - 48) * 10 + (int)(fName[fnOffset - 16] - 48) * 1;
110 minute = (int)(fName[fnOffset - 14] - 48) * 10 + (int)(fName[fnOffset - 13] - 48) * 1;
111 second = (int)(fName[fnOffset - 11] - 48) * 10 + (int)(fName[fnOffset - 10] - 48) * 1;
112 millisecond = (int)(fName[fnOffset - 8] - 48) * 100 + (int)(fName[fnOffset - 7] - 48) * 10 +
113 (int)(fName[fnOffset - 6] - 48) * 1;
114
115 fStartTimeStamp.Set(year, month, day, hour, minute, second, millisecond * 1000000, kTRUE,
116 -8 * 3600); // Offset for Beijing(local) time
117 return kTRUE;
118}
119
120vector<int> fileerrors;
122 // fDataFrame.clear();
123 // fHeaderFrame.clear();
124 // fileerrors.clear();
125
126 // for (int n = 0; n < fInputFiles.size(); n++) {
127 // CoBoHeaderFrame hdrtmp;
128 // fHeaderFrame.push_back(hdrtmp);
129 // fileerrors.push_back(0);
130 //}
131
132 fRunOrigin = fRunInfo->GetRunNumber();
133 fCurrentEvent = -1;
134
135 if (fRunInfo->GetStartTimestamp() != 0) {
136 fStartTimeStamp = TTimeStamp(fRunInfo->GetStartTimestamp());
137 }
138
139 totalbytesRead = 0;
140}
141
142Bool_t TRestRawMultiCoBoAsAdToSignalProcess::AddInputFile(const string& file) {
143 if (file.find(".graw") == string::npos) {
144 return false;
145 }
146 if (TRestRawToSignalProcess::AddInputFile(file)) {
147 CoBoHeaderFrame hdrtmp;
148 fHeaderFrame.push_back(hdrtmp);
149 fileerrors.push_back(0);
150
151 int i = fHeaderFrame.size() - 1;
152 if (fread(fHeaderFrame[i].frameHeader, 256, 1, fInputFiles[i]) != 1 || feof(fInputFiles[i])) {
153 fclose(fInputFiles[i]);
154 fInputFiles[i] = nullptr;
155 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
156 return kFALSE;
157 }
158 totalbytesRead += 256;
159 if (!ReadFrameHeader(fHeaderFrame[i])) {
160 cout << "error when reading frame header in file " << i << " \"" << fInputFileNames[i] << "\""
161 << endl;
162 cout << "event id " << fCurrentEvent + 1 << ". The file will be closed" << endl;
163 fHeaderFrame[i].Show();
164 cout << endl;
165 GetChar();
166 fclose(fInputFiles[i]);
167 fInputFiles[i] = nullptr;
168 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
169 return false;
170 }
171
172 return true;
173 }
174 return false;
175}
176
178 fSignalEvent->Initialize();
179
180 if (EndReading()) {
181 return nullptr;
182 }
183 if (!FillBuffer()) {
184 fSignalEvent->SetOK(false);
185 return fSignalEvent;
186 }
187
188 // Int_t nextId = GetLowestEventId();
189
191 cout << "TRestRawMultiCoBoAsAdToSignalProcess: Generating event with ID: " << fCurrentEvent << endl;
192 }
193
194 TTimeStamp tSt = 0;
195
196 map<int, CoBoDataFrame>::iterator it;
197 it = fDataFrame.begin();
198
199 while (it != fDataFrame.end()) {
200 CoBoDataFrame* data = &(it->second);
201 if (data->evId == fCurrentEvent) {
202 if ((Double_t)tSt == 0) tSt = data->timeStamp;
203
204 for (int m = 0; m < 272; m++) {
205 if (data->chHit[m]) {
206 signal.Initialize();
207 signal.SetSignalID(m + data->asadId * 272);
208
209 for (int j = 0; j < 512; j++) signal.AddPoint((Short_t)data->data[m][j]);
210
211 fSignalEvent->AddSignal(signal);
212
214 cout << "AgetId, chnId, first value, max value: " << m / 68 << ", " << m % 68 << ", "
215 << signal.GetData(0) << ", " << signal.GetMaxValue() << endl;
216 }
217 }
218 }
219 data->evId = -1;
220 for (int m = 0; m < 272; m++) data->chHit[m] = kFALSE;
221 }
222 it++;
223 }
224
226 cout << "TRestRawMultiCoBoAsAdToSignalProcess: event time is : " << tSt << endl;
227 cout << "TRestRawMultiCoBoAsAdToSignalProcess: " << fSignalEvent->GetNumberOfSignals()
228 << " signals added" << endl;
229 cout << "------------------------------------" << endl;
230 }
231 fSignalEvent->SetTimeStamp(tSt);
232 fSignalEvent->SetID(fCurrentEvent);
233 fSignalEvent->SetRunOrigin(0);
234 fSignalEvent->SetSubRunOrigin(0);
235
236 // cout << fSignalEvent->GetNumberOfSignals() << endl;
237 // if( fSignalEvent->GetNumberOfSignals( ) == 0 ) return nullptr;
238
239 return fSignalEvent;
240}
241
243 for (unsigned int i = 0; i < fileerrors.size(); i++) {
244 if (fileerrors[i] > 0) {
245 RESTWarning << "Found " << fileerrors[i] << " error frame headers in file " << i << RESTendl;
246 RESTWarning << "\"" << fInputFileNames[i] << "\"" << RESTendl;
247 }
248 }
249
250 fileerrors.clear();
251}
252
253// true: finish filling
254// false: error when filling
255bool TRestRawMultiCoBoAsAdToSignalProcess::FillBuffer() {
256 // if the file is opened but not read, read header frame
257 for (unsigned int i = 0; i < fInputFiles.size(); i++) {
258 if (fInputFiles[i] && ftell(fInputFiles[i]) == 0) {
259 if (fread(fHeaderFrame[i].frameHeader, 256, 1, fInputFiles[i]) != 1 || feof(fInputFiles[i])) {
260 fclose(fInputFiles[i]);
261 fInputFiles[i] = nullptr;
262 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
263 return kFALSE;
264 }
265 totalbytesRead += 256;
266 if (!ReadFrameHeader(fHeaderFrame[i])) {
267 cout << "error when reading frame header in file " << i << " \"" << fInputFileNames[i] << "\""
268 << endl;
269 cout << "event id " << fCurrentEvent + 1 << ". The file will be closed" << endl;
270 fHeaderFrame[i].Show();
271 cout << endl;
272 GetChar();
273 fclose(fInputFiles[i]);
274 fInputFiles[i] = nullptr;
275 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
276 return false;
277 }
278 }
279 }
280
281 // normally:
282 // 1.use the smallest event id in header frames as current event id
283 unsigned int evt = fHeaderFrame[0].eventIdx;
284 for (unsigned int i = 1; i < fHeaderFrame.size(); i++) {
285 if (fHeaderFrame[i].eventIdx < evt) evt = fHeaderFrame[i].eventIdx;
286 }
287 fCurrentEvent = evt;
288
289 // loop for each file
290 for (unsigned int i = 0; i < fHeaderFrame.size(); i++) {
291 if (fInputFiles[i] == nullptr) {
292 continue;
293 }
294
295 // file position is at the end of last header
296 // if the eventid of last header is same as current, do the following:
297 // a. read the data frame behind
298 // b. read the next frame header
299 // c. if eventid is the same as current, return to a, otherwise break.
300 while (fCurrentEvent >= 0 && fHeaderFrame[i].eventIdx == (unsigned int)fCurrentEvent) {
302 cout << "TRestRawMultiCoBoAsAdToSignalProcess: retrieving frame header in "
303 "file "
304 << i << " (" << fInputFileNames[i] << ")" << endl;
306 fHeaderFrame[i].Show();
307 }
308
309 // reading data according to the header
310 unsigned int type = fHeaderFrame[i].frameType;
311 if (fHeaderFrame[i].frameHeader[0] == 0x08 && type == 1) // partial readout
312 {
313 ReadFrameDataP(fInputFiles[i], fHeaderFrame[i]);
314 } else if (fHeaderFrame[i].frameHeader[0] == 0x08 && type == 2) // full readout
315 {
316 if (fread(frameDataF, 2048, 136, fInputFiles[i]) != 136 || feof(fInputFiles[i])) {
317 fclose(fInputFiles[i]);
318 fInputFiles[i] = nullptr;
319 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
320 break;
321 }
322 totalbytesRead += 278528;
323 ReadFrameDataF(fHeaderFrame[i]);
324 } else {
325 fclose(fInputFiles[i]);
326 fInputFiles[i] = nullptr;
327 fHeaderFrame[i].eventIdx = (unsigned int)4294967295;
328 return false;
329 }
330
331 // reading next header
332 if (fread(fHeaderFrame[i].frameHeader, 256, 1, fInputFiles[i]) != 1 || feof(fInputFiles[i])) {
333 fclose(fInputFiles[i]);
334 fInputFiles[i] = nullptr;
335 fHeaderFrame[i].eventIdx = (unsigned int)4294967295; // maximum of unsigned int
336 break;
337 }
338 totalbytesRead += 256;
339 if (!ReadFrameHeader(fHeaderFrame[i])) {
340 RESTWarning << "Event " << fCurrentEvent << " : error when reading next frame header"
341 << RESTendl;
342 RESTWarning << "in file " << i << " \"" << fInputFileNames[i] << "\"" << RESTendl;
344 RESTWarning << "trying to skip this event and find next header..." << RESTendl;
345 fileerrors[i] += 1;
347 bool found = false;
349 for (int k = 0; k < 1088; k++) // fullreadoutsize(278528)/headersize(256)=1088
350 {
351 if (fread(fHeaderFrame[i].frameHeader, 256, 1, fInputFiles[i]) != 1 ||
352 feof(fInputFiles[i])) {
353 break;
354 }
355 totalbytesRead += 256;
356 if (ReadFrameHeader(fHeaderFrame[i])) {
357 fVerboseLevel = tmp;
358 RESTWarning << "Successfully found next header (EventId : "
359 << fHeaderFrame[i].eventIdx << ")" << RESTendl;
361 fHeaderFrame[i].Show();
362 cout << endl;
363 // GetChar();
364 found = true;
365 fSignalEvent->SetOK(false);
366 break;
367 }
368 }
369 if (!found) {
370 fclose(fInputFiles[i]);
371 fInputFiles[i] = nullptr;
372 fHeaderFrame[i].eventIdx = (unsigned int)4294967295; // maximum of unsigned int
373 }
374 }
375 }
376 }
377
378 return true;
379}
380
381bool TRestRawMultiCoBoAsAdToSignalProcess::ReadFrameHeader(CoBoHeaderFrame& HdrFrame) {
382 UChar_t* Header = &(HdrFrame.frameHeader[0]);
383
384 HdrFrame.frameSize =
385 (unsigned int)Header[1] * 0x10000 + (unsigned int)Header[2] * 0x100 + (unsigned int)Header[3];
386 HdrFrame.frameSize *= 256;
387
388 HdrFrame.frameType = (unsigned int)Header[5] * 0x100 + (unsigned int)Header[6];
389 HdrFrame.revision = (unsigned int)Header[7];
390 HdrFrame.headerSize = (unsigned int)Header[8] * 0x100 + (unsigned int)Header[9];
391 HdrFrame.itemSize = (unsigned int)Header[10] * 0x100 + (unsigned int)Header[11];
392 HdrFrame.nItems = (unsigned int)Header[12] * 0x1000000 + (unsigned int)Header[13] * 0x10000 +
393 (unsigned int)Header[14] * 0x100 + (unsigned int)Header[15];
394 HdrFrame.eventTime = (Long64_t)Header[16] * 0x10000000000 + (Long64_t)Header[17] * 0x100000000 +
395 (Long64_t)Header[18] * 0x1000000 + (Long64_t)Header[19] * 0x10000 +
396 (Long64_t)Header[20] * 0x100 + (Long64_t)Header[21];
397 HdrFrame.eventTime *= 10; // ns at 100MHz experiment clock
398 HdrFrame.eventIdx = (unsigned int)Header[22] * 0x1000000 + (unsigned int)Header[23] * 0x10000 +
399 (unsigned int)Header[24] * 0x100 + (unsigned int)Header[25];
400
401 HdrFrame.asadIdx = (unsigned int)Header[27];
402 HdrFrame.readOffset = (unsigned int)Header[28] * 0x100 + (unsigned int)Header[29];
403 HdrFrame.status = (unsigned int)Header[30];
404
405 // if (fCurrentEvent == -1) { fCurrentEvent = HdrFrame.eventIdx; }
406 // if (fCurrentEvent != HdrFrame.eventIdx) { fNextEvent = HdrFrame.eventIdx;
407 // return 1; }
408
409 // HdrFrame.fEveTimeSec = HdrFrame.eventTime / (Long64_t)1e9;//relative time
410 // HdrFrame.fEveTimeNanoSec = HdrFrame.eventTime % (Long64_t)1e9;
411
412 // HdrFrame.fEveTimeNanoSec = (int)HdrFrame.eventTime;
413 // HdrFrame.fEveTimeStamp.SetSec(HdrFrame.fEveTimeSec);
414 // HdrFrame.fEveTimeStamp.SetNanoSec(HdrFrame.fEveTimeNanoSec);
415 // HdrFrame.fEveTimeStamp.Add(fStartTimeStamp);
416
417 if (HdrFrame.frameType == 1) {
418 if (HdrFrame.itemSize != 4) {
419 RESTWarning << "unsupported item size!" << RESTendl;
420 return false;
421 }
422
423 } else if (HdrFrame.frameType == 2) {
424 if (HdrFrame.itemSize != 2) {
425 RESTWarning << "unsupported item size!" << RESTendl;
426 return false;
427 }
428 if (HdrFrame.nItems != 139264) {
429 RESTWarning << "unsupported nItems!" << RESTendl;
430 return false;
431 }
432 } else {
433 RESTWarning << "unknown frame type" << RESTendl;
434 return false;
435 }
436
437 // warning<<"revision: "<<revision<<endl;
438 if (HdrFrame.revision != 5) {
439 RESTWarning << "unsupported revision!" << RESTendl;
440 return false;
441 }
442
443 // warning<<"frameHeaderSize: "<<frameHeaderSize<<endl;
444 if (HdrFrame.headerSize != 1) {
445 RESTWarning << "unsupported frameHeader size!" << RESTendl;
446 return false;
447 }
448
449 // warning<<"readOffset: "<<readOffset<<endl;
450 if (HdrFrame.readOffset != 0) {
451 RESTWarning << "unsupported readOffset!" << RESTendl;
452 return false;
453 }
454
455 if (HdrFrame.status) {
456 RESTWarning << "bad frame!" << RESTendl;
457 return false;
458 }
459
460 if (HdrFrame.nItems * HdrFrame.itemSize + 256 != HdrFrame.frameSize) {
461 RESTWarning << "Event " << fCurrentEvent << " : item number and frame size unmatch!" << RESTendl;
462
463 // sometimes there is a itemnumber-framesize unmatch problem
464 // nItems*itemSize(=4)+256=frameSize
465 // because nitems/512 should be a integer(signal number), we can make a fix
466 //{
467 // if (HdrFrame.nItems % 512 == 0) {
468 // warning << "...frameSize (" << HdrFrame.frameSize;
469 // HdrFrame.frameSize = HdrFrame.nItems * HdrFrame.itemSize + 256;
470 // warning << ") fixed to " << HdrFrame.frameSize << "..." << endl;
471 // warning << endl;
472 //}
473 // else
474 //{
475 // if(fVerboseLevel>=REST_Info)
476 // HdrFrame.Show();
477 // if (((HdrFrame.frameSize - 256) / HdrFrame.itemSize) % 512 == 0)
478 // {
479 // HdrFrame.nItems = (HdrFrame.frameSize - 256) /
480 // HdrFrame.itemSize; warning << "...nItems fixed to " << HdrFrame.nItems
481 // <<
482 //"..." << endl; warning << endl;
483 // }
484 // else
485 // {
486 // warning << "frame unfixed" << endl;
487 // warning << endl;
488 // }
489 //}
490 }
491
492 return true;
493}
494
495bool TRestRawMultiCoBoAsAdToSignalProcess::ReadFrameDataP(FILE* f, CoBoHeaderFrame& hdr) {
496 unsigned int i;
497 unsigned int agetIdx, chanIdx, buckIdx, sample, chTmp;
498
499 unsigned int asadid = hdr.asadIdx;
500 unsigned int size = hdr.frameSize;
501 // unsigned int items = hdr.nItems;
502 unsigned int eventid = hdr.eventIdx;
503 Long64_t time = hdr.eventTime;
504 TTimeStamp eveTimeStamp;
505 CoBoDataFrame& dataf = fDataFrame[asadid];
506
507 //------------read frame data-----------
508 if (size > 256) {
509 unsigned int NBuckTotal = (size - 256) / 4;
510 for (i = 0; i < NBuckTotal; i++) {
511 if ((fread(frameDataP, 4, 1, f)) != 1 || feof(f)) {
512 fclose(f);
513 f = nullptr;
514 return kFALSE;
515 }
516 totalbytesRead += 4;
517 // total: 4bytes, 32 bits
518 // 11 111111|1 1111111|11 11 1111|11111111
519 // agetIdx chanIdx buckIdx unused samplepoint
520 agetIdx = (frameDataP[0] >> 6); // first 2 bits of the byte
521 chanIdx = ((unsigned int)(frameDataP[0] & 0x3f) * 2 + (frameDataP[1] >> 7));
522 chTmp = agetIdx * 68 + chanIdx;
523 buckIdx = ((unsigned int)(frameDataP[1] & 0x7f) * 4 + (frameDataP[2] >> 6));
524 sample = ((unsigned int)(frameDataP[2] & 0x0f) * 0x100 + frameDataP[3]);
525
526 if (chTmp >= 272) {
527 cout << "channel id error! value: " << chTmp << endl;
528 continue;
529 }
530
531 dataf.chHit[chTmp] = kTRUE;
532 dataf.data[chTmp][buckIdx] = sample;
533 }
534 }
535
536 eveTimeStamp.SetNanoSec(time % ((Long64_t)1e9));
537 eveTimeStamp.SetSec(time / ((Long64_t)1e9));
538 eveTimeStamp.Add(fStartTimeStamp);
539
540 dataf.asadId = asadid;
541 dataf.evId = eventid;
542 dataf.timeStamp = eveTimeStamp;
543
544 return true;
545}
546
547bool TRestRawMultiCoBoAsAdToSignalProcess::ReadFrameDataF(CoBoHeaderFrame& hdr) {
548 int i;
549 int j;
550 unsigned int agetIdx, chanIdx, chanIdx0, chanIdx1, chanIdx2, chanIdx3, sample, chTmp;
551
552 unsigned int asadid = hdr.asadIdx;
553 unsigned int eventid = hdr.eventIdx;
554 Long64_t time = hdr.eventTime;
555 TTimeStamp eveTimeStamp;
556 CoBoDataFrame& dataf = fDataFrame[asadid];
557
558 int tmpP;
559 for (i = 0; i < 512; i++) {
560 chanIdx0 = 0;
561 chanIdx1 = 0;
562 chanIdx2 = 0;
563 chanIdx3 = 0;
564 for (j = 0; j < 272; j++) {
565 // total 8*2= 16 bits
566 // 11 11 1111|11111111
567 // agetIdx unused samplepoint
568
569 tmpP = (i * 272 + j) * 2;
570 agetIdx = (frameDataF[tmpP] >> 6);
571 sample = ((unsigned int)(frameDataF[tmpP] & 0x0f) * 0x100 + frameDataF[tmpP + 1]);
572
573 if (agetIdx == 0) {
574 chanIdx = chanIdx0;
575 chanIdx0++;
576 } else if (agetIdx == 1) {
577 chanIdx = chanIdx1;
578 chanIdx1++;
579 } else if (agetIdx == 2) {
580 chanIdx = chanIdx2;
581 chanIdx2++;
582 } else {
583 chanIdx = chanIdx3;
584 chanIdx3++;
585 }
586 // cout<<"agetIdx: "<<agetIdx<<" chanIdx: "<<chanIdx<<endl;
587
588 if (chanIdx > 67) { /*cout << "buck or channel id error! ChannelId,
589 AgetId: " << chanIdx << ", " << agetIdx << endl;*/
590 continue;
591 }
592 chTmp = agetIdx * 68 + chanIdx;
593 dataf.chHit[chTmp] = kTRUE;
594 dataf.data[chTmp][i] = sample;
595 }
596 }
597
598 eveTimeStamp.SetNanoSec(time % ((Long64_t)1e9));
599 eveTimeStamp.SetSec(time / ((Long64_t)1e9));
600 eveTimeStamp.Add(fStartTimeStamp);
601
602 dataf.asadId = asadid;
603 dataf.evId = eventid;
604 dataf.timeStamp = eveTimeStamp;
605
606 return true;
607}
608
609Bool_t TRestRawMultiCoBoAsAdToSignalProcess::EndReading() {
610 for (auto& m : fDataFrame) {
611 m.second.finished = true;
612 }
613
614 // cout << "header frame: ";
615 // for (int n = 0; n < nFiles; n++) {
616 // cout << fHeaderFrame[n].asadIdx << ":" << fHeaderFrame[n].eventIdx << ", ";
617 //}
618 // cout << endl;
619
620 for (int n = 0; n < nFiles; n++) {
621 // if one header is not 42949..., the asad chain is not finished
622 fDataFrame[fHeaderFrame[n].asadIdx].finished =
623 (fDataFrame[fHeaderFrame[n].asadIdx].finished &&
624 (fHeaderFrame[n].eventIdx == (unsigned int)4294967295));
625 }
626
627 // cout << "data frame: ";
628 // for (auto m : fDataFrame) {
629 // cout << m.second.asadId << ":" << m.second.finished << ", ";
630 //}
631 // cout << endl;
632 // cout << endl;
633
634 for (auto m : fDataFrame) {
635 // if any of the asad chain is finihsed, we ends reading for all the
636 // events
637 if (m.second.asadId == -1) continue;
638
639 if (m.second.finished == true) {
640 return true;
641 }
642 }
643
644 for (const auto& file : fInputFiles) {
645 if (file != nullptr) {
646 return false;
647 }
648 }
649
650 return true;
651}
TRestRun * fRunInfo
< Pointer to TRestRun object where to find metadata.
A base class for any REST event.
Definition: TRestEvent.h:38
endl_t RESTendl
Termination flag object for TRestStringOutput.
TRestStringOutput::REST_Verbose_Level fVerboseLevel
Verbose level used to print debug info.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
void SetSectionName(std::string sName)
set the section name, clear the section content
UChar_t frameDataF[278528]
///for partial readout data frame
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void Initialize() override
Making default settings.
int fCurrentEvent
///reserves a header frame for each file
TTimeStamp fStartTimeStamp
///for full readout data frame
void EndProcess() override
To be executed at the end of the run (outside event loop)
void InitProcess() override
To be executed at the beginning of the run (outside event loop)
std::vector< CoBoHeaderFrame > fHeaderFrame
///asadId, dataframe
Double_t GetMaxValue()
Returns the maximum value found in the data points. It includes baseline correction.
void Initialize()
Initialization of TRestRawSignal members.
void SetSignalID(Int_t sID)
It sets the id number of the signal.
Double_t GetData(Int_t n) const
It returns the data value of point n including baseline correction.
void AddPoint(Short_t)
Adds a new point to the end of the signal data array.
void Initialize() override
Making default settings.
REST_Verbose_Level
Enumerate of verbose level, containing five levels.
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages
@ REST_Silent
show minimum information of the software, as well as error messages
Int_t GetChar(std::string hint="Press a KEY to continue ...")
Helps to pause the program, printing a message before pausing.