Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
Fadc250Module.cxx
Go to the documentation of this file.
1
2//
3// Fadc250Module
4// JLab FADC 250 Module
5//
6// Author: Eric Pooser, pooser@jlab.org
7//
8// For documentation pertaining to the FADC250
9// modules refer to the following link:
10// https://coda.jlab.org/drupal/system/files/pdfs/HardwareManual/fADC250/FADC250_Processing_FPGA_Firmware_ver_0x0C01_Description_Instructions.pdf
11//
12// The following data type summary corresponds to firmware version 0x0C00 (06/09/2016)
13// -----------------------------------------------------------------------------------
14// 0 block header
15// 1 block trailer
16// 2 event header
17// 3 trigger time
18// 4 window raw data
19// 5 (reserved)
20// 6 pulse raw data -> Old firmware
21// 7 pulse integral -> Old firmware
22// 8 pulse time -> Old firmware
23// 9 pulse parameters
24// 10 pulse parameters (pedestal) -> Old firmware
25// 11 (reserved)
26// 12 scaler data
27// 13 (reserved)
28// 14 data not valid (empty module)
29// 15 filler (non-data) word
30//
31// Mode Summary
32// ----------------------------------------------
33// Mode 1: data type 4 -> Old firmware
34// Mode 2: data type 6 -> Old firmware
35// Mode 3: data types 7 & 8 -> Old firmware
36// Mode 4: data types 8 & 10 -> Old firmware
37// Mode 7: data types 7 & 8 & 10 -> Old firmware
38// Mode 8: data types 4 & 8 & 10 -> Old firmware
39// Mode 9: data type 9
40// Mode 10: data types 4 & 9
41//
42// fFirmwareVers == 2
43// Current version (0x0C0D, 9/28/17) of the FADC firmware only
44// reports the pedestal once for up to four identifiable pulses in
45// in the readout window. For now, it is best to tell THaEvData
46// that there are n pedestals associated with n identified pulses
47// and point the GetData method to the first entry in the pedestal
48// data vector. This is now the non firmware specific behavior.
49//
51
52#include "Fadc250Module.h"
53#include "THaSlotData.h"
54#include "TMath.h"
55
56#include <unistd.h>
57#include <iostream>
58#include <iomanip>
59#include <numeric>
60#include <cassert>
61#include <stdexcept>
62#include <map>
63#include <sstream>
64
65using namespace std;
66
67// #define DEBUG
68// #define WITH_DEBUG
69
70#ifdef DEBUG
71#include <fstream>
72#endif
73
74namespace Decoder {
75
76using vsiz_t = vector<int>::size_type;
77
79 DoRegister(ModuleType("Decoder::Fadc250Module", 250));
80
81//_____________________________________________________________________________
83 : PipeliningModule(crate, slot), fadc_data{}, fPulseData(NADCCHAN),
84 data_type_4(false), data_type_6(false), data_type_7(false),
85 data_type_8(false), data_type_9(false), data_type_10(false),
86 block_header_found(false), block_trailer_found(false),
87 event_header_found(false), slots_match(false)
88{
89 IsInit = false;
91}
92
93//_____________________________________________________________________________
94#if defined DEBUG && defined WITH_DEBUG
96{
97 // delete fDebugFile; fDebugFile = 0;
98}
99#else
101#endif
102
103//_____________________________________________________________________________
110
111//_____________________________________________________________________________
112// Clear all data vectors
113inline
115{
116 // Clear all data objects
117 assert(fPulseData.size() == NADCCHAN); // Initialization error in constructor
118 for( uint32_t i = 0; i < NADCCHAN; i++ ) {
119 fPulseData[i].clear();
120 }
121}
122
123//_____________________________________________________________________________
124// Require that slot from base class and slot from
125// data match before populating data vectors
126inline
127void Fadc250Module::PopulateDataVector( vector<uint32_t>& data_vector,
128 uint32_t data ) const
129{
130 if( slots_match )
131 data_vector.push_back(data);
132}
133
134//_____________________________________________________________________________
135// Sum elements contained in data vector
136inline
137uint32_t Fadc250Module::SumVectorElements( const vector<uint32_t>& data_vector )
138{
139 return accumulate(data_vector.begin(), data_vector.end(), 0U);
140}
141
142//_____________________________________________________________________________
144{
145 // Clear event-by-event data
147 fadc_data.clear();
149 // Initialize data_type_def to FILLER and data types to false
150 data_type_def = 15;
151 // Initialize data types to false
154}
155
156//_____________________________________________________________________________
158{
159#if defined DEBUG && defined WITH_DEBUG
160 // This will make a HUGE output
161 // delete fDebugFile; fDebugFile = 0;
162 // fDebugFile = new ofstream;
163 // fDebugFile->open("fadcdebug.dat");
164#endif
165 Clear();
166 IsInit = true;
167 fName = "FADC250 JLab Flash ADC Module";
168}
169
170//_____________________________________________________________________________
172{
173 cout << "FADC250 Decoder has been called" << endl;
174}
175
176//_____________________________________________________________________________
178 UInt_t chan ) const
179{
180 vsiz_t ret = 0;
181 switch( emode ) {
182 case kSampleADC:
183 ret = fPulseData[chan].samples.size();
184 break;
185 case kPulseIntegral:
186 ret = fPulseData[chan].integral.size();
187 break;
188 case kPulseTime:
189 ret = fPulseData[chan].time.size();
190 break;
191 case kPulsePeak:
192 ret = fPulseData[chan].peak.size();
193 break;
194 case kPulsePedestal:
195 if( fFirmwareVers == 2 )
196 ret = fPulseData[chan].pedestal.size();
197 else
198 ret = fPulseData[chan].integral.size();
199 break;
200 case kCoarseTime:
201 ret = fPulseData[chan].coarse_time.size();
202 break;
203 case kFineTime:
204 ret = fPulseData[chan].fine_time.size();
205 break;
206 }
207 if( ret > kMaxUInt )
208 throw overflow_error("ERROR! Fadc250Module::GetNumEvents: "
209 "integer overflow");
210 return ret;
211}
212
213//_____________________________________________________________________________
215 UInt_t ievent ) const
216{
217 switch( emode ) {
218 case kSampleADC:
219 return GetPulseSamplesData(chan, ievent);
220 case kPulseIntegral:
221 return GetPulseIntegralData(chan, ievent);
222 case kPulseTime:
223 return GetPulseTimeData(chan, ievent);
224 case kPulsePeak:
225 return GetPulsePeakData(chan, ievent);
226 case kPulsePedestal:
227 return GetPulsePedestalData(chan, ievent);
228 case kCoarseTime:
229 return GetPulseCoarseTimeData(chan, ievent);
230 case kFineTime:
231 return GetPulseFineTimeData(chan, ievent);
232 }
233 return 0;
234}
235
236//_____________________________________________________________________________
238{
239 vsiz_t nevent = fPulseData[chan].integral.size();
240 if( ievent >= nevent ) {
241 cout << "ERROR:: Fadc250Module:: GetPulseIntegralData:: invalid event number for slot = " << fSlot << ", channel = "
242 << chan << endl;
243 return kMaxUInt;
244 }
245 if( nevent == 0 ) {
246 cout << "ERROR:: Fadc250Module:: GetPulseIntegralData:: data vector empty for slot = " << fSlot << ", channel = "
247 << chan << endl;
248 return kMaxUInt;
249 } else {
250#ifdef WITH_DEBUG
251 if( fDebugFile )
252 *fDebugFile << "Fadc250Module::GetPulseIntegralData channel "
253 << chan << ", event " << ievent << " = "
254 << fPulseData[chan].integral[ievent] << endl;
255#endif
256 return fPulseData[chan].integral[ievent];
257 }
258}
259
260//_____________________________________________________________________________
262{
263 vsiz_t nevent = fPulseData[chan].samples.size();
264 if( nevent == 0 ) {
265 cout << "ERROR:: Fadc250Module:: GetEmulatedPulseIntegralData:: data vector empty for slot = " << fSlot
266 << ", channel = " << chan << "\n"
267 << "Ensure that FADC is operating in mode 1 OR 8" << endl;
268 return kMaxUInt;
269 } else {
270#ifdef WITH_DEBUG
271 if( fDebugFile )
272 *fDebugFile << "Fadc250Module::GetEmulatedPulseIntegralData channel "
273 << chan << " = " << SumVectorElements(fPulseData[chan].samples) << endl;
274#endif
276 }
277}
278
279//_____________________________________________________________________________
281{
282 vsiz_t nevent = fPulseData[chan].time.size();
283 if( ievent >= nevent ) {
284 cout << "ERROR:: Fadc250Module:: GetPulseTimeData:: invalid event number for slot = " << fSlot << ", channel = "
285 << chan << endl;
286 return kMaxUInt;
287 }
288 if( nevent == 0 ) {
289 cout << "ERROR:: Fadc250Module:: GetPulseTimeData:: data vector empty for slot = " << fSlot << ", channel = "
290 << chan << endl;
291 return kMaxUInt;
292 } else {
293#ifdef WITH_DEBUG
294 if( fDebugFile )
295 *fDebugFile << "Fadc250Module::GetPulseTimeData channel "
296 << chan << ", event " << ievent << " = "
297 << fPulseData[chan].time[ievent] << endl;
298#endif
299 return fPulseData[chan].time[ievent];
300 }
301}
302
303//_____________________________________________________________________________
305{
306 vsiz_t nevent = fPulseData[chan].coarse_time.size();
307 if( ievent >= nevent ) {
308 cout << "ERROR:: Fadc250Module:: GetPulseCoarseTimeData:: invalid event number for slot = " << fSlot
309 << ", channel = " << chan << endl;
310 return kMaxUInt;
311 }
312 if( nevent == 0 ) {
313 cout << "ERROR:: Fadc250Module:: GetPulseCoarseTimeData:: data vector empty for slot = " << fSlot << ", channel = "
314 << chan << endl;
315 return kMaxUInt;
316 } else {
317#ifdef WITH_DEBUG
318 if( fDebugFile )
319 *fDebugFile << "Fadc250Module::GetPulseCoarseTimeData channel "
320 << chan << ", event " << ievent << " = "
321 << fPulseData[chan].coarse_time[ievent] << endl;
322#endif
323 return fPulseData[chan].coarse_time[ievent];
324 }
325}
326
327//_____________________________________________________________________________
329{
330 vsiz_t nevent = fPulseData[chan].fine_time.size();
331 if( ievent >= nevent ) {
332 cout << "ERROR:: Fadc250Module:: GetPulseCoarseTimeData:: invalid event number for slot = " << fSlot
333 << ", channel = " << chan << endl;
334 return kMaxUInt;
335 }
336 if( nevent == 0 ) {
337 cout << "ERROR:: Fadc250Module:: GetPulseCoarseTimeData:: data vector empty for slot = " << fSlot << ", channel = "
338 << chan << endl;
339 return kMaxUInt;
340 } else {
341#ifdef WITH_DEBUG
342 if( fDebugFile )
343 *fDebugFile << "Fadc250Module::GetPulseFineTimeData channel "
344 << chan << ", event " << ievent << " = "
345 << fPulseData[chan].fine_time[ievent] << endl;
346#endif
347 return fPulseData[chan].fine_time[ievent];
348 }
349}
350
351//_____________________________________________________________________________
353{
354 vsiz_t nevent = fPulseData[chan].peak.size();
355 if( ievent >= nevent ) {
356 cout << "ERROR:: Fadc250Module:: GetPulsePeakData:: invalid event number for slot = " << fSlot << ", channel = "
357 << chan << endl;
358 return kMaxUInt;
359 }
360 if( nevent == 0 ) {
361 cout << "ERROR:: Fadc250Module:: GetPulsePeakData:: data vector empty for slot = " << fSlot << ", channel = "
362 << chan << endl;
363 return kMaxUInt;
364 } else {
365#ifdef WITH_DEBUG
366 if( fDebugFile )
367 *fDebugFile << "Fadc250Module::GetPulsePeakData channel "
368 << chan << ", event " << ievent << " = "
369 << fPulseData[chan].peak[ievent] << endl;
370#endif
371 return fPulseData[chan].peak[ievent];
372 }
373}
374
375//_____________________________________________________________________________
377{
378 vsiz_t nevent = fPulseData[chan].pedestal.size();
379 if( nevent == 0 ) {
380 cout << "ERROR:: Fadc250Module:: GetPulsePedestalData:: data vector empty for slot = " << fSlot << ", channel = "
381 << chan << endl;
382 return kMaxUInt;
383 }
384 if( fFirmwareVers == 2 ) {
385 if( ievent >= nevent ) {
386 cout << "ERROR:: Fadc250Module:: GetPulsePedestalData:: invalid data vector size = " << fSlot << ", channel = "
387 << chan << endl;
388 return kMaxUInt;
389 } else {
390#ifdef WITH_DEBUG
391 if( fDebugFile )
392 *fDebugFile << "Fadc250Module::GetPulsePedestalData channel "
393 << chan << ", event " << ievent << " = "
394 << fPulseData[chan].pedestal[ievent] << endl;
395#endif
396 return fPulseData[chan].pedestal[ievent];
397 }
398 }
399 if( nevent != 1 ) {
400 cout << "ERROR:: Fadc250Module:: GetPulsePedestalData:: invalid data vector size = " << fSlot << ", channel = "
401 << chan << endl;
402 return kMaxUInt;
403 } else {
404#ifdef WITH_DEBUG
405 if( fDebugFile )
406 *fDebugFile << "Fadc250Module::GetPulsePedestalData channel "
407 << chan << ", event " << ievent << " = "
408 << fPulseData[chan].pedestal[0] << endl;
409#endif
410 return fPulseData[chan].pedestal[0];
411 }
412}
413
414//_____________________________________________________________________________
416{
417 vsiz_t nevent = fPulseData[chan].pedestal_quality.size();
418 if( nevent == 0 ) {
419 cout << "ERROR:: Fadc250Module:: GetPedestalQuality:: data vector empty for slot = " << fSlot << ", channel = "
420 << chan << endl;
421 return kMaxUInt;
422 }
423 if( nevent != 1 ) {
424 cout << "ERROR:: Fadc250Module:: GetPedestalQuality:: invalid data vector size = " << fSlot << ", channel = "
425 << chan << endl;
426 return kMaxUInt;
427 } else {
428#ifdef WITH_DEBUG
429 if( fDebugFile )
430 *fDebugFile << "Fadc250Module::GetPedestalQuality channel "
431 << chan << ", event " << ievent << " = "
432 << fPulseData[chan].pedestal_quality[0] << endl;
433#endif
434 return fPulseData[chan].pedestal_quality[0];
435 }
436}
437
438//_____________________________________________________________________________
440{
441 vsiz_t nevent = fPulseData[chan].overflow.size();
442 if( ievent >= nevent ) {
443 cout << "ERROR:: Fadc250Module:: GetOverflowBit:: invalid event number for slot = " << fSlot << ", channel = "
444 << chan << endl;
445 return kMaxUInt;
446 }
447 if( nevent == 0 ) {
448 cout << "ERROR:: Fadc250Module:: GetOverflowBit:: data vector empty for slot = " << fSlot << ", channel = " << chan
449 << endl;
450 return kMaxUInt;
451 } else {
452#ifdef WITH_DEBUG
453 if( fDebugFile )
454 *fDebugFile << "Fadc250Module::GetOverflowBit channel "
455 << chan << ", event " << ievent << " = "
456 << fPulseData[chan].overflow[ievent] << endl;
457#endif
458 return fPulseData[chan].overflow[ievent];
459 }
460}
461
462//_____________________________________________________________________________
464{
465 vsiz_t nevent = fPulseData[chan].underflow.size();
466 if( ievent >= nevent ) {
467 cout << "ERROR:: Fadc250Module:: GetUnderflowBit:: invalid event number for slot = " << fSlot << ", channel = "
468 << chan << endl;
469 return kMaxUInt;
470 }
471 if( nevent == 0 ) {
472 cout << "ERROR:: Fadc250Module:: GetUnderflowBit:: data vector empty for slot = " << fSlot << ", channel = " << chan
473 << endl;
474 return kMaxUInt;
475 } else {
476#ifdef WITH_DEBUG
477 if( fDebugFile )
478 *fDebugFile << "Fadc250Module::GetUnderflowBit channel "
479 << chan << ", event " << ievent << " = "
480 << fPulseData[chan].underflow[ievent] << endl;
481#endif
482 return fPulseData[chan].underflow[ievent];
483 }
484}
485
486//_____________________________________________________________________________
488{
489 // Truncate to 32 bits
490 UInt_t shorttime = fadc_data.trig_time;
491#ifdef WITH_DEBUG
492 if( fDebugFile )
493 *fDebugFile << "Fadc250Module::GetTriggerTime = "
494 << fadc_data.trig_time << " " << shorttime << endl;
495#endif
496 return shorttime;
497}
498
499//_____________________________________________________________________________
501{
502 vsiz_t nevent = fPulseData[chan].samples.size();
503 if( ievent >= nevent ) {
504 cout << "ERROR:: Fadc250Module:: GetPulseSamplesData:: invalid event number for slot = " << fSlot << ", channel = "
505 << chan << endl;
506 return kMaxUInt;
507 }
508 if( nevent == 0 ) {
509 cout << "ERROR:: Fadc250Module:: GetPulseSamplesData:: data vector empty for slot = " << fSlot << ", channel = "
510 << chan << endl;
511 return kMaxUInt;
512 } else {
513#ifdef WITH_DEBUG
514 if( fDebugFile )
515 *fDebugFile << "Fadc250Module::GetPulseSamplesData channel "
516 << chan << ", event " << ievent << " = "
517 << fPulseData[chan].samples[ievent] << endl;
518#endif
519 return fPulseData[chan].samples[ievent];
520 }
521}
522
523//_____________________________________________________________________________
525{
526 vsiz_t nevent = fPulseData[chan].samples.size();
527 if( nevent == 0 ) {
528 cout << "ERROR:: Fadc250Module:: GetPulseSamplesVector:: data vector empty for slot = " << fSlot << ", channel = "
529 << chan << endl;
530 return {};
531 } else {
532#ifdef WITH_DEBUG
533 if( fDebugFile )
534 *fDebugFile << "Fadc250Module::GetPulseSamplesVector channel "
535 << chan << " = " << &fPulseData[chan].samples << endl;
536#endif
537 return fPulseData[chan].samples;
538 }
539}
540
541//_____________________________________________________________________________
543{
544#ifdef WITH_DEBUG
545 if( !fDebugFile ) return;
546 *fDebugFile << "start, Print Data Type " << endl;
547 if( data_type_4 ) *fDebugFile << "data type 4" << endl;
548 if( data_type_6 ) *fDebugFile << "data type 6" << endl;
549 if( data_type_7 ) *fDebugFile << "data type 7" << endl;
550 if( data_type_8 ) *fDebugFile << "data type 8" << endl;
551 if( data_type_9 ) *fDebugFile << "data type 9" << endl;
552 if( data_type_10 ) *fDebugFile << "data type 10" << endl;
553 *fDebugFile << "end, Print Data Type " << endl;
554#endif
555}
556
557//_____________________________________________________________________________
559{
560 if( fFirmwareVers == 1 ) { // some "older" firmware
561 if( data_type_4 && data_type_6 ) return 8;
562 if( data_type_7 ) return 7;
563 }
564 if( data_type_4 && !(data_type_6) && !(data_type_7) && !(data_type_8) && !(data_type_9) && !(data_type_10) ) return 1;
565 else if( !(data_type_4) && data_type_6 && !(data_type_7) && !(data_type_8) && !(data_type_9) &&
566 !(data_type_10) )
567 return 2;
568 else if( !(data_type_4) && !(data_type_6) && data_type_7 && data_type_8 && !(data_type_9) && !(data_type_10) )
569 return 3;
570 else if( !(data_type_4) && !(data_type_6) && !(data_type_7) && data_type_8 && !(data_type_9) && data_type_10 )
571 return 4;
572 else if( !(data_type_4) && !(data_type_6) && data_type_7 && data_type_8 && !(data_type_9) && data_type_10 ) return 7;
573 else if( data_type_4 && !(data_type_6) && !(data_type_7) && data_type_8 && !(data_type_9) && data_type_10 ) return 8;
574 else if( !(data_type_4) && !(data_type_6) && !(data_type_7) && !(data_type_8) && data_type_9 &&
575 !(data_type_10) )
576 return 9;
577 else if( data_type_4 && !(data_type_6) && !(data_type_7) && !(data_type_8) && data_type_9 && !(data_type_10) )
578 return 10;
579 else {
580 // PrintDataType();
581 cout << "ERROR:: Fadc250Module:: GetFadcMode:: FADC is in invalid mode for slot = " << fSlot << endl;
582 if( fDebugFile )
583 *fDebugFile << "ERROR:: Fadc250Module:: GetFadcMode:: FADC is in invalid mode for slot = " << fSlot << endl;
584 return -1;
585 }
586}
587
588//_____________________________________________________________________________
590{
591 assert(chan < NADCCHAN);
592 UInt_t sz = 0;
595 // For some "old" firmware version
596 if( fFirmwareVers == 1 ) {
597 if( mode == 7 &&
598 (sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size() )
599 return sz;
600 if( mode == 8 )
601 return fPulseData[chan].samples.size();
602 }
603 // The rest for "modern" firmware
604 if( mode == 1 )
605 return 1;
606 else if(
607 (mode == 7 &&
608 ((sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size()) &&
609 (fPulseData[chan].pedestal.size() == sz) &&
610 (fPulseData[chan].peak.size() == sz)) ||
611 (mode == 8 &&
612 (sz = fPulseData[chan].time.size()) == fPulseData[chan].pedestal.size() &&
613 fPulseData[chan].peak.size() == sz) ||
614 (mode == 9 && fFirmwareVers == 2 &&
615 (sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size() &&
616 fPulseData[chan].pedestal.size() == sz &&
617 fPulseData[chan].peak.size() == sz) ||
618 (mode == 10 && fFirmwareVers == 2 &&
619 (sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size() &&
620 fPulseData[chan].pedestal.size() == sz &&
621 fPulseData[chan].peak.size() == sz) ||
622 (mode == 9 &&
623 (sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size() &&
624 fPulseData[chan].pedestal.size() <= 1 &&
625 fPulseData[chan].peak.size() == sz) ||
626 (mode == 10 &&
627 (sz = fPulseData[chan].integral.size()) == fPulseData[chan].time.size() &&
628 fPulseData[chan].pedestal.size() <= 1 &&
629 fPulseData[chan].peak.size() == sz)
630 ) {
631#ifdef WITH_DEBUG
632 if( fDebugFile )
633 *fDebugFile << "Fadc250Module::GetNumFadcEvents channel "
634 << chan << " = " << sz << endl;
635#endif
636 return sz;
637 } else {
638 cout
639 << "ERROR:: Fadc250Module:: GetNumFadcEvents:: FADC not in acceptable mode or data vector sizes do not match for for slot = "
640 << fSlot << ", channel = " << chan << endl;
641 cout << "ERROR:: Fadc250Module:: GetNumFadcEvents:: FADC data potentially corrupt for event " << fadc_data.evt_num
642 << ", PROCEED WITH CAUTION!" << endl;
643 return kMaxUInt;
644 }
645}
646
647//_____________________________________________________________________________
649{
651 if( (mode == 1) || (mode == 8) || (mode == 10) ) {
652 vsiz_t nsamples = fPulseData[chan].samples.size();
653 if( ievent >= nsamples ) {
654 cout << "ERROR:: Fadc250Module:: GetNumFadcSamples:: invalid event number for slot = " << fSlot << ", channel = "
655 << chan << endl;
656 return kMaxUInt;
657 }
658 if( nsamples == 0 ) {
659 //cout << "ERROR:: Fadc250Module:: GetNumFadcSamples:: data vector empty for slot = " << fSlot << ", channel = " << chan << endl;
660 return kMaxUInt;
661 } else {
662#ifdef WITH_DEBUG
663 if( fDebugFile )
664 *fDebugFile << "Fadc250Module::GetNumFadcSamples channel "
665 << chan << ", event " << ievent << " = "
666 << fPulseData[chan].samples.size() << endl;
667#endif
668 return fPulseData[chan].samples.size();
669 }
670 } else {
671 cout << "ERROR:: Fadc250Module:: GetNumFadcSamples:: FADC is not in Mode 1, 8, or 10 for slot = " << fSlot
672 << ", channel = " << chan << endl;
673 return kMaxUInt;
674 }
675}
676
677//_____________________________________________________________________________
678inline
679void Fadc250Module::DecodeBlockHeader( UInt_t pdat, uint32_t data_type_id )
680{
681 static const string here{ "Fadc250Module::DecodeBlockHeader" };
682
683 if( data_type_id ) {
684 block_header_found = true; // Set to true if found
685 fadc_data.slot_blk_hdr = (pdat >> 22) & 0x1F; // Slot number (set by VME64x backplane), mask 5 bits
686 // Debug output
687#ifdef WITH_DEBUG
688 if( fDebugFile )
689 *fDebugFile << "Fadc250Module::Decode:: Slot from FADC block header = " << fadc_data.slot_blk_hdr << endl;
690#endif
691 // Ensure that slots from cratemap and FADC match
692 slots_match = ((uint32_t) fSlot == fadc_data.slot_blk_hdr);
693 if( !slots_match )
694 return;
695 fadc_data.mod_id = (pdat >> 18) & 0xF; // Module ID ('1' for FADC250), mask 4 bits
696 fadc_data.iblock_num = (pdat >> 8) & 0x3FF; // Event block number, mask 10 bits
697 fadc_data.nblock_events = (pdat >> 0) & 0xFF; // Number of events in block, mask 8 bits
698 // Debug output
699#ifdef WITH_DEBUG
700 if( fDebugFile )
701 *fDebugFile << "Fadc250Module::Decode:: FADC BLOCK HEADER"
702 << " >> data = " << hex << pdat << dec
703 << " >> slot = " << fadc_data.slot_blk_hdr
704 << " >> module id = " << fadc_data.mod_id << " >> event block number = " << fadc_data.iblock_num
705 << " >> num events in block = " << fadc_data.nblock_events
706 << endl;
707#endif
708 } else {
709 fadc_data.PL =
710 (pdat >> 18) & 0x7FF; // Number of samples before trigger point for processing to begin, mask 11 bits
711 fadc_data.NSB =
712 (pdat >> 9) & 0x1FF; // Number of samples before threshold crossing to include processing, mask 9 bits
713 fadc_data.NSA =
714 (pdat >> 0) & 0x1FF; // Number of samples after threshold crossing to include processing, mask 9 bits
715 // Debug output
716#ifdef WITH_DEBUG
717 if( fDebugFile )
718 *fDebugFile << "Fadc250Module::Decode:: FADC BLOCK HEADER"
719 << " >> data = " << hex << pdat << dec
720 << " >> NSB trigger point (PL) = " << fadc_data.PL
721 << " >> NSB threshold crossing = " << fadc_data.NSA
722 << " >> NSA threshold crossing = " << fadc_data.NSB
723 << endl;
724#endif
725 }
726}
727
728//_____________________________________________________________________________
730{
731 block_trailer_found = true;
732 fadc_data.slot_blk_trl = (pdat >> 22) & 0x1F; // Slot number (set by VME64x backplane), mask 5 bits
733 fadc_data.nwords_inblock = (pdat >> 0) & 0x3FFFFF; // Total number of words in block of events, mask 22 bits
734 // Debug output
735#ifdef WITH_DEBUG
736 if( fDebugFile )
737 *fDebugFile << "Fadc250Module::Decode:: FADC BLOCK TRAILER"
738 << " >> data = " << hex << pdat << dec
739 << " >> slot = " << fadc_data.slot_blk_trl
740 << " >> nwords in block = " << fadc_data.nwords_inblock
741 << endl;
742#endif
743}
744
745//_____________________________________________________________________________
747{
748 event_header_found = true;
749 // For firmware versions pre 0x0C00 (06/09/2016)
750 fadc_data.slot_evt_hdr = (pdat >> 22) & 0x1F; // Slot number (set by VME64x backplane), mask 5 bits
751 // fadc_data.evt_num = (pdat >> 0) & 0x3FFFFF; // Event number, mask 22 bits
752 // For firmware versions post 0x0C00 (06/09/2016)
753 fadc_data.eh_trig_time = (pdat >> 12) & 0x3FF; // Event header trigger time
754 fadc_data.trig_num = (pdat >> 0) & 0xFFF; // Trigger number
755 // Debug output
756 // #ifdef WITH_DEBUG
757 // if (fDebugFile)
758 // *fDebugFile << "Fadc250Module::Decode:: FADC EVENT HEADER >> data = " << hex
759 // << data << dec << " >> slot = " << fadc_data.slot_evt_hdr
760 // << " >> event number = " << fadc_data.evt_num << endl;
761 // #endif
762#ifdef WITH_DEBUG
763 if( fDebugFile )
764 *fDebugFile << "Fadc250Module::Decode:: FADC EVENT HEADER"
765 << " >> data = " << hex << pdat << dec
766 << " >> event header trigger time = " << fadc_data.eh_trig_time
767 << " >> trigger number = " << fadc_data.trig_num
768 << endl;
769#endif
770}
771
772//_____________________________________________________________________________
773void Fadc250Module::DecodeTriggerTime( UInt_t pdat, uint32_t data_type_id )
774{
775 if( data_type_id ) // Trigger time word 1
776 fadc_data.trig_time_w1 = (pdat >> 0) & 0xFFFFFF; // Time = T_D T_E T_F
777 else // Trigger time word 2
778 fadc_data.trig_time_w2 = (pdat >> 0) & 0xFFFFFF; // Time = T_A T_B T_C
779 // Time = T_A T_B T_C T_D T_E T_F
780 fadc_data.trig_time = (fadc_data.trig_time_w2 << 24) | fadc_data.trig_time_w1;
781 // Debug output
782#ifdef WITH_DEBUG
783 if( fDebugFile )
784 *fDebugFile << "Fadc250Module::Decode:: FADC TRIGGER TIME"
785 << " >> data = " << hex << pdat << dec
786 << " >> trigger time word 1 = " << fadc_data.trig_time_w1
787 << " >> trigger time word 2 = " << fadc_data.trig_time_w2
788 << " >> trigger time = " << fadc_data.trig_time
789 << endl;
790#endif
791}
792
793//_____________________________________________________________________________
794void Fadc250Module::DecodeWindowRawData( UInt_t pdat, uint32_t data_type_id )
795{
796 data_type_4 = true;
797 if( data_type_id == 1 ) {
798 fadc_data.chan = (pdat >> 23) & 0xF; // FADC channel number
799 fadc_data.win_width = (pdat >> 0) & 0xFFF; // Window width
800 // Debug output
801#ifdef WITH_DEBUG
802 if( fDebugFile )
803 *fDebugFile << "Fadc250Module::Decode:: FADC WINDOW RAW DATA"
804 << " >> data = " << hex << pdat << dec
805 << " >> channel = " << fadc_data.chan
806 << " >> window width = " << fadc_data.win_width
807 << endl;
808#endif
809 } else {
810 assert(data_type_id == 0); // Ensure this is a data continuation word
811 Bool_t invalid_1 = (pdat >> 29) & 0x1; // Check if sample x is invalid
812 Bool_t invalid_2 = (pdat >> 13) & 0x1; // Check if sample x+1 is invalid
813 uint32_t sample_1 = 0; // ADC sample x (includes overflow bit)
814 uint32_t sample_2 = 0; // ADC sample x+1 (includes overflow bit)
815 if( !invalid_1 ) sample_1 = (pdat >> 16) & 0x1FFF; // If sample x is valid, assign value
816 if( !invalid_2 ) sample_2 = (pdat >> 0) & 0x1FFF; // If sample x+1 is valid, assign value
817
818 PopulateDataVector(fPulseData[fadc_data.chan].samples, sample_1); // Sample 1
819 fadc_data.invalid_samples |= invalid_1; // Invalid samples
820 fadc_data.overflow = (sample_1 >> 12) & 0x1; // Sample 1 overflow bit
821 if( invalid_2 ) { // Skip last sample if not expected and flagged as invalid
822 if( fPulseData[fadc_data.chan].samples.size() == fadc_data.win_width )
823 return;
824 }
825
826 PopulateDataVector(fPulseData[fadc_data.chan].samples, sample_2); // Sample 2
827 fadc_data.invalid_samples |= invalid_2; // Invalid samples
828 fadc_data.overflow = (sample_2 >> 12) & 0x1; // Sample 2 overflow bit
829 // Debug output
830#ifdef WITH_DEBUG
831 if( fDebugFile )
832 *fDebugFile << "Fadc250Module::Decode:: FADC WINDOW RAW DATA"
833 << " >> data = " << hex << pdat << dec
834 << " >> channel = " << fadc_data.chan
835 << " >> sample 1 = " << sample_1
836 << " >> sample 2 = " << sample_2
837 << " >> size of fPulseSamples = "
838 << fPulseData[fadc_data.chan].samples.size()
839 << endl;
840#endif
841 } // FADC data sample for window raw data
842}
843
844//_____________________________________________________________________________
845void Fadc250Module::DecodePulseRawData( UInt_t pdat, uint32_t data_type_id )
846{
847 data_type_6 = true;
848 if( data_type_id == 1 ) {
849 fadc_data.chan = (pdat >> 23) & 0xF; // FADC channel number
850 fadc_data.pulse_num = (pdat >> 21) & 0x3; // FADC pulse number
851 fadc_data.sample_num_tc = (pdat >> 0) & 0x3FF; // FADC sample number of threshold crossing
852 // Debug output
853#ifdef WITH_DEBUG
854 if( fDebugFile )
855 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE RAW DATA"
856 << " >> data = " << hex << pdat << dec
857 << " >> channel = " << fadc_data.chan
858 << " >> window width = " << fadc_data.win_width
859 << endl;
860#endif
861 } else {
862 assert(data_type_id == 0); // Ensure this is a data continuation word
863 Bool_t invalid_1 = (pdat >> 29) & 0x1; // Check if sample x is invalid
864 Bool_t invalid_2 = (pdat >> 13) & 0x1; // Check if sample x+1 is invalid
865 uint32_t sample_1 = 0; // ADC sample x (includes overflow bit)
866 uint32_t sample_2 = 0; // ADC sample x+1 (includes overflow bit)
867 if( !invalid_1 ) sample_1 = (pdat >> 16) & 0x1FFF; // If sample x is valid, assign value
868 if( !invalid_2 ) sample_2 = (pdat >> 0) & 0x1FFF; // If sample x+1 is valid, assign value
869
870 PopulateDataVector(fPulseData[fadc_data.chan].samples, sample_1); // Sample 1
871 fadc_data.invalid_samples |= invalid_1; // Invalid samples
872 fadc_data.overflow = (sample_1 >> 12) & 0x1; // Sample 1 overflow bit
873 if( invalid_2 ) { // Skip last sample if not expected and flagged as invalid
874 if( fPulseData[fadc_data.chan].samples.size() == fadc_data.win_width )
875 return;
876 }
877
878 PopulateDataVector(fPulseData[fadc_data.chan].samples, sample_2); // Sample 2
879 fadc_data.invalid_samples |= invalid_2; // Invalid samples
880 fadc_data.overflow = (sample_2 >> 12) & 0x1; // Sample 2 overflow bit
881 // Debug output
882#ifdef WITH_DEBUG
883 if( fDebugFile )
884 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE RAW DATA"
885 << " >> data = " << hex << pdat << dec
886 << " >> sample 1 = " << sample_1
887 << " >> sample 2 = " << sample_2
888 << " >> size of fPulseSamples = "
889 << fPulseData[fadc_data.chan].samples.size()
890 << endl;
891#endif
892 } // FADC data sample loop for pulse raw data
893}
894
895//_____________________________________________________________________________
897{
898 data_type_7 = true;
899 fadc_data.chan = (pdat >> 23) & 0xF; // FADC channel number
900 fadc_data.pulse_num = (pdat >> 21) & 0x3; // FADC pulse number
901 fadc_data.qual_factor = (pdat >> 19) & 0x3; // FADC quality factor (0-3)
902 fadc_data.pulse_integral = (pdat >> 0) & 0x7FFFF; // FADC pulse integral
903 // Store data in arrays of vectors
904 PopulateDataVector(fPulseData[fadc_data.chan].integral, fadc_data.pulse_integral);
905 // Debug output
906#ifdef WITH_DEBUG
907 if( fDebugFile )
908 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE INTEGRAL"
909 << " >> data = " << hex << pdat << dec
910 << " >> chan = " << fadc_data.chan
911 << " >> pulse num = " << fadc_data.pulse_num
912 << " >> quality factor = " << fadc_data.qual_factor
913 << " >> pulse integral = " << fadc_data.pulse_integral
914 << endl;
915#endif
916}
917
918//_____________________________________________________________________________
920{
921 data_type_8 = true;
922 fadc_data.chan = (pdat >> 23) & 0xF; // FADC channel number
923 fadc_data.pulse_num = (pdat >> 21) & 0x3; // FADC pulse number
924 fadc_data.qual_factor = (pdat >> 19) & 0x3; // FADC quality factor (0-3)
925 fadc_data.coarse_pulse_time = (pdat >> 6) & 0x1FF; // FADC coarse time (4 ns/count)
926 fadc_data.fine_pulse_time = (pdat >> 0) & 0x3F; // FADC fine time (0.0625 ns/count)
927 fadc_data.time = (pdat >> 0) & 0x7FFF; // FADC time (0.0625 ns/count, bmoffit)
928 // Store data in arrays of vectors
929 PopulateDataVector(fPulseData[fadc_data.chan].coarse_time, fadc_data.coarse_pulse_time);
930 PopulateDataVector(fPulseData[fadc_data.chan].fine_time, fadc_data.fine_pulse_time);
931 PopulateDataVector(fPulseData[fadc_data.chan].time, fadc_data.time);
932 // Debug output
933#ifdef WITH_DEBUG
934 if( fDebugFile )
935 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE TIME"
936 << " >> data = " << hex << pdat << dec
937 << " >> chan = " << fadc_data.chan
938 << " >> pulse num = " << fadc_data.pulse_num
939 << " >> quality factor = " << fadc_data.qual_factor
940 << " >> coarse time = " << fadc_data.coarse_pulse_time
941 << " >> fine time = " << fadc_data.fine_pulse_time
942 << " >> time = " << fadc_data.time
943 << endl;
944#endif
945}
946
947//_____________________________________________________________________________
948void Fadc250Module::DecodePulseParameters( UInt_t pdat, uint32_t data_type_id )
949{
950 data_type_9 = true;
951 if( data_type_id == 1 ) {
952 fadc_data.evnt_of_blk = (pdat >> 19) & 0xFF; // Event number within block
953 fadc_data.chan = (pdat >> 15) & 0xF; // FADC channel number
954 fadc_data.qual_factor = (pdat >> 14) & 0x1; // Pedestal quality
955 fadc_data.pedestal_sum = (pdat >> 0) & 0x3FFF; // Pedestal sum
956 // Populate data vectors
957 PopulateDataVector(fPulseData[fadc_data.chan].pedestal, fadc_data.pedestal_sum);
958 PopulateDataVector(fPulseData[fadc_data.chan].pedestal_quality, fadc_data.qual_factor);
959 // Debug output
960#ifdef WITH_DEBUG
961 if( fDebugFile )
962 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE PARAMETERS"
963 << " >> data = " << hex << pdat << dec
964 << " >> chan = " << fadc_data.chan
965 << " >> event of block = " << fadc_data.evnt_of_blk
966 << " >> quality factor = " << fadc_data.qual_factor
967 << " >> pedestal sum = " << fadc_data.pedestal_sum
968 << endl;
969#endif
970 } else {
971 assert(data_type_id == 0); // Ensure this is a data continuation word
972 if( ((pdat >> 30) & 0x1) == 1 ) { // Ensure that data is the integral of pulse window word
973 fadc_data.sample_sum =
974 (pdat >> 12) & 0x3FFFF; // 18-bit sum of raw samples that constitute the pulse data set
975 fadc_data.nsa_ext = (pdat >> 11) & 0x1; // NSA extended beyond PTW
976 fadc_data.samp_overflow = (pdat >> 10) & 0x1; // One or more samples is overflow
977 fadc_data.samp_underflow = (pdat >> 9) & 0x1; // One or more samples is underflow
978 fadc_data.samp_over_thresh =
979 (pdat >> 0) & 0x1FF; // Number of samples within NSA that the pulse is above threshold
980 // Populate data vectors
981 PopulateDataVector(fPulseData[fadc_data.chan].integral, fadc_data.sample_sum);
982 PopulateDataVector(fPulseData[fadc_data.chan].overflow, fadc_data.samp_overflow);
983 PopulateDataVector(fPulseData[fadc_data.chan].underflow, fadc_data.samp_underflow);
984 // Debug output
985#ifdef WITH_DEBUG
986 if( fDebugFile )
987 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE PARAMETERS"
988 << " >> data = " << hex << pdat << dec
989 << " >> chan = " << fadc_data.chan
990 << " >> integral = " << fadc_data.sample_sum
991 << " >> nsa extended = " << fadc_data.nsa_ext
992 << " >> sample overflow = " << fadc_data.samp_overflow
993 << " >> sample underflow = " << fadc_data.samp_underflow
994 << " >> sample over threshold = " << fadc_data.samp_over_thresh
995 << endl;
996#endif
997 }
998 if( ((pdat >> 30) & 0x1) == 0 ) { // Ensure that data is the time of pulse window word
999 fadc_data.coarse_pulse_time = (pdat >> 21) & 0x1FF; // Coarse time (4 ns/count)
1000 fadc_data.fine_pulse_time = (pdat >> 15) & 0x3F; // Fine time (0.0625 ns/count)
1001 fadc_data.time = (pdat >> 15) & 0x7FFF; // FADC time (0.0625 ns/count)
1002 fadc_data.pulse_peak = (pdat >> 3) & 0xFFF; // Pulse peak
1003 fadc_data.peak_beyond_nsa =
1004 (pdat >> 2) & 0x1; // Pulse peak is beyond NSA or could be beyond window end
1005 fadc_data.peak_not_found = (pdat >> 1) & 0x1; // Pulse peak cannot be found
1006 fadc_data.peak_above_maxped =
1007 (pdat >> 0) & 0x1; // 1 or more of first four samples is above either MaxPed or TET
1008 // Populate data vectors
1009 PopulateDataVector(fPulseData[fadc_data.chan].coarse_time, fadc_data.coarse_pulse_time);
1010 PopulateDataVector(fPulseData[fadc_data.chan].fine_time, fadc_data.fine_pulse_time);
1011 PopulateDataVector(fPulseData[fadc_data.chan].time, fadc_data.time);
1012 PopulateDataVector(fPulseData[fadc_data.chan].peak, fadc_data.pulse_peak);
1013 // Debug output
1014#ifdef WITH_DEBUG
1015 if( fDebugFile )
1016 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE PARAMETERS"
1017 << " >> data = " << hex << pdat << dec
1018 << " >> chan = " << fadc_data.chan
1019 << " >> coarse time = " << fadc_data.coarse_pulse_time
1020 << " >> fine time = " << fadc_data.fine_pulse_time
1021 << " >> time = " << fadc_data.time
1022 << " >> pulse peak = " << fadc_data.pulse_peak
1023 << " >> peak beyond nsa = " << fadc_data.peak_beyond_nsa
1024 << " >> peak not found = " << fadc_data.peak_not_found
1025 << " >> peak above maxped = " << fadc_data.peak_above_maxped
1026 << endl;
1027#endif
1028 }
1029 }
1030}
1031
1032//_____________________________________________________________________________
1034{
1035 data_type_10 = true;
1036 fadc_data.chan = (pdat >> 23) & 0xF; // FADC channel number
1037 fadc_data.pulse_num = (pdat >> 21) & 0x3; // FADC pulse number
1038 fadc_data.pedestal = (pdat >> 12) & 0x1FF; // FADC pulse pedestal
1039 fadc_data.pulse_peak = (pdat >> 0) & 0xFFF; // FADC pulse peak
1040 // Store data in arrays of vectors
1041 PopulateDataVector(fPulseData[fadc_data.chan].pedestal, fadc_data.pedestal);
1042 PopulateDataVector(fPulseData[fadc_data.chan].peak, fadc_data.pulse_peak);
1043 // Debug output
1044#ifdef WITH_DEBUG
1045 if( fDebugFile )
1046 *fDebugFile << "Fadc250Module::Decode:: FADC PULSE PEDESTAL"
1047 << " >> data = " << hex << pdat << dec
1048 << " >> chan = " << fadc_data.chan
1049 << " >> pulse num = " << fadc_data.pulse_num
1050 << " >> pedestal = " << fadc_data.pedestal
1051 << " >> pulse peak = " << fadc_data.pulse_peak
1052 << endl;
1053#endif
1054}
1055
1056//_____________________________________________________________________________
1058{
1059 fadc_data.scaler_words = (pdat >> 0) & 0x3F; // FADC scaler words
1060 // Debug output
1061#ifdef WITH_DEBUG
1062 if( fDebugFile )
1063 *fDebugFile << "Fadc250Module::Decode:: FADC SCALER HEADER"
1064 << " >> data = " << hex << pdat << dec
1065 << " >> data words = " << hex << fadc_data.scaler_words << dec
1066 << endl;
1067#endif
1068}
1069
1070//_____________________________________________________________________________
1071void Fadc250Module::UnsupportedType( UInt_t pdat, uint32_t data_type_id )
1072{
1073 // Handle unsupported, invalid, or irrelevant/non-decodable data types
1074#ifdef WITH_DEBUG
1075 // Data type descriptions
1076 static const vector<string> what_text{ "UNDEFINED TYPE",
1077 "DATA NOT VALID",
1078 "FILLER WORD",
1079 "INCORRECT DECODING" };
1080 // Lookup table data_type -> message number
1081 static const map<uint32_t, uint32_t> what_map = {
1082 // undefined type
1083 { 5, 0 },
1084 { 11, 0 },
1085 { 13, 0 },
1086 // data not valid
1087 { 14, 1 },
1088 // filler word
1089 { 15, 2 }
1090 };
1091 auto idx = what_map.find(data_type_def);
1092 // Message index. The last message means this function was called when
1093 // it shouldn't have been called, i.e. coding error in DecodeOneWord
1094 size_t i = (idx == what_map.end()) ? what_text.size() - 1 : idx->second;
1095 const string& what = what_text[i];
1096 ostringstream str;
1097 str << "Fadc250Module::Decode:: " << what
1098 << " >> data = " << hex << pdat << dec
1099 << " >> data type id = " << data_type_id
1100 << endl;
1101 if( fDebugFile )
1102 *fDebugFile << str.str();
1103 if( idx == what_map.end() )
1104 cerr << str.str();
1105#endif
1106}
1107
1108//_____________________________________________________________________________
1110{
1111 assert(pdat);
1112 uint32_t data = *pdat;
1113 uint32_t data_type_id = (data >> 31) & 0x1; // Data type identification, mask 1 bit
1114 if( data_type_id == 1 )
1115 data_type_def = (data >> 27) & 0xF; // Data type defining words, mask 4 bits
1116
1117 // Debug output
1118#ifdef WITH_DEBUG
1119 if( fDebugFile )
1120 *fDebugFile << "Fadc250Module::Decode:: FADC DATA TYPES"
1121 << " >> data = " << hex << data << dec
1122 << " >> data word id = " << data_type_id
1123 << " >> data type = " << data_type_def
1124 << endl;
1125
1126#endif
1127
1128 // Ensure that slots match and do not decode if they differ.
1129 // This should never happen if PipeliningModule::LoadBank selected the
1130 // correct bank.
1131 if( !slots_match && data_type_def != 0 ) {
1132#ifdef WITH_DEBUG
1133 if( fDebugFile )
1134 *fDebugFile << "Fadc250Module::Decode:: "
1135 << "fSlot & FADC slot do not match AND data type != 0"
1136 << endl;
1137#endif
1138 return -1;
1139 }
1140
1141 // Acquire data objects depending on the data type defining word
1142 switch( data_type_def ) {
1143 case 0: // Block header, indicates the beginning of a block of events
1144 DecodeBlockHeader(data, data_type_id);
1145 break;
1146 case 1: // Block trailer, indicates the end of a block of events
1148 break;
1149 case 2: // Event header, indicates start of an event, includes the trigger number
1151 break;
1152 case 3: // Trigger time, time of trigger occurrence relative to the most recent global reset
1153 DecodeTriggerTime(data, data_type_id);
1154 break;
1155 case 4: // Window raw data
1156 DecodeWindowRawData(data, data_type_id);
1157 break;
1158 case 6: // Pulse raw data
1159 DecodePulseRawData(data, data_type_id);
1160 break;
1161 case 7: // Pulse integral
1163 break;
1164 case 8: // Pulse time
1166 break;
1167 case 9: // Pulse Parameters
1168 DecodePulseParameters(data, data_type_id);
1169 break;
1170 case 10: // Pulse Pedestal
1172 break;
1173 case 12: // Scaler header
1175 break;
1176 case 5: // Undefined type
1177 case 11: // Undefined type
1178 case 13: // Undefined type
1179 case 14: // Data not valid
1180 case 15: // Filler Word, should be ignored
1181 UnsupportedType(data, data_type_id);
1182 break;
1183 default:
1184 throw logic_error("Fadc250Module: incorrect masking of data_type_def");
1185 } // data_type_def switch
1186
1187#ifdef WITH_DEBUG
1188 if( fDebugFile )
1189 *fDebugFile << "**********************************************************************"
1190 << "\n" << endl;
1191#endif
1192
1193 return block_trailer_found;
1194}
1195
1196//_____________________________________________________________________________
1198{
1199 // Load THaSlotData
1200 for( uint32_t chan = 0; chan < NADCCHAN; chan++ ) {
1201 // Pulse Integral
1202 for( vsiz_t ievent = 0; ievent < fPulseData[chan].integral.size(); ievent++ )
1203 sldat->loadData("adc", chan, fPulseData[chan].integral[ievent], fPulseData[chan].integral[ievent]);
1204 // Pulse Time
1205 for( vsiz_t ievent = 0; ievent < fPulseData[chan].time.size(); ievent++ )
1206 sldat->loadData("adc", chan, fPulseData[chan].time[ievent], fPulseData[chan].time[ievent]);
1207 // Pulse Peak
1208 for( vsiz_t ievent = 0; ievent < fPulseData[chan].peak.size(); ievent++ )
1209 sldat->loadData("adc", chan, fPulseData[chan].peak[ievent], fPulseData[chan].peak[ievent]);
1210 // Pulse Pedestal
1211 for( vsiz_t ievent = 0; ievent < fPulseData[chan].pedestal.size(); ievent++ )
1212 sldat->loadData("adc", chan, fPulseData[chan].pedestal[ievent], fPulseData[chan].pedestal[ievent]);
1213 // Pulse Samples
1214 for( vsiz_t ievent = 0; ievent < fPulseData[chan].samples.size(); ievent++ )
1215 sldat->loadData("adc", chan, fPulseData[chan].samples[ievent], fPulseData[chan].samples[ievent]);
1216 } // Channel loop
1217}
1218
1219//_____________________________________________________________________________
1221 const UInt_t* pstop )
1222{
1223 // Load from evbuffer between [evbuffer,pstop]
1224
1225 return LoadSlot(sldat, evbuffer, 0, pstop + 1 - evbuffer);
1226}
1227
1228//_____________________________________________________________________________
1230 UInt_t pos, UInt_t len)
1231{
1232 // Load from bank data in evbuffer between [pos,pos+len)
1233
1234
1235 const auto* p = evbuffer + pos;
1236 const auto* q = p + len;
1237 while( p != q ) {
1238 if( Decode(p++) == 1 )
1239 break; // block trailer found
1240 }
1241
1242 LoadTHaSlotDataObj(sldat);
1243
1244 return fWordsSeen = p - (evbuffer + pos);
1245}
1246
1247//_____________________________________________________________________________
1248
1249} //namespace Decoder
1250
int Int_t
unsigned int UInt_t
uint32_t samples
std::vector< uint32_t > integral
uint32_t chan
bool Bool_t
const UInt_t kMaxUInt
const char Option_t
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
static const char *const here
Definition THaVar.cxx:64
float * q
#define NADCCHAN
virtual UInt_t GetData(Decoder::EModuleType mtype, UInt_t chan, UInt_t ievent) const
void DecodeWindowRawData(UInt_t pdat, uint32_t data_type_id)
void DecodePulseParameters(UInt_t pdat, uint32_t data_type_id)
std::vector< fadc_pulse_data > fPulseData
virtual UInt_t GetPulseCoarseTimeData(UInt_t chan, UInt_t ievent) const
virtual void CheckDecoderStatus() const
virtual Int_t GetFadcMode() const
virtual UInt_t GetNumFadcEvents(UInt_t chan) const
virtual UInt_t GetEmulatedPulseIntegralData(UInt_t chan) const
virtual std::vector< uint32_t > GetPulseSamplesVector(UInt_t chan) const
void DecodeEventHeader(UInt_t pdat)
void DecodeTriggerTime(UInt_t pdat, uint32_t data_type_id)
virtual void Clear(Option_t *opt="")
void DecodePulseTime(UInt_t pdat)
void DecodeBlockHeader(UInt_t pdat, uint32_t data_type_id)
virtual UInt_t GetPulsePedestalData(UInt_t chan, UInt_t ievent) const
virtual UInt_t GetOverflowBit(UInt_t chan, UInt_t ievent) const
void LoadTHaSlotDataObj(THaSlotData *sldat)
virtual UInt_t GetNumFadcSamples(UInt_t chan, UInt_t ievent) const
void PopulateDataVector(std::vector< uint32_t > &data_vector, uint32_t data) const
void DecodeBlockTrailer(UInt_t pdat)
virtual Int_t Decode(const UInt_t *data)
static const size_t NADCCHAN
void DecodePulseIntegral(UInt_t pdat)
virtual UInt_t GetNumEvents() const
void DecodePulsePedestal(UInt_t pdat)
void DecodeScalerHeader(UInt_t pdat)
virtual UInt_t GetPulseSamplesData(UInt_t chan, UInt_t ievent) const
virtual UInt_t GetPulseFineTimeData(UInt_t chan, UInt_t ievent) const
virtual UInt_t LoadSlot(THaSlotData *sldat, const UInt_t *evbuffer, const UInt_t *pstop)
virtual UInt_t GetPulseIntegralData(UInt_t chan, UInt_t ievent) const
static uint32_t SumVectorElements(const std::vector< uint32_t > &data_vector)
void DecodePulseRawData(UInt_t pdat, uint32_t data_type_id)
virtual Bool_t HasCapability(Decoder::EModuleType type)
void UnsupportedType(UInt_t pdat, uint32_t data_type_id)
virtual UInt_t GetTriggerTime() const
virtual UInt_t GetPedestalQuality(UInt_t chan, UInt_t ievent) const
virtual UInt_t GetPulseTimeData(UInt_t chan, UInt_t ievent) const
virtual UInt_t GetPulsePeakData(UInt_t chan, UInt_t ievent) const
static TypeIter_t fgThisType
virtual UInt_t GetUnderflowBit(UInt_t chan, UInt_t ievent) const
std::ofstream * fDebugFile
Definition Module.h:156
Int_t fFirmwareVers
Definition Module.h:153
UInt_t fSlot
Definition Module.h:142
UInt_t fWordsSeen
Definition Module.h:145
TypeSet_t::iterator TypeIter_t
Definition Module.h:40
Bool_t IsInit
Definition Module.h:151
virtual void Clear(Option_t *opt="")
Int_t loadData(const char *type, UInt_t chan, UInt_t dat, UInt_t raw)
TString fName
EModuleType
Definition Decoder.h:57
@ kPulseTime
Definition Decoder.h:57
@ kPulsePedestal
Definition Decoder.h:58
@ kCoarseTime
Definition Decoder.h:58
@ kFineTime
Definition Decoder.h:58
@ kSampleADC
Definition Decoder.h:57
@ kPulsePeak
Definition Decoder.h:58
@ kPulseIntegral
Definition Decoder.h:57
vector< int >::size_type vsiz_t
STL namespace.
ClassImp(TPyArg)