DUNE-DAQ
DUNE Trigger and Data Acquisition software
Loading...
Searching...
No Matches
utils.py
Go to the documentation of this file.
1#general imports
2import sys
3
4#dunedaq imports
5import daqdataformats
6import detdataformats
7import fddetdataformats
8import trgdataformats
9import detchannelmaps
10
11#unpacker imports
16import h5py
17
18#analysis imports
19import numpy as np
20import numpy.fft
21
23
24 def get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, ch):
25 return (ch + (stream << 10 ) + (slot << 18 ) + (crate << 22))
27 return -1
29 return 0xffff
31 return "Null"
32
34
35 is_fragment_unpacker = False
36
37 is_detector_unpacker = False
38 is_trigger_unpacker = False
39
40 def __init__(self, index=None):
41 self.index = index
42
43 def get_all_data(self,in_data=None):
44 return None
45
47
48 is_fragment_unpacker = False
49
50 is_detector_unpacker = False
51 is_trigger_unpacker = False
52
53 def get_srcid_data(self,sid):
54 return [ SourceIDData(run=self.index.run,
55 trigger=self.index.trigger,
56 sequence=self.index.sequence,
57 src_id=sid.id,
58 subsystem=sid.subsystem,
59 subsystem_str=daqdataformats.SourceID.subsystem_to_string(sid.subsystem),
60 version=sid.version) ]
61
62 def get_all_data(self,in_data):
63 #in_data = sid
64 return { "sid": self.get_srcid_data(in_data) }
65
67
68 is_fragment_unpacker = False
69
70 def get_trh_data(self,trh,n_fragments):
71 return [ TriggerRecordData(run=trh.get_run_number(),
72 trigger=trh.get_trigger_number(),
73 sequence=trh.get_sequence_number(),
74 trigger_timestamp_dts=trh.get_trigger_timestamp(),
75 n_fragments=n_fragments,
76 n_requested_components=trh.get_num_requested_components(),
77 error_bits=trh.get_header().error_bits,
78 trigger_type=trh.get_trigger_type(),
79 max_sequence_number=trh.get_max_sequence_number(),
80 total_size_bytes=trh.get_total_size_bytes()) ]
81
82 def get_all_data(self,in_data):
83 #in_data[0]=trh
84 #in_data[1]=n_fragments
85 return { "trh": self.get_trh_data(in_data[0],in_data[1]) }
86
88
89 is_fragment_unpacker = True
90
91 is_detector_unpacker = False
92 is_trigger_unpacker = False
93
94 def get_n_obj(self,frag):
95 return None
96
97 def get_trg_data(self,in_data):
98 return None, None
99
100 def get_det_data(self,in_data):
101 return None, None, None, None
102
103 def get_frh_data(self,frag):
104 frh = frag.get_header()
105 return [ FragmentHeaderData(run=frh.run_number,
106 trigger=frh.trigger_number,
107 sequence=frh.sequence_number,
108 src_id=frh.element_id.id,
109 trigger_timestamp_dts=frh.trigger_timestamp,
110 window_begin_dts=frh.window_begin,
111 window_end_dts=frh.window_end,
112 det_id=frh.detector_id,
113 error_bits=frh.error_bits,
114 fragment_type=frh.fragment_type,
115 total_size_bytes=frh.size,
116 data_size_bytes=frag.get_data_size()) ]
117
118 def get_all_data(self,in_data):
119 #in_data = fragment
120
121 data_dict = { "frh": self.get_frh_data(in_data) }
122
123 #if no data, nothing to unpack further
124 if in_data.get_data_size()==0:
125 return data_dict
126
127 type_string = f'{detdataformats.DetID.Subdetector(in_data.get_detector_id()).name}_{in_data.get_fragment_type().name}'
128
130 trgh, trgd = self.get_trg_data(in_data)
131 if trgh is not None: data_dict[f"trgh_{type_string}"] = trgh
132 if trgd[0] is not None: data_dict[f"trgd_{type_string}"] = trgd[0]
133 if trgd[1] is not None: data_dict[f"trgd_{type_string}_inputs"] = trgd[1]
134
136 daqh, deth, detd, detw = self.get_det_data(in_data)
137 if daqh is not None: data_dict["daqh"] = daqh
138 if deth is not None: data_dict[f"deth_{type_string}"] = deth
139 if detd is not None: data_dict[f"detd_{type_string}"] = detd
140 if detw is not None: data_dict[f"detw_{type_string}"] = detw
141
142 return data_dict
143
145
146 is_trigger_unpacker = True
147
148 def get_trg_data(self,frag):
149 return self.get_trg_header_data(frag),self.get_trg_obj_data(frag)
150
151 def get_trg_header_data(self,frag):
152 frh = frag.get_header()
153 return [ TriggerHeaderData(run=frh.run_number,
154 trigger=frh.trigger_number,
155 sequence=frh.sequence_number,
156 src_id=frh.element_id.id,
157 n_obj=self.get_n_obj(frag),
158 version=self.get_trg_data_version(frag)) ]
159
160
162
163 trg_obj = trgdataformats.TriggerPrimitive
164
165 def __init__(self,channel_map=None):
166 super().__init__()
167 if 'TPC' in channel_map:
168 self.channel_map = detchannelmaps.make_tpc_map(channel_map)
169 elif 'PDS' in channel_map:
170 self.channel_map = detchannelmaps.make_pds_map(channel_map)
171 else:
172 self.channel_map = NullChannelMap
173
174 def get_n_obj(self,frag):
175 return int(frag.get_data_size()/self.trg_obj.sizeof())
176
177 def get_trg_data_version(self,frag):
178 if self.get_n_objget_n_obj(frag)==0:
179 return None
180 trig_obj = self.trg_obj(frag.get_data())
181 return trig_obj.version
182
183 def get_trg_obj_data(self,frag):
184 frh = frag.get_header()
185 tpd_list = []
186 for i_tp in range(self.get_n_objget_n_obj(frag)):
187 tp = self.trg_obj(frag.get_data(i_tp*self.trg_obj.sizeof()))
188 ch_info = self.channel_map.get_channel_info_from_offline_channel(tp.channel)
189 tpd_list.append( TriggerPrimitiveData(run=frh.run_number,
190 trigger=frh.trigger_number,
191 sequence=frh.sequence_number,
192 src_id=frh.element_id.id,
193 time_start=tp.time_start,
194 samples_to_peak=tp.samples_to_peak,
195 samples_over_threshold=tp.samples_over_threshold,
196 channel=tp.channel,
197 plane=self.channel_map.get_plane_from_offline_channel(tp.channel),
198 element=ch_info.element,
199 adc_integral=tp.adc_integral,
200 adc_peak=tp.adc_peak,
201 detid=tp.detid,
202 flag=tp.flag,
203 id_ta=-1) )
204 return tpd_list, None
205
207
208 trg_obj = trgdataformats.TriggerActivity
209
210 def __init__(self,channel_map=None):
211 super().__init__()
212 if 'TPC' in channel_map:
213 self.channel_map = detchannelmaps.make_tpc_map(channel_map)
214 elif 'PDS' in channel_map:
215 self.channel_map = detchannelmaps.make_pds_map(channel_map)
216 else:
217 self.channel_map = NullChannelMap
218
219 def get_n_obj(self,frag):
220 frag_data_size = frag.get_data_size()
221 size_so_far = 0
222 n_tobjs = 0
223 while size_so_far < frag_data_size:
224 tobj = self.trg_obj(frag.get_data(size_so_far))
225 n_tobjs = n_tobjs + 1
226 size_so_far = size_so_far + tobj.sizeof()
227 return n_tobjs
228
229 def get_trg_data_version(self,frag):
230 if self.get_n_objget_n_obj(frag)==0:
231 return None
232 trig_obj = self.trg_obj(frag.get_data())
233 return trig_obj.data.version
234
235 def get_trg_obj_data(self,frag):
236 frh = frag.get_header()
237 ta_list = []
238 tpd_list = []
239 size_so_far = 0
240 for i_ta in range(self.get_n_objget_n_obj(frag)):
241 ta= self.trg_obj(frag.get_data(size_so_far))
242 size_so_far = size_so_far + ta.sizeof()
243 ch_info_ta = self.channel_map.get_channel_info_from_offline_channel(ta.data.channel_peak)
244 ta_list.append( TriggerActivityData(run=frh.run_number,
245 trigger=frh.trigger_number,
246 sequence=frh.sequence_number,
247 src_id=frh.element_id.id,
248 id=i_ta,
249 time_start=ta.data.time_start,
250 time_end=ta.data.time_end,
251 time_peak=ta.data.time_peak,
252 time_activity=ta.data.time_activity,
253 channel_start=ta.data.channel_start,
254 channel_end=ta.data.channel_end,
255 channel_peak=ta.data.channel_peak,
256 plane=self.channel_map.get_plane_from_offline_channel(ta.data.channel_peak),
257 element=ch_info_ta.element,
258 adc_integral=ta.data.adc_integral,
259 adc_peak=ta.data.adc_peak,
260 detid=ta.data.detid,
261 ta_type=ta.data.type,
262 algorithm=ta.data.algorithm,
263 n_tps=len(ta),
264 id_tc=-1) )
265 for i_tp in range(len(ta)):
266 tp = ta[i_tp]
267 ch_info_tp = self.channel_map.get_channel_info_from_offline_channel(tp.channel)
268 tpd_list.append( TriggerPrimitiveData(run=frh.run_number,
269 trigger=frh.trigger_number,
270 sequence=frh.sequence_number,
271 src_id=frh.element_id.id,
272 time_start=tp.time_start,
273 samples_to_peak=tp.samples_to_peak,
274 samples_over_threshold=tp.samples_over_threshold,
275 channel=tp.channel,
276 plane=self.channel_map.get_plane_from_offline_channel(tp.channel),
277 element=ch_info_tp.element,
278 adc_integral=tp.adc_integral,
279 adc_peak=tp.adc_peak,
280 detid=tp.detid,
281 flag=tp.flag,
282 id_ta=i_ta) )
283
284 if len(tpd_list)==0:
285 tpd_list=None
286 return ta_list, tpd_list
287
289
290 trg_obj = trgdataformats.TriggerCandidate
291
292 def __init__(self,channel_map=None):
293 super().__init__()
294 if not channel_map:
295 self.channel_map = NullChannelMap
296 elif 'TPC' in channel_map:
297 self.channel_map = detchannelmaps.make_tpc_map(channel_map)
298 elif 'PDS' in channel_map:
299 self.channel_map = detchannelmaps.make_pds_map(channel_map)
300 else:
301 self.channel_map = NullChannelMap
302
303 def get_n_obj(self,frag):
304 frag_data_size = frag.get_data_size()
305 size_so_far = 0
306 n_tobjs = 0
307 while size_so_far < frag_data_size:
308 tobj = self.trg_obj(frag.get_data(size_so_far))
309 n_tobjs = n_tobjs + 1
310 size_so_far = size_so_far + tobj.sizeof()
311 return n_tobjs
312
313 def get_trg_data_version(self,frag):
314 if self.get_n_objget_n_obj(frag)==0:
315 return None
316 trig_obj = self.trg_obj(frag.get_data())
317 return trig_obj.data.version
318
319 def get_trg_obj_data(self,frag):
320 frh = frag.get_header()
321 tc_list = []
322 ta_list = []
323 size_so_far = 0
324 for i_tc in range(self.get_n_objget_n_obj(frag)):
325 tc = self.trg_obj(frag.get_data(size_so_far))
326 size_so_far = size_so_far + tc.sizeof()
327 tc_list.append( TriggerCandidateData(run=frh.run_number,
328 trigger=frh.trigger_number,
329 sequence=frh.sequence_number,
330 src_id=frh.element_id.id,
331 id=i_tc,
332 time_start=tc.data.time_start,
333 time_end=tc.data.time_end,
334 time_candidate=tc.data.time_candidate,
335 detid=tc.data.detid,
336 tc_type=tc.data.type,
337 algorithm=tc.data.algorithm,
338 n_tas=len(tc) ) )
339 for i_ta in range(len(tc)):
340 ta = tc[i_ta]
341 ch_info_ta = self.channel_map.get_channel_info_from_offline_channel(ta.channel_peak)
342 ta_list.append( TriggerActivityData(run=frh.run_number,
343 trigger=frh.trigger_number,
344 sequence=frh.sequence_number,
345 src_id=frh.element_id.id,
346 id=i_ta,
347 time_start=ta.time_start,
348 time_end=ta.time_end,
349 time_peak=ta.time_peak,
350 time_activity=ta.time_activity,
351 channel_start=ta.channel_start,
352 channel_end=ta.channel_end,
353 channel_peak=ta.channel_peak,
354 plane=self.channel_map.get_plane_from_offline_channel(ta.channel_peak),
355 element=ch_info_ta.element,
356 adc_integral=ta.adc_integral,
357 adc_peak=ta.adc_peak,
358 detid=ta.detid,
359 ta_type=ta.type,
360 algorithm=ta.algorithm,
361 n_tps=-1,
362 id_tc=i_tc) )
363 if len(ta_list)==0:
364 ta_list=None
365 return tc_list, ta_list
366
367
369
370 is_detector_unpacker = True
371
372 def __init__(self,channel_map=None,ana_data_prescale=1,wvfm_data_prescale=None):
373 super().__init__()
374 self.ana_data_prescale = None if not ana_data_prescale else int(ana_data_prescale)
375 self.wvfm_data_prescale = None if not wvfm_data_prescale else int(wvfm_data_prescale)
376 if not channel_map:
377 self.channel_map = NullChannelMap
378 elif 'TPC' in channel_map:
379 self.channel_map = detchannelmaps.make_tpc_map(channel_map)
380 elif 'PDS' in channel_map:
381 self.channel_map = detchannelmaps.make_pds_map(channel_map)
382 else:
383 self.channel_map = NullChannelMap
384
386 return None
387
388 def get_det_data_version(self,frag):
389 return None
390
391 def get_timestamp_first(self,frag):
392 return None
393
395 return None, None, None, None
396
397 def get_daq_header_data(self,frag):
398 frh = frag.get_header()
399 det_id, crate_id, slot_id, stream_id = self.get_det_crate_slot_stream(frag)
400 return [ DAQHeaderData(run=frh.run_number,
401 trigger=frh.trigger_number,
402 sequence=frh.sequence_number,
403 src_id=frh.element_id.id,
404 n_obj=self.get_n_obj(frag),
405 daq_header_version=self.get_daq_header_version(frag),
406 det_data_version=self.get_det_data_version(frag),
407 det_id=det_id,
408 crate_id=crate_id,
409 slot_id=slot_id,
410 stream_id=stream_id,
411 timestamp_first_dts=self.get_timestamp_first(frag)) ]
412
413 def get_det_header_data(self,frag):
414 return None
415
416 def get_det_data_all(self,frag):
417 return None, None
418
419 def get_det_data(self,frag):
420 det_ana_data, det_wvfm_data = self.get_det_data_all(frag)
421 return self.get_daq_header_data(frag), self.get_det_header_data(frag), det_ana_data, det_wvfm_data
422
423
425
427 frame_obj = fddetdataformats.WIBEthFrame
428
429 SAMPLING_PERIOD = 32
430 N_CHANNELS_PER_FRAME = 64
431
432 def get_n_obj(self,frag):
433 return self.unpacker.get_n_frames(frag)
434
436 return self.frame_obj(frag.get_data()).get_daqheader().version
437
438 def get_timestamp_first(self,frag):
439 return self.frame_obj(frag.get_data()).get_timestamp()
440
441 def get_det_data_version(self,frag):
442 return self.frame_obj(frag.get_data()).get_wibheader().version
443
445 dh = self.frame_obj(frag.get_data()).get_daqheader()
446 return dh.det_id, dh.crate_id, dh.slot_id, dh.stream_id
447
448 def get_det_header_data(self,frag):
449 frh = frag.get_header()
450
451 n_frames = self.get_n_objget_n_obj(frag)
452
453 pulser_arr = np.empty(n_frames)
454 calibration_arr = np.empty(n_frames)
455 ready_arr = np.empty(n_frames)
456 context_arr = np.empty(n_frames)
457 wib_sync_arr = np.empty(n_frames)
458 femb_sync_arr = np.empty(n_frames)
459 cd_arr = np.empty(n_frames)
460 crc_err_arr = np.empty(n_frames)
461 link_valid_arr = np.empty(n_frames)
462 lol_arr = np.empty(n_frames)
463 colddata_ts0_arr = np.empty(n_frames)
464 colddata_ts1_arr = np.empty(n_frames)
465
466 for i in range(n_frames):
467 wh = self.frame_obj(frag.get_data(i*self.frame_obj.sizeof())).get_wibheader()
468
469 pulser_arr[i] = wh.pulser
470 calibration_arr[i] = wh.calibration
471 ready_arr[i] = wh.ready
472 context_arr[i] = wh.context
473
474 wib_sync_arr[i] = wh.wib_sync
475 femb_sync_arr[i] = wh.femb_sync
476
477 cd_arr[i] = wh.cd
478 crc_err_arr[i] = wh.crc_err
479 link_valid_arr[i] = wh.link_valid
480 lol_arr[i] = wh.lol
481
482 colddata_ts0_arr[i] = wh.colddata_timestamp_0
483 colddata_ts1_arr[i] = wh.colddata_timestamp_1
484
485 pulser_change_idx, pulser_change_val, _ = sparsify_array_diff_locs_and_vals(pulser_arr)
486 calibration_change_idx, calibration_change_val, _ = sparsify_array_diff_locs_and_vals(calibration_arr)
487 ready_change_idx, ready_change_val, _ = sparsify_array_diff_locs_and_vals(context_arr)
488 context_change_idx, context_change_val, _ = sparsify_array_diff_locs_and_vals(context_arr)
489
490 wib_sync_change_idx, wib_sync_change_val, _ = sparsify_array_diff_locs_and_vals(wib_sync_arr)
491 femb_sync_change_idx, femb_sync_change_val, _ = sparsify_array_diff_locs_and_vals(femb_sync_arr)
492
493 cd_change_idx, cd_change_val, _ = sparsify_array_diff_locs_and_vals(cd_arr)
494 crc_err_change_idx, crc_err_change_val, _ = sparsify_array_diff_locs_and_vals(crc_err_arr)
495 link_valid_change_idx, link_valid_change_val, _ = sparsify_array_diff_locs_and_vals(link_valid_arr)
496 lol_change_idx, lol_change_val, _ = sparsify_array_diff_locs_and_vals(context_arr)
497
498 colddata_ts0_diff = np.diff(colddata_ts0_arr)
499 colddata_ts0_diff[colddata_ts0_diff<0] = colddata_ts0_diff[colddata_ts0_diff<0]+0x8000
500 colddata_ts0_diff_change_idx, colddata_ts0_diff_change_val, _ = sparsify_array_diff_locs_and_vals(colddata_ts0_diff)
501
502 colddata_ts1_diff = np.diff(colddata_ts1_arr)
503 colddata_ts1_diff[colddata_ts1_diff<0] = colddata_ts1_diff[colddata_ts1_diff<0]+0x8000
504 colddata_ts1_diff_change_idx, colddata_ts1_diff_change_val, _ = sparsify_array_diff_locs_and_vals(colddata_ts1_diff)
505
506 ts_arr = self.unpacker.np_array_timestamp(frag)
507 ts_diff_change_idx, ts_diff_change_val, _ = sparsify_array_diff_locs_and_vals(np.diff(ts_arr))
508
509 wh = self.frame_obj(frag.get_data()).get_wibheader()
510 ts_diff_vals, ts_diff_counts = np.unique(np.diff(self.unpacker.np_array_timestamp(frag)),return_counts=True)
511 return [ WIBEthHeaderData(run=frh.run_number,
512 trigger=frh.trigger_number,
513 sequence=frh.sequence_number,
514 src_id=frh.element_id.id,
515 femb_id=(wh.channel>>1)&0x3,
516 colddata_id=wh.channel&0x1,
517 version=wh.version,
518 pulser_vals=pulser_change_val, pulser_idx=pulser_change_idx,
519 calibration_vals=calibration_change_val, calibration_idx=calibration_change_idx,
520 ready_vals=ready_change_val, ready_idx=ready_change_idx,
521 context_vals=context_change_val, context_idx=context_change_idx,
522 wib_sync_vals=wib_sync_change_val, wib_sync_idx=wib_sync_change_idx,
523 femb_sync_vals=femb_sync_change_val, femb_sync_idx=femb_sync_change_idx,
524 cd_vals=cd_change_val, cd_idx=cd_change_idx,
525 crc_err_vals=crc_err_change_val, crc_err_idx=crc_err_change_idx,
526 link_valid_vals=link_valid_change_val, link_valid_idx=link_valid_change_idx,
527 lol_vals=lol_change_val, lol_idx=lol_change_idx,
528 colddata_timestamp_0_diff_vals=colddata_ts0_diff_change_val,
529 colddata_timestamp_0_diff_idx=colddata_ts0_diff_change_idx,
530 colddata_timestamp_0_first=colddata_ts0_arr[0],
531 colddata_timestamp_1_diff_vals=colddata_ts1_diff_change_val,
532 colddata_timestamp_1_diff_idx=colddata_ts1_diff_change_idx,
533 colddata_timestamp_1_first=colddata_ts1_arr[0],
534 timestamp_dts_diff_vals=ts_diff_change_val, timestamp_dts_diff_idx=ts_diff_change_idx,
535 timestamp_dts_first=ts_arr[0],
536 n_frames=n_frames,
537 n_channels=self.N_CHANNELS_PER_FRAME,
538 sampling_period=self.SAMPLING_PERIOD) ]
539
540 def get_det_data_all(self,frag):
541 frh = frag.get_header()
542 trigger_number = frh.trigger_number
543
544 get_ana_data = (self.ana_data_prescale is not None and (trigger_number % self.ana_data_prescale)==0)
545 get_wvfm_data = (self.wvfm_data_prescale is not None and (trigger_number % self.wvfm_data_prescale)==0)
546
547 #print(f'\t\tTrigger number {trigger_number}: get_ana_data? {get_ana_data} get_wvfm_data? {get_wvfm_data}')
548
549 if not (get_ana_data or get_wvfm_data):
550 return None,None
551
552 ana_data = None
553 wvfm_data = None
554
555 adcs = self.unpacker.np_array_adc(frag)
556 det, crate, slot, stream = self.get_det_crate_slot_streamget_det_crate_slot_stream(frag)
557 channels = [ self.channel_map.get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, c) for c in range(self.N_CHANNELS_PER_FRAME) ]
558 planes = [ self.channel_map.get_plane_from_offline_channel(uc) for uc in channels ]
559 elements = [ self.channel_map.get_element_id_from_offline_channel(uc) for uc in channels ]
560 wib_chans = range(self.N_CHANNELS_PER_FRAME)
561
562 if get_ana_data:
563 adc_mean = np.mean(adcs,axis=0)
564 adc_rms = np.std(adcs,axis=0)
565 adc_max = np.max(adcs,axis=0)
566 adc_min = np.min(adcs,axis=0)
567 adc_median = np.median(adcs,axis=0)
568 ana_data = [ WIBEthAnalysisData(run=frh.run_number,
569 trigger=frh.trigger_number,
570 sequence=frh.sequence_number,
571 src_id=frh.element_id.id,
572 channel=channels[i_ch],
573 plane=planes[i_ch],
574 element=elements[i_ch],
575 wib_chan=wib_chans[i_ch],
576 adc_mean=adc_mean[i_ch],
577 adc_rms=adc_rms[i_ch],
578 adc_max=adc_max[i_ch],
579 adc_min=adc_min[i_ch],
580 adc_median=adc_median[i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
581 if get_wvfm_data:
582 timestamps = self.unpacker.np_array_timestamp(frag)
583 ffts = np.abs(np.fft.rfft(adcs,axis=0))
584 wvfm_data = [ WIBEthWaveformData(run=frh.run_number,
585 trigger=frh.trigger_number,
586 sequence=frh.sequence_number,
587 src_id=frh.element_id.id,
588 channel=channels[i_ch],
589 plane=planes[i_ch],
590 element=elements[i_ch],
591 wib_chan=wib_chans[i_ch],
592 timestamps=timestamps,
593 adcs=adcs[:,i_ch],
594 fft_mag=ffts[:,i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
595
596 return ana_data, wvfm_data
597
598
600
602 frame_obj = fddetdataformats.TDEEthFrame
603
604 SAMPLING_PERIOD = 31.25
605 N_CHANNELS_PER_FRAME = 64
606
607 def get_n_obj(self,frag):
608 return self.unpacker.get_n_frames(frag)
609
611 return self.frame_obj(frag.get_data()).get_daqheader().version
612
613 def get_timestamp_first(self,frag):
614 return self.frame_obj(frag.get_data()).get_timestamp()
615
616 def get_det_data_version(self,frag):
617 return self.frame_obj(frag.get_data()).get_tdeheader().version
618
620 dh = self.frame_obj(frag.get_data()).get_daqheader()
621 return dh.det_id, dh.crate_id, dh.slot_id, dh.stream_id
622
623 def get_det_header_data(self,frag):
624 frh = frag.get_header()
625
626 n_frames = self.get_n_objget_n_obj(frag)
627
628 errors_arr = np.empty(n_frames)
629 tai_time_arr = np.empty(n_frames)
630
631 for i in range(n_frames):
632 tdeh = self.frame_obj(frag.get_data(i*self.frame_obj.sizeof())).get_tdeheader()
633
634 errors_arr[i] = tdeh.tde_errors
635 tai_time_arr[i] = tdeh.TAItime
636
637 errors_change_idx, errors_change_val, _ = sparsify_array_diff_locs_and_vals(errors_arr)
638
639 tai_time_diff = np.diff(tai_time_arr)
640 tai_time_diff_change_idx, tai_time_diff_change_val, _ = sparsify_array_diff_locs_and_vals(tai_time_diff)
641
642 ts_arr = self.unpacker.np_array_timestamp(frag)
643 ts_diff_change_idx, ts_diff_change_val, _ = sparsify_array_diff_locs_and_vals(np.diff(ts_arr))
644
645 tdeh = self.frame_obj(frag.get_data()).get_tdeheader()
646 ts_diff_vals, ts_diff_counts = np.unique(np.diff(self.unpacker.np_array_timestamp(frag)),return_counts=True)
647 return [ TDEEthHeaderData(run=frh.run_number,
648 trigger=frh.trigger_number,
649 sequence=frh.sequence_number,
650 src_id=frh.element_id.id,
651 channel_id=tdeh.channel,
652 tde_header=tdeh.tde_header,
653 version=tdeh.version,
654 errors_vals=errors_change_val, errors_idx=errors_change_idx,
655 tai_time_diff_vals=tai_time_diff_change_val,
656 tai_time_diff_idx=tai_time_diff_change_idx,
657 tai_time_first=tai_time_arr[0],
658 timestamp_dts_diff_vals=ts_diff_change_val,
659 timestamp_dts_diff_idx=ts_diff_change_idx,
660 timestamp_dts_first=ts_arr[0],
661 n_frames=n_frames,
662 n_channels=self.N_CHANNELS_PER_FRAME,
663 sampling_period=self.SAMPLING_PERIOD) ]
664
665 def get_det_data_all(self,frag):
666 frh = frag.get_header()
667 trigger_number = frh.trigger_number
668
669 get_ana_data = (self.ana_data_prescale is not None and (trigger_number % self.ana_data_prescale)==0)
670 get_wvfm_data = (self.wvfm_data_prescale is not None and (trigger_number % self.wvfm_data_prescale)==0)
671
672 #print(f'\t\tTrigger number {trigger_number}: get_ana_data? {get_ana_data} get_wvfm_data? {get_wvfm_data}')
673
674 if not (get_ana_data or get_wvfm_data):
675 return None,None
676
677 ana_data = None
678 wvfm_data = None
679
680 adcs = self.unpacker.np_array_adc(frag)
681 det, crate, slot, stream = self.get_det_crate_slot_streamget_det_crate_slot_stream(frag)
682 channels = [ self.channel_map.get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, c) for c in range(self.N_CHANNELS_PER_FRAME) ]
683 planes = [ self.channel_map.get_plane_from_offline_channel(uc) for uc in channels ]
684 elements = [ self.channel_map.get_element_id_from_offline_channel(uc) for uc in channels ]
685 tde_chans = range(self.N_CHANNELS_PER_FRAME)
686
687 if get_ana_data:
688 adc_mean = np.mean(adcs,axis=0)
689 adc_rms = np.std(adcs,axis=0)
690 adc_max = np.max(adcs,axis=0)
691 adc_min = np.min(adcs,axis=0)
692 adc_median = np.median(adcs,axis=0)
693 ana_data = [ TDEEthAnalysisData(run=frh.run_number,
694 trigger=frh.trigger_number,
695 sequence=frh.sequence_number,
696 src_id=frh.element_id.id,
697 channel=channels[i_ch],
698 plane=planes[i_ch],
699 element=elements[i_ch],
700 tde_chan=tde_chans[i_ch],
701 adc_mean=adc_mean[i_ch],
702 adc_rms=adc_rms[i_ch],
703 adc_max=adc_max[i_ch],
704 adc_min=adc_min[i_ch],
705 adc_median=adc_median[i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
706 if get_wvfm_data:
707 timestamps = self.unpacker.np_array_timestamp(frag)
708 ffts = np.abs(np.fft.rfft(adcs,axis=0))
709 wvfm_data = [ TDEEthWaveformData(run=frh.run_number,
710 trigger=frh.trigger_number,
711 sequence=frh.sequence_number,
712 src_id=frh.element_id.id,
713 channel=channels[i_ch],
714 plane=planes[i_ch],
715 element=elements[i_ch],
716 tde_chan=tde_chans[i_ch],
717 timestamps=timestamps,
718 adcs=adcs[:,i_ch],
719 fft_mag=ffts[:,i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
720
721 return ana_data, wvfm_data
722
724
726 frame_obj = fddetdataformats.DAPHNEStreamFrame
727
728 SAMPLING_PERIOD = 1
729 N_CHANNELS_PER_FRAME = 4
730
731 def get_n_obj(self,frag):
732 return self.unpacker.get_n_frames_stream(frag)
733
735 return self.frame_obj(frag.get_data()).get_daqheader().version
736
737 def get_timestamp_first(self,frag):
738 return self.frame_obj(frag.get_data()).get_timestamp()
739
740 def get_det_data_version(self,frag):
741 return 0
742
744 dh = self.frame_obj(frag.get_data()).get_daqheader()
745 return dh.det_id, dh.crate_id, dh.slot_id, dh.link_id
746
747 def get_det_header_data(self,frag):
748 frh = frag.get_header()
749 #dh = self.frame_obj(frag.get_data()).get_header()
750 ts_diffs_vals, ts_diffs_counts = np.unique(np.diff( np.array(self.unpacker.np_array_timestamp_stream(frag), dtype=np.int64)), return_counts=True)
751 return [ DAPHNEStreamHeaderData(run=frh.run_number,
752 trigger=frh.trigger_number,
753 sequence=frh.sequence_number,
754 src_id=frh.element_id.id,
755 n_channels=self.N_CHANNELS_PER_FRAME,
756 sampling_period=self.SAMPLING_PERIOD,
757 ts_diffs_vals=ts_diffs_vals,
758 ts_diffs_counts=ts_diffs_counts) ]
759
760 def get_det_data_all(self,frag):
761 frh = frag.get_header()
762 trigger_number = frh.trigger_number
763
764 get_ana_data = (self.ana_data_prescale is not None and (trigger_number % self.ana_data_prescale)==0)
765 get_wvfm_data = (self.wvfm_data_prescale is not None and (trigger_number % self.wvfm_data_prescale)==0)
766
767 if not (get_ana_data or get_wvfm_data):
768 return None,None
769
770 ana_data = None
771 wvfm_data = None
772
773 adcs = self.unpacker.np_array_adc_stream(frag)
774 dh = self.frame_obj(frag.get_data()).get_header()
775 det, crate, slot, stream = self.get_det_crate_slot_streamget_det_crate_slot_stream(frag)
776 daphne_chans = [ dh.channel_0, dh.channel_1, dh.channel_2, dh.channel_3 ]
777 channels = [ self.channel_map.get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, c) for c in daphne_chans ]
778
779 if get_ana_data:
780 adc_mean = np.mean(adcs,axis=0)
781 adc_rms = np.std(adcs,axis=0)
782 adc_max = np.max(adcs,axis=0)
783 adc_min = np.min(adcs,axis=0)
784 adc_median = np.median(adcs,axis=0)
785 ana_data = [ DAPHNEStreamAnalysisData(run=frh.run_number,
786 trigger=frh.trigger_number,
787 sequence=frh.sequence_number,
788 src_id=frh.element_id.id,
789 channel=channels[i_ch],
790 daphne_chan=daphne_chans[i_ch],
791 adc_mean=adc_mean[i_ch],
792 adc_rms=adc_rms[i_ch],
793 adc_max=adc_max[i_ch],
794 adc_min=adc_min[i_ch],
795 adc_median=adc_median[i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
796 if get_wvfm_data:
797 timestamps = self.unpacker.np_array_timestamp_stream(frag)
798 ffts = np.abs(np.fft.rfft(adcs,axis=0))
799 wvfm_data = [ DAPHNEStreamWaveformData(run=frh.run_number,
800 trigger=frh.trigger_number,
801 sequence=frh.sequence_number,
802 src_id=frh.element_id.id,
803 channel=channels[i_ch],
804 daphne_chan=channels[i_ch],
805 adcs=adcs[:,i_ch],
806 timestamps=timestamps,
807 fft_mag=ffts[:,i_ch]) for i_ch in range(self.N_CHANNELS_PER_FRAME) ]
808 return ana_data, wvfm_data
809
812 frame_obj = fddetdataformats.DAPHNEStreamFrame
813
815 dh = self.frame_objframe_obj(frag.get_data()).get_daqheader()
816 return dh.det_id, dh.crate_id, dh.slot_id, dh.stream_id
817
819
821 frame_obj = fddetdataformats.DAPHNEFrame
822
823 SAMPLING_PERIOD = 1
824 N_CHANNELS_PER_FRAME = 1
825
826 def get_n_obj(self,frag):
827 return self.unpacker.get_n_frames(frag)
828
830 return self.frame_obj(frag.get_data()).get_daqheader().version
831
832 def get_timestamp_first(self,frag):
833 return self.frame_obj(frag.get_data()).get_timestamp()
834
835 def get_det_data_version(self,frag):
836 return 0
837
839 dh = self.frame_obj(frag.get_data()).get_daqheader()
840 return dh.det_id, dh.crate_id, dh.slot_id, dh.link_id
841
842 def get_det_header_data(self,frag):
843 return None
844
845 def get_det_data_all(self,frag):
846 frh = frag.get_header()
847 trigger_number = frh.trigger_number
848 wvfm_data = None
849 ana_data = None
850
851 get_ana_data = (self.ana_data_prescale is not None and (trigger_number % self.ana_data_prescale)==0)
852 get_wvfm_data = (self.wvfm_data_prescale is not None and (trigger_number % self.wvfm_data_prescale)==0)
853
854 if not (get_ana_data or get_wvfm_data):
855 return None,None
856
857 n_frames = self.get_n_objget_n_obj(frag)
858 adcs = self.unpacker.np_array_adc(frag)
859
860 daphne_headers = [ self.frame_obj(frag.get_data(iframe*self.frame_obj.sizeof())).get_header() for iframe in range(n_frames) ]
861 timestamp = self.unpacker.np_array_timestamp(frag)
862
863 if (len(adcs)) == 0:
864 return None, None
865
866 det, crate, slot, stream = self.get_det_crate_slot_streamget_det_crate_slot_stream(frag)
867
868 if get_ana_data:
869 ax = 1
870 adc_mean = np.mean(adcs,axis=ax)
871 adc_rms = np.std(adcs,axis=ax)
872 adc_max = np.max(adcs,axis=ax)
873 adc_min = np.min(adcs,axis=ax)
874 adc_median = np.median(adcs,axis=ax)
875 ts_max = np.argmax(adcs,axis=ax)*self.SAMPLING_PERIOD + timestamp
876 ts_min = np.argmin(adcs,axis=ax)*self.SAMPLING_PERIOD + timestamp
877
878 ana_data = [ DAPHNEAnalysisData(run=frh.run_number,
879 trigger=frh.trigger_number,
880 sequence=frh.sequence_number,
881 src_id=frh.element_id.id,
882 channel=self.channel_map.get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, daphne_headers[iframe].channel),
883 daphne_chan=daphne_headers[iframe].channel,
884 timestamp_dts=timestamp[iframe],
885 trigger_sample_value=daphne_headers[iframe].trigger_sample_value,
886 threshold=daphne_headers[iframe].threshold,
887 baseline=daphne_headers[iframe].baseline,
888 adc_mean=adc_mean[iframe],
889 adc_rms=adc_rms[iframe],
890 adc_max=adc_max[iframe],
891 adc_min=adc_min[iframe],
892 adc_median=adc_median[iframe],
893 timestamp_max_dts=ts_max[iframe],
894 timestamp_min_dts=ts_min[iframe]) for iframe in range(n_frames) ]
895
896
897 if get_wvfm_data:
898
899 wvfm_data = [ DAPHNEWaveformData(run=frh.run_number,
900 trigger=frh.trigger_number,
901 sequence=frh.sequence_number,
902 src_id=frh.element_id.id,
903 channel=self.channel_map.get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, daphne_headers[iframe].channel),
904 daphne_chan=daphne_headers[iframe].channel,
905 timestamp_dts=timestamp[iframe],
906 timestamps=np.arange(np.size(adcs[iframe,:]))*self.SAMPLING_PERIOD+timestamp[iframe],
907 adcs=adcs[iframe,:]) for iframe in range(n_frames) ]
908
909 return ana_data, wvfm_data
910
913 frame_obj = fddetdataformats.DAPHNEEthFrame
914
916 dh = self.frame_objframe_obj(frag.get_data()).get_daqheader()
917 return dh.det_id, dh.crate_id, dh.slot_id, dh.stream_id
918
__init__(self, channel_map=None, ana_data_prescale=1, wvfm_data_prescale=None)
Definition utils.py:372
get_offline_channel_from_det_crate_slot_stream_chan(det, crate, slot, stream, ch)
Definition utils.py:24
__init__(self, index=None)
Definition utils.py:40
get_all_data(self, in_data=None)
Definition utils.py:43
sparsify_array_diff_locs_and_vals(arr)
Sparsification and desparsifications for arrays.