REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawFEUDreamToSignalProcess.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
24// The TRestRawFEUDreamToSignalProcess ...
25//
26// DOCUMENTATION TO BE WRITTEN (main description, methods, data members)
27//
28// TODO: This process requires optimization to improve the data processing
29// rate.
30//
31// \warning This process might be obsolete today. It may need additional
32// revision, validation, and documentation. Use it under your own risk. If you
33// find this process useful for your work feel free to use it, improve it,
34// validate and/or document this process. If all those points are addressed
35// these lines can be removed.
36//
37// <hr>
38//
39// \warning **⚠ REST is under continous development.** This
40// documentation
41// is offered to you by the REST community. Your HELP is needed to keep this
42// code
43// up to date. Your feedback will be worth to support this software, please
44// report
45// any problems/suggestions you may find will using it at [The REST Framework
46// forum](http://ezpc10.unizar.es). You are welcome to contribute fixing typos,
47// updating
48// information or adding/proposing new contributions. See also our
49// <a href="https://github.com/rest-for-physics/framework/blob/master/CONTRIBUTING.md">Contribution
50// Guide</a>.
51//
52//
53//--------------------------------------------------------------------------
54//
55// RESTsoft - Software for Rare Event Searches with TPCs
56//
57// History of developments:
58//
59// 2019-May: First implementation
60// Damien Neyret
61//
62// 2021-May: Readapted to compile in REST v2.3.X
63// Damien Neyret
64//
65// \class TRestRawFEUDreamToSignalProcess
66// \author Damien Neyret
67// \author Javier Galan
68//
69// <hr>
71#include "TRestRawFEUDreamToSignalProcess.h"
72
73using namespace std;
74
75#include "TTimeStamp.h"
76
78
79TRestRawFEUDreamToSignalProcess::TRestRawFEUDreamToSignalProcess() { Initialize(); }
80
81TRestRawFEUDreamToSignalProcess::TRestRawFEUDreamToSignalProcess(const char* configFilename)
82 : TRestRawToSignalProcess(configFilename) {
83 Initialize();
84}
85
86TRestRawFEUDreamToSignalProcess::~TRestRawFEUDreamToSignalProcess() {
87 // TRestRawFEUDreamToSignalProcess destructor
88}
89
92
93 // MaxThreshold = 4000;
94 bad_event = false;
95 line = 0; // line number
96 Nevent = 0, Nbadevent = 0; // current event number
97 IDEvent = 0;
98}
99
101 tStart = 0; // timeStamp of the run initially set to 0
102 RESTInfo << "TRestRawFEUDreamToSignalProcess::InitProcess" << RESTendl;
103
104 totalBytesReaded = 0;
105}
106
108 FeuReadOut Feu;
109 bool badreadfg = false;
110
111 RESTDebug << "---------------Start of TRestRawFEUDreamToSignalProcess::ProcessEvent------------"
112 << RESTendl;
113
114 fSignalEvent->Initialize();
115
116 while (true) { // loop on events
117
118 Feu.NewEvent(); // reset Feu structure
119
120 // Check header and fill Event IDs and TimeStamps
121 badreadfg = ReadFeuHeaders(Feu);
122 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: header read, badreadfg " << badreadfg
123 << RESTendl;
124 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: event to read EventID " << Feu.EventID
125 << " Time " << Feu.TimeStamp << " isample " << Feu.isample << RESTendl;
126
127 if (badreadfg) {
128 RESTWarning
129 << "TRestRawFEUDreamToSignalProcess::ProcessEvent: Error in reading feu header (bad file or "
130 "end of file), trying to go to the next file"
131 << RESTendl;
132 if (GoToNextFile()) {
133 badreadfg = ReadFeuHeaders(Feu); // reading event from the next file
134 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: header read, badreadfg "
135 << badreadfg << RESTendl;
136 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: event to read EventID "
137 << Feu.EventID << " Time " << Feu.TimeStamp << " isample " << Feu.isample
138 << RESTendl;
139 } else {
140 return nullptr;
141 }
142 }
143
144 // Read event
145 badreadfg = ReadEvent(Feu);
146 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: event read, badreadfg " << badreadfg
147 << RESTendl;
148 if (badreadfg) {
149 RESTError << "TRestRawFEUDreamToSignalProcess::ProcessEvent: Error in event reading at event "
150 << Nevent << RESTendl;
151 break;
152 }
153
154 Nevent++;
155 if (badreadfg) break; // break from loop
156 if (bad_event) Nbadevent++;
157 if ((Nevent % 100) == 0)
158 cout << "TRestRawFEUDreamToSignalProcess::ProcessEvent: " << Nevent
159 << " events processed in file, and " << Nbadevent << " bad events skipped " << endl;
160
162 RESTInfo << "-- TRestRawFEUDreamToSignalProcess::ProcessEvent ---" << RESTendl;
163 RESTInfo << "Event ID : " << fSignalEvent->GetID() << RESTendl;
164 RESTInfo << "Time stamp : " << fSignalEvent->GetTimeStamp() << RESTendl;
165 RESTInfo << "Number of Signals : " << fSignalEvent->GetNumberOfSignals() << RESTendl;
166 RESTInfo << "Number of Samples : " << (Feu.isample + 1) << RESTendl;
167 RESTInfo << "-------------------------------------------------" << RESTendl;
168 }
169
170 Feu.isample = -1;
171 Feu.isample_prev = -2;
172 Feu.NewEvent();
173 bad_event = false;
174
175 if (fSignalEvent->GetNumberOfSignals() == 0) {
176 RESTError << "TRestRawFEUDreamToSignalProcess::ProcessEvent: no signal in event" << RESTendl;
177 return nullptr;
178 }
179
180 RESTDebug << "TRestRawFEUDreamToSignalProcess::ProcessEvent: returning signal event fSignalEvent "
181 << fSignalEvent << RESTendl;
182 if (GetVerboseLevel() > TRestStringOutput::REST_Verbose_Level::REST_Debug) fSignalEvent->PrintEvent();
183 return fSignalEvent;
184 }
185
186 return nullptr; // can't read data
187}
188
189// Definition of decoding methods
190bool TRestRawFEUDreamToSignalProcess::ReadEvent(FeuReadOut& Feu) {
191 bool badreadfg = false;
192
193 while (!Feu.event_completed) {
194 badreadfg = ReadFeuHeaders(Feu); // read feu header if not done
195 RESTDebug << "TRestRawFEUDreamToSignalProcess::ReadEvent: header read, badreadfg " << badreadfg
196 << RESTendl;
197 if (badreadfg) {
198 RESTWarning << "TRestRawFEUDreamToSignalProcess::ReadEvent: error in reading FEU headers "
199 << RESTendl; // failed
200 return true;
201 }
202
203 badreadfg = ReadDreamData(Feu); // read dream data
204 RESTDebug << "TRestRawFEUDreamToSignalProcess::ReadEvent: data read, badreadfg " << badreadfg
205 << RESTendl;
206 if (badreadfg) {
207 RESTError << "TRestRawFEUDreamToSignalProcess::ReadEvent: error in reading Dream data "
208 << RESTendl;
209 return true;
210 }
211
212 badreadfg = ReadFeuTrailer(Feu); // read feu trailer
213 RESTDebug << "TRestRawFEUDreamToSignalProcess::ReadEvent: trailer read, badreadfg " << badreadfg
214 << RESTendl;
215 if (badreadfg) {
216 RESTError << "TRestRawFEUDreamToSignalProcess::ReadEvent: error in reading FEU trailer"
217 << RESTendl;
218 return true;
219 }
220
221 } // end loop
222 RESTInfo << "TRestRawFEUDreamToSignalProcess::ReadEvent: Event ID " << Feu.EventID
223 << " processed successfully, Time " << Feu.TimeStamp << " isample " << Feu.isample << RESTendl;
224
225 return false;
226}
227
228bool TRestRawFEUDreamToSignalProcess::ReadFeuHeaders(FeuReadOut& Feu) {
229 if (Feu.FeuHeaderLoaded) return false; // already done
230
231 if (!Feu.data_to_treat) { // data not loaded
232
233 int nbytes = fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile);
234 totalBytesReaded += sizeof(Feu.current_data);
235 if (nbytes == 0) {
236 // perror("TRestRawFEUDreamToSignalProcess::ReadFeuHeaders: Error in reading FeuHeaders !");
237 RESTWarning
238 << "TRestRawFEUDreamToSignalProcess::ReadFeuHeaders: Problem in reading raw file, ferror "
239 << ferror(fInputBinFile) << " feof " << feof(fInputBinFile) << " fInputBinFile "
240 << fInputBinFile << RESTendl;
241 // fclose(fInputBinFile);
242 return true; // failed
243 }
244 // debug<<" Reading FeuHeaders ok, nbytes "<<nbytes<<endl;
245 Feu.current_data.ntohs_();
246 Feu.data_to_treat = true;
247 }
248
249 while (true) { // loop on words of the header
250 if (Feu.FeuHeaderLine < 8 && Feu.current_data.is_Feu_header()) {
251 if (Feu.FeuHeaderLine == 0) {
252 Feu.zs_mode = Feu.current_data.get_zs_mode();
253 Feu.Id = Feu.current_data.get_Feu_ID();
254 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " ZS mode "
255 << Feu.zs_mode << " Id " << Feu.Id << RESTendl;
256
257 } else if (Feu.FeuHeaderLine == 1) {
258 Feu.EventID = Feu.current_data.get_data();
259 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " EventID "
260 << Feu.EventID << RESTendl;
261
262 } else if (Feu.FeuHeaderLine == 2) {
263 Feu.TimeStamp = Feu.current_data.get_data();
264 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " TimeStamp "
265 << Feu.TimeStamp << RESTendl;
266
267 } else if (Feu.FeuHeaderLine == 3) {
268 Feu.FineTimeStamp = Feu.current_data.get_finetstp();
269 Feu.isample_prev = Feu.isample;
270 Feu.isample = Feu.current_data.get_sample_ID();
271 // fprintf(stderr, "sampleIDprevious %d, sampleID %d \n", Feu.isample_prev,Feu.isample);
272 Feu.FeuHeaderLoaded = true;
273 if (Feu.isample != Feu.isample_prev + 1) {
274 if (Feu.isample_prev ==
275 fMinPoints - 2) { // finishing the current event and starting the next one
276 // fprintf(stderr, "Event ID %d, processed \n", Feu.EventID-1);
277 } else {
278 RESTError
279 << "TRestRawFEUDreamToSignalProcess::ReadFeuHeaders: non continuous sample index "
280 "number, isample = "
281 << Feu.isample << " prev_isample = " << Feu.isample_prev << RESTendl;
282 bad_event = true;
283 }
284 }
285 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " FineTimeStamp "
286 << Feu.FineTimeStamp << " isample " << Feu.isample << RESTendl;
287
288 // Reading optionals
289 } else if (Feu.FeuHeaderLine == 4) {
290 Feu.EventID_Op = Feu.current_data.get_data();
291 Feu.EventID += (1 << 12) * Feu.EventID_Op;
292 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " EventID_Op "
293 << Feu.EventID_Op << RESTendl;
294
295 } else if (Feu.FeuHeaderLine == 5) {
296 Feu.TimeStamp_Op1 = Feu.current_data.get_data();
297 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " TimeStamp_Op1 "
298 << Feu.TimeStamp_Op1 << RESTendl;
299
300 } else if (Feu.FeuHeaderLine == 6) {
301 Feu.TimeStamp_Op2 = Feu.current_data.get_data();
302 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " TimeStamp_Op2 "
303 << Feu.TimeStamp_Op2 << RESTendl;
304
305 } else if (Feu.FeuHeaderLine == 7) {
306 Feu.TimeStamp_Op3 = Feu.current_data.get_TimeStamp_Op();
307 Feu.TimeStamp += ((long long)1 << 36) * Feu.TimeStamp_Op3 + (1 << 24) * Feu.TimeStamp_Op2 +
308 (1 << 12) * Feu.TimeStamp_Op1;
309 RESTDebug << "ReadFeuHeaders: header FeuHeaderLine " << Feu.FeuHeaderLine << " TimeStamp_Op3 "
310 << Feu.TimeStamp_Op3 << " TimeStamp " << Feu.TimeStamp << RESTendl;
311 }
312
313 Feu.data_to_treat = false;
314 Feu.FeuHeaderLine++;
315
316 } else if (Feu.FeuHeaderLine > 8 && Feu.current_data.is_Feu_header()) {
317 RESTError << "TRestRawFEUDreamToSignalProcess::ReadFeuHeaders: too long Feu header part "
318 << Feu.FeuHeaderLine << RESTendl;
319 bad_event = true;
320 } else if (Feu.FeuHeaderLine > 3 && !Feu.current_data.is_Feu_header())
321 break; // header finished
322
323 if (fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile) == 0) return true;
324 totalBytesReaded += sizeof(Feu.current_data);
325 Feu.current_data.ntohs_();
326 Feu.data_to_treat = true;
327
328 } // end while
329
330 fSignalEvent->SetID(Feu.EventID);
331 fSignalEvent->SetTime(tStart + Feu.TimeStamp * 8.E-9); // timeStamp in seconds
332
333 return false;
334}
335
336bool TRestRawFEUDreamToSignalProcess::ReadDreamData(FeuReadOut& Feu) {
337 bool got_raw_data_header = false;
338 bool got_channel_id = false;
339 int ichannel = 0;
340
341 if (!Feu.FeuHeaderLoaded) { // already loaded
342 RESTError
343 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: error in ReadDreamData, Feu header not loaded"
344 << RESTendl;
345 return true;
346 }
347
348 if (!Feu.data_to_treat) { // no data to treat
349 int nbytes = fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile);
350 totalBytesReaded += sizeof(Feu.current_data);
351 if (nbytes == 0) {
352 perror("TRestRawFEUDreamToSignalProcess::ReadDreamData: no Dream data to read in file");
353 RESTError << "TRestRawFEUDreamToSignalProcess::ReadDreamData: problem in reading raw data file, "
354 "ferror "
355 << ferror(fInputBinFile) << " feof " << feof(fInputBinFile) << " fInputBinFile "
356 << fInputBinFile << RESTendl;
357 fclose(fInputBinFile);
358 return true; // failed
359 }
360 // debug<<" Reading DreamData ok, nbytes "<<nbytes<<endl;
361 Feu.current_data.ntohs_();
362 Feu.data_to_treat = true;
363 }
364
365 while (true) { // loop on words in the Dream data main structure (not header or trailer ones)
366 if (Feu.FeuHeaderLine > 3 && !Feu.current_data.is_Feu_header()) {
367 if (Feu.DataHeaderLine < 4 && Feu.current_data.is_data_header()) { // data header treatment
368 if (Feu.DataHeaderLine == 0) {
369 Feu.TriggerID = Feu.current_data.get_data(); // trigger Id MSB
370 RESTDebug << "ReadDreamData: header DataHeaderLine " << Feu.DataHeaderLine
371 << " TriggerID MSB " << Feu.TriggerID << RESTendl;
372 } else if (Feu.DataHeaderLine == 1) {
373 Feu.TriggerID_ISB = Feu.current_data.get_data();
374 RESTDebug << "ReadDreamData: header DataHeaderLine " << Feu.DataHeaderLine
375 << " TriggerID_ISB " << Feu.TriggerID_ISB << RESTendl;
376 } else if (Feu.DataHeaderLine == 2) {
377 Feu.TriggerID_LSB = Feu.current_data.get_data();
378 Feu.TriggerID *= (1 << 24);
379 Feu.TriggerID += Feu.TriggerID_LSB + (1 << 12) * Feu.TriggerID_ISB;
380 RESTDebug << "ReadDreamData: header DataHeaderLine " << Feu.DataHeaderLine
381 << " TriggerID_LSB " << Feu.TriggerID_LSB << " TriggerID " << Feu.TriggerID
382 << RESTendl;
383 } else if (Feu.DataHeaderLine == 3) {
384 Feu.asicN = Feu.current_data.get_dream_ID(); // Dream_ID
385 ichannel = 0; // reset counting of channels for the current Dream chip
386 Feu.channelN = 0;
387 Feu.DataTrailerLine = 0;
388 // fprintf(stderr, " asic N %d \n", Feu.asicN);
389 got_raw_data_header = true;
390 RESTDebug << "ReadDreamData: header DataHeaderLine " << Feu.DataHeaderLine << " asicN "
391 << Feu.asicN << RESTendl;
392 }
393 Feu.DataHeaderLine++;
394 Feu.data_to_treat = false;
395
396 } else if (Feu.DataHeaderLine > 3 && Feu.current_data.is_data_header()) {
397 bad_event = true;
398 RESTError << "TRestRawFEUDreamToSignalProcess::ReadDreamData: too many data header lines, "
399 "DataHeaderLine "
400 << Feu.DataHeaderLine << RESTendl;
401 return true;
402
403 } else if (Feu.current_data.is_data() &&
404 !Feu.zs_mode) { // data lines treatment, non-zero suppression mode
405 if (!got_raw_data_header) {
406 bad_event = true;
407 RESTError
408 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: data lines without header in "
409 "non-ZS mode "
410 << RESTendl;
411 }
412 Feu.channelN = ichannel;
413 if (!bad_event && Feu.channelN > -1 && Feu.channelN < NstripMax &&
414 Feu.isample < fMinPoints) { // fMinPoints=128
415 Feu.physChannel = Feu.asicN * NstripMax + Feu.channelN; // channel's number on the DREAM
416
417 // loop on samples
418 if (Feu.physChannel < MaxPhysChannel) {
419 Int_t sgnlIndex = fSignalEvent->GetSignalIndex(Feu.physChannel);
420 if (sgnlIndex == -1) {
421 sgnlIndex = fSignalEvent->GetNumberOfSignals();
422 TRestRawSignal sgnl(fMinPoints);
423 sgnl.SetSignalID(Feu.physChannel);
424 fSignalEvent->AddSignal(sgnl);
425 }
426 fSignalEvent->AddChargeToSignal(Feu.physChannel, Feu.isample,
427 Feu.current_data.get_data());
428 } else
429 RESTError
430 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: too large physical Channel "
431 "in "
432 "non-ZS mode , Feu.physChannel= "
433 << Feu.physChannel << " > MaxPhysChannel " << MaxPhysChannel << RESTendl;
434 RESTExtreme << "ReadDreamData: nonZS physChannel " << Feu.physChannel << " get_data "
435 << Feu.current_data.get_data() << RESTendl;
436 }
437 ichannel++;
438 Feu.data_to_treat = false;
439
440 } else if (Feu.current_data.is_data_zs() && Feu.zs_mode) { // zero-suppression mode
441 if (got_raw_data_header) {
442 bad_event = true;
443 RESTError
444 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: data lines with header in ZS "
445 "mode "
446 << RESTendl;
447 }
448 if (!got_channel_id && Feu.current_data.is_channel_ID()) { // get channelID and dreamID
449 ichannel = Feu.current_data.get_channel_ID();
450 Feu.channelN = ichannel;
451 Feu.asicN = Feu.current_data.get_dream_ID(); // Dream_ID
452 Feu.physChannel = Feu.asicN * NstripMax + Feu.channelN; // channel's number on the DREAM
453 got_channel_id = true;
454 if (Feu.channelN < 0 || Feu.channelN >= NstripMax) {
455 RESTError
456 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: too large channel number in "
457 "ZS "
458 "mode , Feu.channelN= "
459 << Feu.channelN << " > MaxPhysChannel " << MaxPhysChannel << RESTendl;
460 bad_event = true;
461 }
462 RESTExtreme << "ReadDreamData: ZS header physChannel " << Feu.physChannel << " asicN "
463 << Feu.asicN << " ichannel " << ichannel << RESTendl;
464 } else { // get channel data
465 got_channel_id = false;
466 if (!bad_event && Feu.channelN > -1 && Feu.channelN < NstripMax &&
467 Feu.isample < fMinPoints) {
468 Feu.channel_data = Feu.current_data.get_data();
469 }
470 RESTExtreme << "ReadDreamData: ZS data get_data " << Feu.current_data.get_data()
471 << RESTendl;
472 }
473
474 if (Feu.physChannel < MaxPhysChannel) {
475 Int_t sgnlIndex = fSignalEvent->GetSignalIndex(Feu.physChannel);
476 if (sgnlIndex == -1) {
477 sgnlIndex = fSignalEvent->GetNumberOfSignals();
478 TRestRawSignal sgnl(fMinPoints);
479 sgnl.SetSignalID(Feu.physChannel);
480 fSignalEvent->AddSignal(sgnl);
481 }
482 fSignalEvent->AddChargeToSignal(Feu.physChannel, Feu.isample, Feu.channel_data);
483 } else
484 RESTError
485 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: too large physical Channel in ZS "
486 "mode , Feu.physChannel= "
487 << Feu.physChannel << " > MaxPhysChannel " << MaxPhysChannel << RESTendl;
488 Feu.data_to_treat = false;
489
490 } else if (Feu.current_data.is_data_trailer()) { // data trailer treatment
491
492 if (ichannel != NstripMax && !Feu.zs_mode) {
493 bad_event = true;
494 RESTError
495 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: trailer with missing channel "
496 "numbers in non-ZS mode, ichannel "
497 << ichannel << RESTendl;
498 return true;
499 }
500 if (got_channel_id && Feu.zs_mode) {
501 bad_event = true;
502 RESTError
503 << "TRestRawFEUDreamToSignalProcess::ReadDreamData: trailer with channel Id without "
504 "channel data in ZS data, got_channel_id true"
505 << RESTendl;
506 return true;
507 }
508
509 if (Feu.DataTrailerLine == 0) { // CMN
510 Feu.CMN = Feu.current_data.get_data();
511 RESTDebug << "ReadDreamData: trailer DataTrailerLine " << Feu.DataTrailerLine << " CMN "
512 << Feu.CMN << RESTendl;
513 } else if (Feu.DataTrailerLine == 1) {
514 Feu.CMN_rest = Feu.current_data.get_data();
515 RESTDebug << "ReadDreamData: trailer DataTrailerLine " << Feu.DataTrailerLine
516 << " CMN_rest " << Feu.CMN_rest << RESTendl;
517 } else if (Feu.DataTrailerLine == 2) { // Cell_ID
518 Feu.Cell_ID_MSB = Feu.current_data.get_data();
519 RESTDebug << "ReadDreamData: trailer DataTrailerLine " << Feu.DataTrailerLine
520 << " Cell_ID_MSB " << Feu.Cell_ID_MSB << RESTendl;
521 } else if (Feu.DataTrailerLine == 3) {
522 Feu.Cell_ID_ISB = Feu.current_data.get_data();
523 RESTDebug << "ReadDreamData: trailer DataTrailerLine " << Feu.DataTrailerLine
524 << " Cell_ID_ISB " << Feu.Cell_ID_ISB << RESTendl;
525 } else if (Feu.DataTrailerLine == 4) {
526 Feu.Cell_ID_LSB = Feu.current_data.get_data();
527 Feu.Cell_ID = Feu.Cell_ID_LSB + (1 << 12) * Feu.Cell_ID_ISB +
528 ((long long)1 << 24) * Feu.Cell_ID_MSB;
529 RESTDebug << "ReadDreamData: trailer DataTrailerLine " << Feu.DataTrailerLine
530 << " Cell_ID_LSB " << Feu.Cell_ID_LSB << " Cell_ID " << Feu.Cell_ID << RESTendl;
531 Feu.DataHeaderLine = 0;
532 got_raw_data_header = false;
533 Feu.channelN = 0;
534 }
535 Feu.DataTrailerLine++;
536 Feu.data_to_treat = false;
537
538 } else if (Feu.DataTrailerLine > 4 && Feu.current_data.is_data_trailer()) {
539 bad_event = true;
540 RESTError << "TRestRawFEUDreamToSignalProcess::ReadDreamData: too many data trailer lines, "
541 "DataTrailerLine "
542 << Feu.DataTrailerLine << RESTendl;
543 return true;
544
545 } else if (Feu.current_data.is_final_trailer())
546 break; // Dream raw data finished
547 }
548
549 if (fread((char*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile) == 0) return true;
550 totalBytesReaded += sizeof(Feu.current_data);
551 Feu.current_data.ntohs_();
552 Feu.data_to_treat = true;
553
554 } // end while
555 return false;
556}
557
558bool TRestRawFEUDreamToSignalProcess::ReadFeuTrailer(FeuReadOut& Feu) {
559 if (!Feu.data_to_treat) {
560 int nbytes = fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile);
561 totalBytesReaded += sizeof(Feu.current_data);
562 if (nbytes == 0) {
563 perror("TRestRawFEUDreamToSignalProcess::ReadFeuTrailer: can't read new data from file");
564 RESTError
565 << "TRestRawFEUDreamToSignalProcess::ReadFeuTrailer: can't read new data from file, ferror "
566 << ferror(fInputBinFile) << " feof " << feof(fInputBinFile) << " fInputBinFile "
567 << fInputBinFile << RESTendl;
568 fclose(fInputBinFile);
569 return true; // failed
570 }
571 RESTDebug << "TRestRawFEUDreamToSignalProcess::ReadFeuTrailer: Reading FeuTrailer ok, nbytes "
572 << nbytes << RESTendl;
573 Feu.current_data.ntohs_();
574 Feu.data_to_treat = true;
575 }
576
577 while (true) {
578 if (Feu.current_data.is_final_trailer()) {
579 if (Feu.channelN != 0) {
580 bad_event = true;
581 RESTError << "TRestRawFEUDreamToSignalProcess::ReadFeuTrailer: channel number not NULL in "
582 "trailer, "
583 "Feu.channelN "
584 << Feu.channelN << RESTendl;
585 return true;
586 }
587 if (Feu.current_data.is_end_of_event()) {
588 if (Feu.isample != (fMinPoints - 1)) {
589 RESTWarning
590 << "TRestRawFEUDreamToSignalProcess::ReadFeuTrailer: not all samples read at end of "
591 "event, isample "
592 << Feu.isample << " MinPoints " << fMinPoints << RESTendl;
593 }
594 // Feu.isample=-1; Feu.isample_prev=-2;
595 Feu.event_completed = true;
596 }
597
598 Feu.FeuHeaderLine = 0;
599 Feu.FeuHeaderLoaded = false;
600 Feu.zs_mode = false;
601 Feu.data_to_treat = false;
602
603 // Reading VEP, not used
604 int z = fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile);
605 if (z == 0)
606 RESTError << "TRestRawFEUDreamToSignalProcess::ReadFeuTrailer. Error reading file"
607 << RESTendl;
608 totalBytesReaded += sizeof(Feu.current_data);
609 break;
610 }
611
612 if (fread((void*)&(Feu.current_data), sizeof(Feu.current_data), 1, fInputBinFile) == 0) return true;
613 totalBytesReaded += sizeof(Feu.current_data);
614 Feu.current_data.ntohs_();
615 Feu.data_to_treat = true;
616
617 } // end while
618 return false;
619}
A base class for any REST event.
Definition: TRestEvent.h:38
void SetTime(Double_t time)
Definition: TRestEvent.cxx:85
endl_t RESTendl
Termination flag object for TRestStringOutput.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
An process to read binary data from FEUDream electronics.
void Initialize() override
Making default settings.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void InitProcess() override
To be executed at the beginning of the run (outside event loop)
void AddChargeToSignal(Int_t sgnlID, Int_t bin, Short_t value)
It defines a Short_t array with a physical parameter that evolves in time using a fixed time bin.
A base class for any process reading a binary external file as input to REST.
void Initialize() override
Making default settings.
@ REST_Info
+show most of the information for each steps
@ REST_Debug
+show the defined debug messages