Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaSlotData.h
Go to the documentation of this file.
1#ifndef Podd_THaSlotData_h_
2#define Podd_THaSlotData_h_
3
5//
6// THaSlotData
7// Data in one slot of one crate from DAQ.
8//
9// THaSlotData contains data from one slot of one crate
10// from a CODA event. Public methods allow to obtain
11// raw data for this crate & slot, or to obtain TDC, ADC,
12// or scaler data for each channel in the slot. Methods
13// clearEvent() and loadData() should only be used by
14// the decoder. WARNING: For efficiency, only the
15// hit counters are zero'd each event, not the data
16// arrays, see below.
17//
18// author Robert Michaels (rom@jlab.org)
19//
21
22#include "Decoder.h"
23#include "Module.h"
24#include "CustomAlloc.h"
25#include <cassert>
26#include <fstream>
27#include <string>
28#include <vector>
29#include <memory>
30
31const int SD_WARN = -2;
32const int SD_ERR = -1;
33const int SD_OK = 1;
34
35namespace Decoder {
36
38
39public:
40
41 static const UInt_t DEFNCHAN; // Default number of channels
42 static const UInt_t DEFNDATA; // Default number of data words
43 static const UInt_t DEFNHITCHAN; // Default number of hits per channel
44
47 virtual ~THaSlotData();
48 const char* devType() const; // "adc", "tdc", "scaler"
49 Int_t loadModule( const THaCrateMap* map );
50 UInt_t getNumRaw() const { return numraw; }; // Amount of raw CODA data
51 UInt_t getRawData(UInt_t ihit) const; // Returns raw data words
53 UInt_t getNumHits(UInt_t chan) const; // Num hits on a channel
54 UInt_t getNumChan() const; // Num unique channels hit
55 UInt_t getNextChan(UInt_t index) const; // List of unique channels hit
56 UInt_t getData(UInt_t chan, UInt_t hit) const;// Data (adc,tdc,scaler) on 1 chan
57 UInt_t getCrate() const { return crate; }
58 UInt_t getSlot() const { return slot; }
59 UInt_t getNchan() const { return fNchan; }
60 void clearEvent(); // clear event counters
61 Int_t loadData( const char* type, UInt_t chan, UInt_t dat, UInt_t raw );
63
64 // new
65 UInt_t LoadIfSlot( const UInt_t* evbuffer, const UInt_t* pstop );
66 UInt_t LoadBank( const UInt_t* p, UInt_t pos, UInt_t len );
68 Bool_t IsMultiBlockMode() { if (fModule) return fModule->IsMultiBlockMode(); return false; };
69 Bool_t BlockIsDone() { if (fModule) return fModule->BlockIsDone(); return false; };
70
71 void SetDebugFile(std::ofstream *file) { fDebugFile = file; };
72 Module* GetModule() { return fModule.get(); };
73
74 // Define crate, slot
76 UInt_t ndata = DEFNDATA, UInt_t nhitperchan = DEFNHITCHAN );
77 void print() const;
78 void print_to_file() const;
79 void compressdataindex(UInt_t numidx);
80
81private:
82
85 std::string device;
86 std::unique_ptr<Module> fModule;
87 UInt_t numhitperchan; // expected number of hits per channel
88 UInt_t numraw; // Hit counters (numraw, numHits, numchanhit)
89 UInt_t numchanhit; // can be zero'd by clearEvent each event.
90 UInt_t firstfreedataidx; // pointer to first free space in dataindex array
92 VectorUInt numHits; // numHits[channel]
93 VectorUIntNI chanlist; // chanlist[hitindex]
94 VectorUIntNI idxlist; // [channel] pointer to 1st entry in dataindex
95 VectorUIntNI chanindex; // [channel] gives hitindex
96 VectorUIntNI dataindex; // [idxlist] pointer to rawdata and data
97 VectorUIntNI numMaxHits; // [channel] current maximum number of hits
98 VectorUIntNI rawData; // rawData[hit] (all bits)
99 VectorUIntNI data; // data[hit] (only data bits)
100 std::ofstream *fDebugFile; // debug output to this file, if nonzero
101 bool didini; // true if object initialized via define()
102 UInt_t fNchan; // Number of channels for this device
103
104 void compressdataindexImpl(UInt_t numidx);
105
106 ClassDef(THaSlotData,0) // Data in one slot of fastbus, vme, camac
107};
108
109//______________ inline functions _____________________________________________
111 // Returns the raw data (all bits)
112 assert( hit < numraw );
113 if (hit < numraw) return rawData[hit];
114 return 0;
115}
116
117//_____________________________________________________________________________
118// Data (words on 1 chan)
119inline
121 assert(chan < fNchan && hit < numHits[chan] );
122 if ( chan >= fNchan || numHits[chan] <= hit)
123 return 0;
125 assert(index < numraw);
126 if (index < numraw) return rawData[index];
127 return 0;
128}
129
130//_____________________________________________________________________________
131inline
133 // Num hits on a channel
134 assert(chan < fNchan );
135 return (chan < fNchan) ? numHits[chan] : 0;
136}
137
138//_____________________________________________________________________________
139inline
141 // Num unique channels # hit
142 return numchanhit;
143}
144
145//_____________________________________________________________________________
146inline
148 // List of unique channels hit
149 assert(index < numchanhit);
150 if (index < numchanhit)
151 return chanlist[index];
152 return 0;
153}
154
155//_____________________________________________________________________________
156// Data (words on 1 chan)
157inline
159 assert(chan < fNchan && hit < numHits[chan] );
160 if ( chan >= fNchan || numHits[chan] <= hit)
161 return 0;
163 assert(index < numraw);
164 if (index < numraw) return data[index];
165 return 0;
166}
167
168//_____________________________________________________________________________
169// Device type (adc, tdc, scaler)
170inline
171const char* THaSlotData::devType() const {
172 return device.c_str();
173}
174
175//_____________________________________________________________________________
176inline
178 // Only the minimum is cleared; e.g. data array is not cleared.
179 // CAUTION: this code is critical for performance
180 numraw = 0;
183 while( numchanhit>0 ) numHits[chanlist[--numchanhit]] = 0;
184}
185
186//_____________________________________________________________________________
187inline
189
190 if( firstfreedataidx+numidx >= dataindex.size() )
191 compressdataindexImpl(numidx);
192}
193
194}
195
196#endif
int Int_t
unsigned int UInt_t
uint32_t chan
bool Bool_t
#define ClassDef(name, id)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
const int SD_WARN
Definition THaSlotData.h:31
const int SD_OK
Definition THaSlotData.h:33
const int SD_ERR
Definition THaSlotData.h:32
UInt_t LoadBank(const UInt_t *p, UInt_t pos, UInt_t len)
UInt_t getSlot() const
Definition THaSlotData.h:58
UInt_t getNumChan() const
UInt_t getNumRaw() const
Definition THaSlotData.h:50
UInt_t getCrate() const
Definition THaSlotData.h:57
std::unique_ptr< Module > fModule
Definition THaSlotData.h:86
VectorUIntNI dataindex
Definition THaSlotData.h:96
Int_t loadData(const char *type, UInt_t chan, UInt_t dat, UInt_t raw)
UInt_t getNumHits(UInt_t chan) const
void print_to_file() const
VectorUIntNI chanlist
Definition THaSlotData.h:93
static const UInt_t DEFNHITCHAN
Definition THaSlotData.h:43
VectorUIntNI numMaxHits
Definition THaSlotData.h:97
UInt_t getNchan() const
Definition THaSlotData.h:59
VectorUIntNI data
Definition THaSlotData.h:99
void SetDebugFile(std::ofstream *file)
Definition THaSlotData.h:71
VectorUIntNI idxlist
Definition THaSlotData.h:94
UInt_t getRawData(UInt_t ihit) const
UInt_t getNextChan(UInt_t index) const
void compressdataindex(UInt_t numidx)
Bool_t IsMultiBlockMode()
Definition THaSlotData.h:68
void define(UInt_t crate, UInt_t slot, UInt_t nchan=DEFNCHAN, UInt_t ndata=DEFNDATA, UInt_t nhitperchan=DEFNHITCHAN)
std::ofstream * fDebugFile
VectorUIntNI chanindex
Definition THaSlotData.h:95
static const UInt_t DEFNDATA
Definition THaSlotData.h:42
UInt_t getData(UInt_t chan, UInt_t hit) const
Int_t loadModule(const THaCrateMap *map)
void compressdataindexImpl(UInt_t numidx)
VectorUIntNI rawData
Definition THaSlotData.h:98
const char * devType() const
UInt_t LoadIfSlot(const UInt_t *evbuffer, const UInt_t *pstop)
static const UInt_t DEFNCHAN
Definition THaSlotData.h:41
std::vector< UInt_t > VectorUInt
Definition CustomAlloc.h:36
std::vector< UInt_t, default_init_allocator< UInt_t > > VectorUIntNI
Definition CustomAlloc.h:38