JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwADC18_Channel.h
Go to the documentation of this file.
1/**********************************************************\
2* File: QwADC18_Channel.h *
3* *
4* Author: P. M. King, W. Deconinck, B. Michaels *
5\**********************************************************/
6
7#pragma once
8
9// System headers
10#include <vector>
11
12// Qweak headers
13#include "VQwHardwareChannel.h"
14#include "MQwMockable.h"
15
16// Forward declarations
17class TTree;
18class QwBlinder;
19class QwParameterFile;
20#ifdef __USE_DATABASE__
22#endif
23
24/**
25 * \class QwADC18_Channel
26 * \ingroup QwAnalysis_ADC
27 * \brief Concrete hardware channel for HAPPEX 18-bit ADC modules
28 *
29 * Decodes and processes data from HAPPEX 18-bit ADC channels, providing
30 * access to raw and calibrated values, statistical moments, single-event
31 * cuts, and running statistics. Implements the dual-operator pattern for
32 * both type-specific and polymorphic operations.
33 */
35/****************************************************************//**
36 * Class: QwADC18_Channel
37 * Base class containing decoding functions for the HAPPEX 18-bit ADC
38 * The functions in this class will decode a single channel
39 * worth of ADC18_Channel data and provide the components
40 * through member functions.
41 ******************************************************************/
42 public:
43 static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex);
44 static void PrintErrorCounterHead();
45 static void PrintErrorCounterTail();
46
47 static const Double_t kTimePerSample;
48
50
56
59
60 public:
62 InitializeChannel("","");
63 SetADC18SaturationLimt(8.5);//FIXME set the default saturation limit
64 };
65 QwADC18_Channel(TString name, TString datatosave = "raw"): MQwMockable() {
66 InitializeChannel(name, datatosave);
67 SetADC18SaturationLimt(8.5);//FIXME set the default saturation limit
68 };
77 VQwHardwareChannel(value,datatosave), MQwMockable(value),
80 {
81 *this = value;
82 };
83 ~QwADC18_Channel() override { };
84
86
88 return new QwADC18_Channel(*this,datatosave);
89 };
90
91 /// \brief Initialize the fields in this object
92 void InitializeChannel(TString name, TString datatosave) override;
93
94 /// \brief Initialize the fields in this object
95 void InitializeChannel(TString subsystem, TString instrumenttype, TString name, TString datatosave) override;
96
97 void LoadChannelParameters(QwParameterFile &paramfile) override;
98
99 // Will update the default sample size for the module.
100 void SetDefaultSampleSize(size_t num_samples_map) {
101 // This will be checked against the no.of samples read by the module
102 fNumberOfSamples_map = num_samples_map;
103 };
104
105 void ClearEventData() override;
106
107 /// Internally generate random event data
108 void RandomizeEventData(int helicity = 0.0, double time = 0.0) override;
109
110 /// Forces the event "number of samples" variable to be what was expected from the mapfile.
111 /// NOTE: this should only be used in mock data generation!
113
114//------------------------------------------------------------------------------------------
115 void SmearByResolution(double resolution) override;
116//------------------------------------------------------------------------------------------
117
118 void SetEventData(Double_t value);
119 void SetRawEventData() override;
120
121 /// Encode the event data into a CODA buffer
122 void EncodeEventData(std::vector<UInt_t> &buffer) override;
123
124 /// Decode the event data from a CODA buffer
125 Bool_t IsHeaderWord(UInt_t word) const;
126 Int_t ProcessDataWord(UInt_t word);
127 Int_t ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t index = 0) override;
128
129 /// Process the event data according to pedestal and calibration factor
130 void ProcessEvent() override;
131
133 void AssignScaledValue(const QwADC18_Channel &value, Double_t scale);
134 void AssignValueFrom(const VQwDataElement* valueptr) override;
135 void AddValueFrom(const VQwHardwareChannel* valueptr) override;
136 void SubtractValueFrom(const VQwHardwareChannel* valueptr) override;
137 void MultiplyBy(const VQwHardwareChannel* valueptr) override;
138 void DivideBy(const VQwHardwareChannel* valueptr) override;
139 void ArcTan(const QwADC18_Channel &value);
140
144
145 VQwHardwareChannel& operator+=(const VQwHardwareChannel& input) override;
146 VQwHardwareChannel& operator-=(const VQwHardwareChannel& input) override;
147 VQwHardwareChannel& operator*=(const VQwHardwareChannel& input) override;
148 VQwHardwareChannel& operator/=(const VQwHardwareChannel& input) override;
149
150 const QwADC18_Channel operator+ (const QwADC18_Channel &value) const;
151 const QwADC18_Channel operator- (const QwADC18_Channel &value) const;
152 const QwADC18_Channel operator* (const QwADC18_Channel &value) const;
153 void Sum(const QwADC18_Channel &value1, const QwADC18_Channel &value2);
154 void Difference(const QwADC18_Channel &value1, const QwADC18_Channel &value2);
155 void Ratio(const QwADC18_Channel &numer, const QwADC18_Channel &denom);
156 void Product(const QwADC18_Channel &value1, const QwADC18_Channel &value2);
157 void DivideBy(const QwADC18_Channel& denom);
158 void AddChannelOffset(Double_t Offset);
159 void Scale(Double_t Offset) override;
160
161 /**
162 * Accumulate event values into the running sum with optional scaling.
163 * @param value Source channel to accumulate from.
164 * @param count Event count scaling (0 means use value.fGoodEventCount).
165 * @param ErrorMask Bit mask of error flags to exclude when accumulating.
166 */
167 void AccumulateRunningSum(const QwADC18_Channel& value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF);
168 void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override{
169 const QwADC18_Channel *tmp_ptr = dynamic_cast<const QwADC18_Channel*>(value);
170 if (tmp_ptr != NULL) {
171 AccumulateRunningSum(*tmp_ptr, count, ErrorMask);
172 } else {
173 throw std::invalid_argument("Standard exception from QwADC18_Channel::AccumulateRunningSum: incompatible hardware channel type");
174 }
175 };
176 ////deaccumulate one value from the running sum
177 inline void DeaccumulateRunningSum(const QwADC18_Channel& value, Int_t ErrorMask=0xFFFFFFF){
178 AccumulateRunningSum(value, -1, ErrorMask);
179 };
180 /*
181 void DeaccumulateRunningSum(VQwHardwareChannel *value){
182 const QwADC18_Channel *tmp_ptr = dynamic_cast<const QwADC18_Channel*>(value);
183 if (tmp_ptr != NULL) {
184 DeaccumulateRunningSum(*tmp_ptr);
185 } else {
186 throw std::invalid_argument("Standard exception from QwADC18_Channel::DeaccumulateRunningSum: incompatible hardware channel type");
187 }
188 };
189 */
190
191 void CalculateRunningAverage() override;
192
193 Bool_t MatchSequenceNumber(size_t seqnum);
194 Bool_t MatchNumberOfSamples(size_t numsamp);
195
196 /*Event cut related routines*/
197 //check values read from modules are at desired level
198 Bool_t ApplySingleEventCuts(Double_t LL, Double_t UL);
199 //check values read from modules are at desired level by comparing upper and lower limits (fULimit and fLLimit) set on this channel
200 Bool_t ApplySingleEventCuts() override;
201 // report number of events failed due to HW and event cut failure
202 void PrintErrorCounters() const override;
203
204 // FIXME Set the absolute staturation limit in volts
205 void SetADC18SaturationLimt(Double_t sat_volts = 8.5){
206 fSaturationABSLimit = sat_volts;
207 }
208
209 // Get the absolute staturation limit in volts
211 return fSaturationABSLimit;
212 }
213
214
215 // Check for hardware errors in the devices. This will return the device error code.
216 Int_t ApplyHWChecks() override;
217
218 // Update the error counters based on the internal fErrorFlag
219 void IncrementErrorCounters() override;
220
221 /*End*/
222
223 Int_t GetRawValue(size_t element) const override { return fValue_Raw; };
224 Double_t GetValue(size_t element) const override { return fValue; };
225 Double_t GetValueM2(size_t element) const override { return fValueM2; };
226 Double_t GetValueError(size_t element) const override { return fValueError; };
227
228 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
229 void FillHistograms() override;
230
231 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override;
232 void ConstructBranch(TTree *tree, TString &prefix) override;
233 void FillTreeVector(QwRootTreeBranchVector &values) const override;
234#ifdef HAS_RNTUPLE_SUPPORT
235 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) override;
236 void FillNTupleVector(std::vector<Double_t>& values) const override;
237#endif // HAS_RNTUPLE_SUPPORT
238
239 Double_t GetAverageVolts() const;
240
241 size_t GetSequenceNumber() const {return (fSequenceNumber);};
242 size_t GetNumberOfSamples() const {return (fNumberOfSamples);};
243
245
246 friend std::ostream& operator<< (std::ostream& stream, const QwADC18_Channel& channel);
247 void PrintValue() const override;
248 void PrintInfo() const override;
249
250 /// \brief Blind this channel as an asymmetry
251 void Blind(const QwBlinder *blinder);
252 /// \brief Blind this channel as a difference
253 void Blind(const QwBlinder *blinder, const QwADC18_Channel& yield);
254
255 void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override;
256
257#ifdef __USE_DATABASE__
258 // Error Counters exist in QwADC18_Channel, not in VQwHardwareChannel
259 //
260 void AddErrEntriesToList(std::vector<QwErrDBInterface> &row_list);
261#endif
262
263 protected:
265
266 private:
267 UInt_t fDiff_Raw;
268 UInt_t fBase_Raw;
269 UInt_t fPeak_Raw;
271
272 Double_t fValue;
273 Double_t fValueM2;
274 Double_t fValueError;
275
276 private:
277 static const Bool_t kDEBUG;
278 static const Int_t kHeaderWordsPerBank; //no.of header words per bank in the CODA buffer
279 static const Int_t kFooterWordsPerBank; //no.of footer words per bank in the CODA buffer
280 static const Int_t kHeaderWordsPerModule; //no.of header words per module in the CODA buffer
281 static const Int_t kFooterWordsPerModule; //no.of footer words per module in the CODA buffer
282 static const Int_t kDataWordsPerChannel; //no.of data words per channel in the CODA buffer
283 static const Int_t kMaxChannels; //no.of channels per module
284
285 private:
286 static const UInt_t mask31x; // = 0x80000000; // Header bit mask
287 static const UInt_t mask3029x; // = 0x60000000; // Channel number mask
288 static const UInt_t mask2625x; // = 0x06000000; // Divider value mask
289 static const UInt_t mask2422x; // = 0x01c00000; // Data type mask
290 static const UInt_t mask21x; // = 0x00200000; // Data type 0 value sign mask
291 static const UInt_t mask200x; // = 0x001fffff; // Data type 0 value field mask
292 static const UInt_t mask2118x; // = 0x003c0000; // Data types 1-2 sample number mask
293 static const UInt_t mask170x; // = 0x0003ffff; // Data types 1-2 value field mask
294 static const UInt_t mask150x; // = 0x0000ffff; // Data type 4 value field mask
295
296
297 /// Pointer to the running sum for this channel
299
300 /// Pointer to the DAC channel for this channel
302
303 /*! \name ADC Calibration */
304 // @{
305 static const Double_t kADC18_VoltsPerBit;
306 //@}
307
308
309 /*! \name Channel information data members */
310
311 /*! \name Channel configuration data members */
312 // @{
313
314 UInt_t fSequenceNumber; ///< Event sequence number for this channel
315 UInt_t fPreviousSequenceNumber; ///< Previous event sequence number for this channel
316 UInt_t fNumberOfSamples; ///< Number of samples read through the module
317 UInt_t fNumberOfSamples_map; ///< Number of samples in the expected to read through the module. This value is set in the QwBeamline map file
318
319 // Set of error counters for each HW test.
320 Int_t fErrorCount_HWSat; ///< check to see ADC channel is saturated
321 Int_t fErrorCount_sample; ///< for sample size check
322 Int_t fErrorCount_SW_HW; ///< HW_sum==SW_sum check
323 Int_t fErrorCount_Sequence; ///< sequence number check
324 Int_t fErrorCount_SameHW; ///< check to see ADC returning same HW value
325 Int_t fErrorCount_ZeroHW; ///< check to see ADC returning zero
326
327 Int_t fNumEvtsWithEventCutsRejected; ///< Counts the Event cut rejected events
328
329
330
331 Int_t fADC_Same_NumEvt; ///< Keep track of how many events with same ADC value returned
332 Int_t fSequenceNo_Prev; ///< Keep the sequence number of the last event
333 Int_t fSequenceNo_Counter; ///< Internal counter to keep track of the sequence number
334 Double_t fPrev_HardwareBlockSum; ///< Previous Module-based sum of the four sub-blocks
335
336
337
338 Double_t fSaturationABSLimit;///<absolute value of the ADC18 saturation volt
339
340
341 const static Bool_t bDEBUG=kFALSE;///<debugging display purposes
342
343 ///<For ADC18 data element trimming uses
344 Bool_t bHw_sum;
346 Bool_t bBlock;
351
352};
virtual void LoadMockDataParameters(QwParameterFile &paramfile)
Load the mock data parameters from the current line in the param file.
Concrete hardware channel for HAPPEX 18-bit ADC modules.
static const Int_t kFooterWordsPerModule
static const Int_t kMaxChannels
Double_t GetValue(size_t element) const override
static const Int_t kHeaderWordsPerBank
Double_t fPrev_HardwareBlockSum
Previous Module-based sum of the four sub-blocks.
UInt_t fNumberOfSamples
Number of samples read through the module.
static const Int_t kHeaderWordsPerModule
static const Bool_t kDEBUG
Int_t fErrorCount_HWSat
check to see ADC channel is saturated
void PrintInfo() const override
Print multiple lines of information about this data element.
QwADC18_Channel & operator*=(const QwADC18_Channel &value)
static const Double_t kADC18_VoltsPerBit
void DeaccumulateRunningSum(const QwADC18_Channel &value, Int_t ErrorMask=0xFFFFFFF)
void AddValueFrom(const VQwHardwareChannel *valueptr) override
Int_t ApplyHWChecks() override
Int_t fErrorCount_SW_HW
HW_sum==SW_sum check.
void RandomizeEventData(int helicity=0.0, double time=0.0) override
Internally generate random event data.
void ProcessEvent() override
Process the event data according to pedestal and calibration factor.
static void PrintErrorCounterHead()
Double_t fSaturationABSLimit
absolute value of the ADC18 saturation volt
static const UInt_t mask3029x
static const Int_t kFooterWordsPerBank
const QwADC18_Channel operator+(const QwADC18_Channel &value) const
Int_t fSequenceNo_Counter
Internal counter to keep track of the sequence number.
Int_t fErrorCount_SameHW
check to see ADC returning same HW value
void FillHistograms() override
Fill the histograms for this data element.
void Product(const QwADC18_Channel &value1, const QwADC18_Channel &value2)
void SubtractValueFrom(const VQwHardwareChannel *valueptr) override
Int_t fErrorCount_sample
for sample size check
Int_t fSequenceNo_Prev
Keep the sequence number of the last event.
QwADC18_Channel(TString name, TString datatosave="raw")
void ArcTan(const QwADC18_Channel &value)
Double_t GetADC18SaturationLimt()
void Sum(const QwADC18_Channel &value1, const QwADC18_Channel &value2)
void ConstructHistograms(TDirectory *folder, TString &prefix) override
Construct the histograms for this data element.
void Scale(Double_t Offset) override
void CalculateRunningAverage() override
size_t GetSequenceNumber() const
void Difference(const QwADC18_Channel &value1, const QwADC18_Channel &value2)
void PrintValue() const override
Print single line of value and error of this data element.
void FillTreeVector(QwRootTreeBranchVector &values) const override
Int_t fErrorCount_ZeroHW
check to see ADC returning zero
void ConstructBranch(TTree *tree, TString &prefix) override
static const UInt_t mask2422x
void AddChannelOffset(Double_t Offset)
void PrintErrorCounters() const override
report number of events failed due to HW and event cut failure
QwADC18_Channel * fDAC
Pointer to the DAC channel for this channel.
void AccumulateRunningSum(const QwADC18_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
static const UInt_t mask2625x
void Ratio(const QwADC18_Channel &numer, const QwADC18_Channel &denom)
static const UInt_t mask170x
Int_t fADC_Same_NumEvt
Keep track of how many events with same ADC value returned.
void SetRawEventData() override
void ClearEventData() override
Clear the event data in this element.
void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override
void AssignScaledValue(const QwADC18_Channel &value, Double_t scale)
const QwADC18_Channel operator-(const QwADC18_Channel &value) const
QwADC18_Channel & operator=(const QwADC18_Channel &value)
Bool_t MatchNumberOfSamples(size_t numsamp)
void ForceMapfileSampleSize()
static const Bool_t bDEBUG
debugging display purposes
QwADC18_Channel * fRunningSum
Pointer to the running sum for this channel.
QwADC18_Channel(const QwADC18_Channel &value)
Double_t GetValueError(size_t element) const override
void SmearByResolution(double resolution) override
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
Int_t GetRawValue(size_t element) const override
void SetADC18SaturationLimt(Double_t sat_volts=8.5)
void SetEventData(Double_t value)
void LoadChannelParameters(QwParameterFile &paramfile) override
UInt_t fPreviousSequenceNumber
Previous event sequence number for this channel.
Int_t fNumEvtsWithEventCutsRejected
Counts the Event cut rejected events.
QwADC18_Channel & operator-=(const QwADC18_Channel &value)
static const UInt_t mask200x
static const UInt_t mask31x
VQwHardwareChannel & operator/=(const VQwHardwareChannel &input) override
~QwADC18_Channel() override
Double_t GetAverageVolts() const
void IncrementErrorCounters() override
void Blind(const QwBlinder *blinder)
Blind this channel as an asymmetry.
static const Double_t kTimePerSample
void MultiplyBy(const VQwHardwareChannel *valueptr) override
Bool_t ApplySingleEventCuts() override
QwADC18_Channel & operator+=(const QwADC18_Channel &value)
static const UInt_t mask2118x
static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex)
Bool_t IsHeaderWord(UInt_t word) const
Decode the event data from a CODA buffer.
void InitializeChannel(TString name, TString datatosave) override
Initialize the fields in this object.
static void PrintErrorCounterTail()
QwADC18_Channel(const QwADC18_Channel &value, VQwDataElement::EDataToSave datatosave)
static const UInt_t mask21x
static const Int_t kDataWordsPerChannel
VQwHardwareChannel * Clone(VQwDataElement::EDataToSave datatosave) const override
Bool_t MatchSequenceNumber(size_t seqnum)
UInt_t fNumberOfSamples_map
Number of samples in the expected to read through the module. This value is set in the QwBeamline map...
size_t GetNumberOfSamples() const
Double_t GetValueM2(size_t element) const override
void AssignValueFrom(const VQwDataElement *valueptr) override
Int_t fErrorCount_Sequence
sequence number check
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t index=0) override
Process the CODA event buffer for this element.
void EncodeEventData(std::vector< UInt_t > &buffer) override
Encode the event data into a CODA buffer.
void DivideBy(const VQwHardwareChannel *valueptr) override
friend std::ostream & operator<<(std::ostream &stream, const QwADC18_Channel &channel)
void SetDefaultSampleSize(size_t num_samples_map)
const QwADC18_Channel operator*(const QwADC18_Channel &value) const
static const UInt_t mask150x
Int_t ProcessDataWord(UInt_t word)
void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
UInt_t fSequenceNumber
Event sequence number for this channel.
Configuration file parser with flexible tokenization and search capabilities.
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:53
The pure virtual base class of all data elements.
Abstract base for concrete hardware channels implementing dual-operator pattern.
Double_t GetValueError() const
Double_t GetValueM2() const
void SetCalibrationFactor(Double_t factor)
virtual VQwHardwareChannel * Clone() const
virtual void AddErrEntriesToList(std::vector< QwErrDBInterface > &)
Double_t GetValueWidth() const
Double_t GetValue() const
virtual void DeaccumulateRunningSum(const VQwHardwareChannel *value, Int_t ErrorMask=0xFFFFFFF)
virtual void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
Data blinding utilities for parity violation analysis.
Definition QwBlinder.h:57