REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawAFTERToSignalProcess.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
59
60#include "TRestRawAFTERToSignalProcess.h"
61
62#include <bitset>
63
64#include "TTimeStamp.h"
65#ifdef WIN32
66#include <Windows.h>
67#else
68#include <arpa/inet.h>
69#endif
70#include "mygblink.h"
71
72using namespace std;
73
75
80
93
98
108
109 prevTime = 0;
110 reducedTime = 0;
111}
112
121 // Binary file header
122
123 // The binary starts here
124 char runUid[21], initTime[21];
125 int z = fread(runUid, 1, 20, fInputBinFile);
126 if (z == 0) RESTError << "TRestRawAFTERToSignalProcess. Problems reading input file." << RESTendl;
127 runUid[20] = '\0';
128 sprintf(initTime, "%s", runUid);
129 printf("File UID is %s \n", initTime);
130 totalbytesRead = sizeof(runUid);
131
132 int year, day, month, hour, minute, second;
133 sscanf(runUid, "R%d.%02d.%02d-%02d:%02d:%02d", &year, &month, &day, &hour, &minute, &second);
134 printf("R%d_%02d_%02d-%02d_%02d_%02d\n", year, month, day, hour, minute, second);
135 TTimeStamp tS(year, month, day, hour, minute, second);
136 tStart = tS.AsDouble();
137 cout << tStart << endl;
138 // Timestamp of the run
139}
140
145 EventHeader head;
146 DataPacketHeader pHeader;
147 DataPacketEnd pEnd;
148
149 fSignalEvent->Initialize();
150
151 // Read next header or quit of end of file
152 if (fread(&head, sizeof(EventHeader), 1, fInputBinFile) != 1) {
153 fclose(fInputBinFile);
154 cout << "Error reading event header :-(" << endl;
155 cout << "... or end of file found :-)" << endl;
156 return nullptr;
157 }
158
159 head.eventSize = ntohl(head.eventSize);
160 head.eventNumb = ntohl(head.eventNumb);
161
162 payload = head.eventSize;
163 frameBits = sizeof(head);
164
165 RESTDebug << "Event number from header --> 0x" << std::hex << head.eventNumb << std::dec << RESTendl;
166 RESTDebug << " event header size " << sizeof(head) << RESTendl;
167 RESTDebug << " total rawdata size 0x" << std::hex << head.eventSize << std::dec << RESTendl;
168 RESTDebug << "Payload " << payload << RESTendl;
169
170 fSignalEvent->SetID(head.eventNumb);
171
172 int timeBin = 0;
173
174 int fecN;
175 int channel;
176 int asicN;
177 int physChannel;
178
179 uint32_t eventTime, deltaTime;
180 uint32_t th, tl;
181 int tempAsic1, tempAsic2, sampleCountRead, pay;
182 uint16_t data, dat;
183
184 bool isData = false;
185
186 bool first = true;
187
188 // Bucle till it finds the readed bits equals the payload
189 while (frameBits < payload) {
190 int x = fread(&pHeader, sizeof(DataPacketHeader), 1, fInputBinFile);
191 if (x == 0)
192 RESTError << "TRestRawAFTERToSignalProcess::ProcessEvent. Problems reading input file."
193 << RESTendl;
194 frameBits += sizeof(DataPacketHeader);
195
196 if (first) // Timestamping (A. Tomas, 30th June 2011)
197 {
198 th = ntohs(pHeader.ts_h);
199 tl = ntohs(pHeader.ts_l);
200 eventTime = th << 16 | tl; // built time from MSB and LSB
201
202 if (eventTime > prevTime)
203 deltaTime = eventTime - prevTime;
204 else
205 deltaTime = (0xFFFFFFFF - prevTime) + eventTime;
206
207 reducedTime += deltaTime;
208
209 // Set timestamp and event ID
210 fSignalEvent->SetTime(tStart + reducedTime * 2.E-8); // Set TimeStamp
211
212 prevTime = eventTime;
213 first = false;
214
215 RESTDebug << "Timestamp: " << eventTime << RESTendl;
216 }
217
218 RESTDebug << "******Event data packet header:******" << RESTendl;
219
220 RESTDebug << "Size " << ntohs(pHeader.size) << RESTendl;
221
222 RESTDebug << "Event data packet header: " << RESTendl;
223 RESTDebug << std::hex << "Size 0x" << ntohs(pHeader.size) << RESTendl;
224#ifdef NEW_DAQ_T2K_2_X
225 RESTDebug << "DCC 0x" << ntohs(pHeader.dcc) << RESTendl;
226#endif
227 RESTDebug << "Hdr word 0x" << ntohs(pHeader.hdr) << RESTendl;
228 RESTDebug << "Args 0x" << ntohs(pHeader.args) << RESTendl;
229 RESTDebug << "TS_H 0x" << ntohs(pHeader.ts_h) << RESTendl;
230 RESTDebug << "TS_L 0x" << ntohs(pHeader.ts_l) << RESTendl;
231 RESTDebug << "Ecnt 0x" << ntohs(pHeader.ecnt) << RESTendl;
232 RESTDebug << "Scnt 0x" << ntohs(pHeader.scnt) << std::dec << RESTendl;
233
234#ifdef NEW_DAQ_T2K_2_X
235 RESTDebug << "RawDCC Head 0x" << std::hex << ntohs(pHeader.dcc) << std::dec << " Version "
236 << GET_EVENT_TYPE(ntohs(pHeader.dcc));
237 RESTDebug << " Flag " << ((ntohs(pHeader.dcc) & 0x3000) >> 12);
238 RESTDebug << " RT " << ((ntohs(pHeader.dcc) & 0x0C00) >> 10) << " DCCInd "
239 << ((ntohs(pHeader.dcc) & 0x03F0) >> 4);
240 RESTDebug << " FEMInd " << (ntohs(pHeader.dcc) & 0x000F) << RESTendl;
241
242 RESTDebug << "FEM0Ind " << ntohs(pHeader.hdr) << " Type " << ((ntohs(pHeader.hdr) & 0xF000) >> 12);
243 RESTDebug << " L " << ((ntohs(pHeader.hdr) & 0x0800) >> 11);
244 RESTDebug << " U " << ((ntohs(pHeader.hdr) & 0x0800) >> 10) << " FECFlags "
245 << ((ntohs(pHeader.hdr) & 0x03F0) >> 4);
246 RESTDebug << " Index " << (ntohs(pHeader.hdr) & 0x000F) << RESTendl;
247
248 RESTDebug << "RawFEM 0x" << std::hex << ntohs(pHeader.args) << std::dec << " M "
249 << ((ntohs(pHeader.args) & 0x8000) >> 15);
250 RESTDebug << " N " << ((ntohs(pHeader.args) & 0x4000) >> 14) << " Zero "
251 << ((ntohs(pHeader.args) & 0x1000) >> 13);
252 RESTDebug << " Arg2 " << GET_RB_ARG2(ntohs(pHeader.args)) << " Arg2 "
253 << GET_RB_ARG1(ntohs(pHeader.args)) << RESTendl;
254 RESTDebug << "TimeStampH " << ntohs(pHeader.ts_h) << RESTendl;
255 RESTDebug << "TimeStampL " << ntohs(pHeader.ts_l) << RESTendl;
256 RESTDebug << "RawEvType 0x" << std::hex << ntohs(pHeader.ecnt) << std::dec << " EvTy "
257 << GET_EVENT_TYPE(ntohs(pHeader.ecnt));
258 RESTDebug << " EventCount " << GET_EVENT_COUNT(ntohs(pHeader.ecnt)) << RESTendl;
259 RESTDebug << "Samples " << ntohs(pHeader.scnt) << RESTendl;
260#endif
261
262 tempAsic1 = GET_RB_ARG1(ntohs(pHeader.args));
263 tempAsic2 = GET_RB_ARG2(ntohs(pHeader.args));
264 channel = tempAsic1 / 6;
265 asicN = (10 * (tempAsic1 % 6) / 2 + tempAsic2) % 4;
266 fecN = (10 * (tempAsic1 % 6) / 2 + tempAsic2) / 4;
267
268 RESTDebug << " channel " << channel << " asic " << asicN << " fec " << fecN << RESTendl;
269
270 sampleCountRead = ntohs(pHeader.scnt);
271 pay = sampleCountRead % 2;
272
273 physChannel = -10;
274 if (channel > 2 && channel < 15) {
275 physChannel = channel - 3;
276 } else if (channel > 15 && channel < 28) {
277 physChannel = channel - 4;
278 } else if (channel > 28 && channel < 53) {
279 physChannel = channel - 5;
280 } else if (channel > 53 && channel < 66) {
281 physChannel = channel - 6;
282 } else if (channel > 66) {
283 physChannel = channel - 7;
284 }
285
286 if (physChannel >= 0) // isThisAphysChannel?
287 // in this case we'd hold it in hEvent
288 // but we need to redefine the physChannel number
289 // so as not to overwrite the information
290 // and to know which ASIC and FEC it belongs to.
291 {
292 isData = true;
293 physChannel = fecN * 72 * 4 + asicN * 72 + physChannel;
294
295 } else
296 isData = false;
297
298 timeBin = 0;
299
300 if (sampleCountRead < 9) isData = false;
301 for (int i = 0; i < sampleCountRead; i++) {
302 int y = fread(&dat, sizeof(uint16_t), 1, fInputBinFile);
303 if (y == 0)
304 RESTError << "TRestRawAFTERToSignalProcess::ProcessEvent. Problems reading input file."
305 << RESTendl;
306 frameBits += sizeof(dat);
307 data = ntohs(dat);
308
309 std::bitset<16> bs(data);
310 RESTDebug << bs << RESTendl;
311
312 if (((data & 0xFE00) >> 9) == 8) {
313 timeBin = GET_CELL_INDEX(data);
314 if (timeBin == 511) isData = false;
315 RESTDebug << data << " Time bin " << timeBin << RESTendl;
316 } else if ((((data & 0xF000) >> 12) == 0) && isData) {
317 fSignalEvent->AddChargeToSignal(physChannel, timeBin, data);
318 RESTDebug << "Time bin " << timeBin << " ADC: " << data << RESTendl;
319 timeBin++;
320 }
321 }
322
323 RESTDebug << pay << RESTendl;
324 if (pay) {
325 int z = fread(&dat, sizeof(uint16_t), 1, fInputBinFile);
326 if (z == 0)
327 RESTError << "TRestRawAFTERToSignalProcess::ProcessEvent. Problems reading input file."
328 << RESTendl;
329 frameBits += sizeof(uint16_t);
330 }
331
332 int w = fread(&pEnd, sizeof(DataPacketEnd), 1, fInputBinFile);
333 if (w == 0)
334 RESTError << "TRestRawAFTERToSignalProcess::ProcessEvent. Problems reading input file."
335 << RESTendl;
336 frameBits += sizeof(DataPacketEnd);
337
338 RESTDebug << "Read "
339 << sampleCountRead * sizeof(uint16_t) + sizeof(DataPacketHeader) + sizeof(DataPacketEnd) +
340 sampleCountRead % 2 * sizeof(uint16_t)
341 << " vs HeadSize " << ntohs(pHeader.size) << " Diff "
342 << ntohs(pHeader.size) - (sampleCountRead + sizeof(DataPacketHeader) +
343 sizeof(DataPacketEnd) + sampleCountRead % 2)
344 << RESTendl;
345 RESTDebug << "Trailer_H " << ntohs(pEnd.crc1) << " Trailer_L " << ntohs(pEnd.crc2) << RESTendl;
346 RESTDebug << "Trailer " << eventTime << "\n" << RESTendl;
347
348 } // end while
349 totalbytesRead += frameBits;
350
351 // printf("Event ID %d time stored
352 // %.3lf\n",fSignalEvent->GetID(),fSignalEvent->GetTime());
353
354 RESTDebug << "End of event " << RESTendl;
355
356 return fSignalEvent;
357}
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.
A process to read binary files produced with AFTER electronics.
~TRestRawAFTERToSignalProcess()
Constructor loading data from a config file.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
void InitProcess() override
Process initialization.
TRestRawAFTERToSignalProcess()
Default constructor.
void Initialize() override
Function to initialize input/output event members and define the section name.
void AddChargeToSignal(Int_t sgnlID, Int_t bin, Short_t value)
void Initialize() override
Making default settings.