JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwVQWK_Channel.h
Go to the documentation of this file.
1
2/*!
3 * \file QwVQWK_Channel.h
4 * \brief Decoding and management for VQWK ADC channels (6x32-bit datawords)
5 * \author P. M. King
6 * \date 2007-05-08
7 */
8
9#pragma once
10
11// System headers
12#include <vector>
13
14// ROOT headers
15#include "TTree.h"
16
17// RNTuple headers
18#ifdef HAS_RNTUPLE_SUPPORT
19#include "ROOT/RNTupleModel.hxx"
20#endif // HAS_RNTUPLE_SUPPORT
21
22// Qweak headers
23#include "VQwHardwareChannel.h"
24#include "MQwMockable.h"
25
26// Forward declarations
27class QwBlinder;
28class QwParameterFile;
29#ifdef __USE_DATABASE__
31#endif
32
33///
34/// \ingroup QwAnalysis_ADC
35///
36/// \ingroup QwAnalysis_BL
37/**
38 * \class QwVQWK_Channel
39 * \ingroup QwAnalysis_ADC
40 * \brief Concrete hardware channel for VQWK ADC modules (6x32-bit words)
41 *
42 * Decodes and processes the data for a single VQWK channel, providing access
43 * to sub-blocks and hardware sums, statistical moments, single-event cuts,
44 * and running statistics. Implements the dual-operator pattern by offering
45 * efficient type-specific operators and delegating polymorphic operators from
46 * the VQwHardwareChannel interface using dynamic_cast.
47 */
49/****************************************************************//**
50 * Class: QwVQWK_Channel
51 * Base class containing decoding functions for the VQWK_Channel
52 * 6 32-bit datawords.
53 * The functions in this class will decode a single channel
54 * worth of VQWK_Channel data and provide the components
55 * through member functions.
56 ******************************************************************/
57 public:
58 static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex);
59 static void PrintErrorCounterHead();
60 static void PrintErrorCounterTail();
61
62 static const Double_t kTimePerSample;
63
65
71
74
75 public:
77 InitializeChannel("","");
78 SetVQWKSaturationLimt(8.5);//set the default saturation limit
79 };
80 QwVQWK_Channel(TString name, TString datatosave = "raw"): MQwMockable() {
81 InitializeChannel(name, datatosave);
82 SetVQWKSaturationLimt(8.5);//set the default saturation limit
83 };
93 VQwHardwareChannel(value,datatosave), MQwMockable(value),
97 {
98 *this = value;
99 };
100 ~QwVQWK_Channel() override { };
101
109
111
113 return new QwVQWK_Channel(*this,datatosave);
114 };
115
116 /// \brief Initialize the fields in this object
117 void InitializeChannel(TString name, TString datatosave) override;
118
119 /// \brief Initialize the fields in this object
120 void InitializeChannel(TString subsystem, TString instrumenttype, TString name, TString datatosave) override;
121
122 void LoadChannelParameters(QwParameterFile &paramfile) override;
123
124 // Will update the default sample size for the module.
125 void SetDefaultSampleSize(size_t num_samples_map) {
126 // This will be checked against the no.of samples read by the module
127 fNumberOfSamples_map = num_samples_map;
128 };
129
130 void ClearEventData() override;
131
132 /// Internally generate random event data
133 void RandomizeEventData(int helicity = 0.0, double time = 0.0) override;
134
135 /// Forces the event "number of samples" variable to be what was expected from the mapfile.
136 /// NOTE: this should only be used in mock data generation!
138
139//------------------------------------------------------------------------------------------
140 void SmearByResolution(double resolution) override;
141//------------------------------------------------------------------------------------------
142
143 /// TODO: SetHardwareSum should be removed, and SetEventData
144 /// should be made protected.
145 void SetHardwareSum(Double_t hwsum, UInt_t sequencenumber = 0);
146 void SetEventData(Double_t* block, UInt_t sequencenumber = 0);
147 void SetRawEventData() override;
148
149 /// Encode the event data into a CODA buffer
150 void EncodeEventData(std::vector<UInt_t> &buffer) override;
151 /// Decode the event data from a CODA buffer
152 Int_t ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t index = 0) override;
153 /// Process the event data according to pedestal and calibration factor
154 void ProcessEvent() override;
155
156
158 void AssignScaledValue(const QwVQWK_Channel &value, Double_t scale);
159 void AssignValueFrom(const VQwDataElement* valueptr) override;
160 void AddValueFrom(const VQwHardwareChannel* valueptr) override;
161 void SubtractValueFrom(const VQwHardwareChannel* valueptr) override;
162 void MultiplyBy(const VQwHardwareChannel* valueptr) override;
163 void DivideBy(const VQwHardwareChannel* valueptr) override;
164 void ArcTan(const QwVQWK_Channel &value);
165
169
170 VQwHardwareChannel& operator+=(const VQwHardwareChannel& input) override;
171 VQwHardwareChannel& operator-=(const VQwHardwareChannel& input) override;
172 VQwHardwareChannel& operator*=(const VQwHardwareChannel& input) override;
173 VQwHardwareChannel& operator/=(const VQwHardwareChannel& input) override;
174
175 const QwVQWK_Channel operator+ (const QwVQWK_Channel &value) const;
176 const QwVQWK_Channel operator- (const QwVQWK_Channel &value) const;
177 const QwVQWK_Channel operator* (const QwVQWK_Channel &value) const;
178 void Sum(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2);
179 void Difference(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2);
180 void Ratio(const QwVQWK_Channel &numer, const QwVQWK_Channel &denom);
181 void Product(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2);
182 void DivideBy(const QwVQWK_Channel& denom);
183 void AddChannelOffset(Double_t Offset);
184 void Scale(Double_t Offset) override;
185
186
187 /**
188 * Accumulate event values into the running sum with optional scaling.
189 * @param value Source channel to accumulate from.
190 * @param count Event count scaling (0 means use value.fGoodEventCount).
191 * @param ErrorMask Bit mask of error flags to exclude when accumulating.
192 */
193 void AccumulateRunningSum(const QwVQWK_Channel& value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF);
194 void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override{
195 const QwVQWK_Channel *tmp_ptr = dynamic_cast<const QwVQWK_Channel*>(value);
196 if (tmp_ptr != NULL) {
197 AccumulateRunningSum(*tmp_ptr, count, ErrorMask);
198 } else {
199 throw std::invalid_argument("Standard exception from QwVQWK_Channel::AccumulateRunningSum: incompatible hardware channel type");
200 }
201 };
202 ////deaccumulate one value from the running sum
203 inline void DeaccumulateRunningSum(const QwVQWK_Channel& value, Int_t ErrorMask=0xFFFFFFF){
204 AccumulateRunningSum(value, -1, ErrorMask);
205 };
206 /*
207 void DeaccumulateRunningSum(VQwHardwareChannel *value){
208 const QwVQWK_Channel *tmp_ptr = dynamic_cast<const QwVQWK_Channel*>(value);
209 if (tmp_ptr != NULL) DeaccumulateRunningSum(*tmp_ptr);
210 };
211 */
212
213 void CalculateRunningAverage() override;
214
215 Bool_t MatchSequenceNumber(size_t seqnum);
216 Bool_t MatchNumberOfSamples(size_t numsamp);
217
218 /*Event cut related routines*/
219 Bool_t ApplySingleEventCuts(Double_t LL,Double_t UL);//check values read from modules are at desired level
220 Bool_t ApplySingleEventCuts() override;//check values read from modules are at desired level by comparing upper and lower limits (fULimit and fLLimit) set on this channel
221 void PrintErrorCounters() const override;// report number of events failed due to HW and event cut failure
222
223 void SetVQWKSaturationLimt(Double_t sat_volts=8.5){//Set the absolute staturation limit in volts.
224 fSaturationABSLimit=sat_volts;
225 }
226
227 Double_t GetVQWKSaturationLimt(){//Get the absolute staturation limit in volts.
228 return fSaturationABSLimit;
229 }
230
231
232 Int_t ApplyHWChecks() override; //Check for hardware errors in the devices. This will return the device error code.
233
234 void IncrementErrorCounters() override;//update the error counters based on the internal fErrorFlag
235
236 /*End*/
237
238 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
239 void FillHistograms() override;
240
241 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override;
242 void ConstructBranch(TTree *tree, TString &prefix) override;
243 void FillTreeVector(QwRootTreeBranchVector &values) const override;
244
245#ifdef HAS_RNTUPLE_SUPPORT
246 // RNTuple support methods
247 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString &prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>> &fieldPtrs) override;
248 void FillNTupleVector(std::vector<Double_t>& values) const override;
249#endif // HAS_RNTUPLE_SUPPORT
250
251 Int_t GetRawValue(size_t element) const override {
252 RangeCheck(element);
253 if (element==0) return fHardwareBlockSum_raw;
254 return fBlock_raw[element-1];
255 }
256 Double_t GetValue(size_t element) const override {
257 RangeCheck(element);
258 if (element==0) return fHardwareBlockSum;
259 return fBlock[element-1];
260 }
261 Double_t GetValueM2(size_t element) const override {
262 RangeCheck(element);
263 if (element==0) return fHardwareBlockSumM2;
264 return fBlockM2[element-1];
265 }
266 Double_t GetValueError(size_t element) const override {
267 RangeCheck(element);
268 if (element==0) return fHardwareBlockSumError;
269 return fBlockError[element-1];
270 }
271
272
273 Double_t GetAverageVolts() const;
274
275 size_t GetSequenceNumber() const {return (fSequenceNumber);};
276 size_t GetNumberOfSamples() const {return (fNumberOfSamples);};
277
279
280 friend std::ostream& operator<< (std::ostream& stream, const QwVQWK_Channel& channel);
281 void PrintValue() const override;
282 void PrintInfo() const override;
283
284 /// \brief Blind this channel as an asymmetry
285 void Blind(const QwBlinder *blinder);
286 /// \brief Blind this channel as a difference
287 void Blind(const QwBlinder *blinder, const QwVQWK_Channel& yield);
288
289 void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override;
290
291#ifdef __USE_DATABASE__
292 // Error Counters exist in QwVQWK_Channel, not in VQwHardwareChannel
293 //
294 void AddErrEntriesToList(std::vector<QwErrDBInterface> &row_list);
295#endif
296
297 protected:
299
300private:
301 // The following specific access methods should only be used internally,
302 // if at all.
303 Double_t GetBlockValue(size_t blocknum) const { return GetValue(blocknum+1);};
304 Double_t GetBlockErrorValue(size_t blocknum) const { return GetValueError(blocknum+1);};
305
306 Double_t GetHardwareSum() const { return GetValue(0);};
307 Double_t GetHardwareSumM2() const { return GetValueM2(0); };
308 Double_t GetHardwareSumWidth() const { return GetValueWidth(0); };
309 Double_t GetHardwareSumError() const { return GetValueError(0); };
310 // Double_t GetSoftwareSum() const {return fSoftwareBlockSum;};
311
312 Int_t GetRawBlockValue(size_t blocknum) const {return GetRawValue(blocknum+1);};
313 Int_t GetRawHardwareSum() const { return GetRawValue(0);};
315
316 private:
317 static const Bool_t kDEBUG;
318 static const Int_t kWordsPerChannel; //no.of words per channel in the CODA buffer
319 static const Int_t kMaxChannels; //no.of channels per module
320
321 /*! \name ADC Calibration */
322 // @{
323 static const Double_t kVQWK_VoltsPerBit;
324 //@}
325
326
327 /*! \name Channel information data members */
328
329 /*! \name Channel configuration data members */
330 // @{
331
332 //UInt_t fBlocksPerEvent;
334 // @}
335
336
337 /*! \name Event data members---Raw values */
338 // @{
339 Int_t fBlock_raw[4]; ///< Array of the sub-block data as read from the module
340 Int_t fHardwareBlockSum_raw; ///< Module-based sum of the four sub-blocks as read from the module
341 Int_t fSoftwareBlockSum_raw; ///< Sum of the data in the four sub-blocks raw
342 /*! \name Event data members---Potentially calibrated values*/
343 // @{
344 // The following values potentially have pedestal removed and calibration applied
345 Double_t fBlock[4]; ///< Array of the sub-block data
346 Double_t fHardwareBlockSum; ///< Module-based sum of the four sub-blocks
347 // @}
348
349
350 /// \name Calculation of the statistical moments
351 // @{
352 // Moments of the separate blocks
353 Double_t fBlockM2[4]; ///< Second moment of the sub-block
354 Double_t fBlockError[4]; ///< Uncertainty on the sub-block
355 // Moments of the hardware sum
356 Double_t fHardwareBlockSumM2; ///< Second moment of the hardware sum
357 Double_t fHardwareBlockSumError; ///< Uncertainty on the hardware sum
358 // @}
359
360
361 UInt_t fSequenceNumber; ///< Event sequence number for this channel
362 UInt_t fPreviousSequenceNumber; ///< Previous event sequence number for this channel
363 UInt_t fNumberOfSamples; ///< Number of samples read through the module
364 UInt_t fNumberOfSamples_map; ///< Number of samples in the expected to read through the module. This value is set in the QwBeamline map file
365
366
367 // Set of error counters for each HW test.
368 Int_t fErrorCount_HWSat; ///< check to see ADC channel is saturated
369 Int_t fErrorCount_sample; ///< for sample size check
370 Int_t fErrorCount_SW_HW; ///< HW_sum==SW_sum check
371 Int_t fErrorCount_Sequence; ///< sequence number check
372 Int_t fErrorCount_SameHW; ///< check to see ADC returning same HW value
373 Int_t fErrorCount_ZeroHW; ///< check to see ADC returning zero
374
375 Int_t fNumEvtsWithEventCutsRejected; ///< Counts the Event cut rejected events
376
377
378
379
380
381
382 Int_t fADC_Same_NumEvt; ///< Keep track of how many events with same ADC value returned
383 Int_t fSequenceNo_Prev; ///< Keep the sequence number of the last event
384 Int_t fSequenceNo_Counter; ///< Internal counter to keep track of the sequence number
385 Double_t fPrev_HardwareBlockSum; ///< Previous Module-based sum of the four sub-blocks
386
387
388
389 Double_t fSaturationABSLimit;///<absolute value of the VQWK saturation volt
390
391
392 const static Bool_t bDEBUG=kFALSE;///<debugging display purposes
393
394 ///<For VQWK data element trimming uses
395 Bool_t bHw_sum;
397 Bool_t bBlock;
402
403private:
404
405
406
407
408
409
410};
virtual void LoadMockDataParameters(QwParameterFile &paramfile)
Load the mock data parameters from the current line in the param file.
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
Concrete hardware channel for VQWK ADC modules (6x32-bit words)
Int_t fADC_Same_NumEvt
Keep track of how many events with same ADC value returned.
Short_t fBlocksPerEvent
void SubtractValueFrom(const VQwHardwareChannel *valueptr) override
Double_t GetValueError() const
Double_t fBlockError[4]
Uncertainty on the sub-block.
QwVQWK_Channel & operator=(const QwVQWK_Channel &value)
void InitializeChannel(TString name, TString datatosave) override
Initialize the fields in this object.
Double_t GetValueError(size_t element) const override
void SetCalibrationToVolts()
void SetVQWKSaturationLimt(Double_t sat_volts=8.5)
Double_t GetHardwareSumM2() const
~QwVQWK_Channel() override
void Sum(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
void SetRawEventData() override
Double_t GetValueM2() const
Bool_t MatchNumberOfSamples(size_t numsamp)
Double_t GetVQWKSaturationLimt()
UInt_t fNumberOfSamples
Number of samples read through the module.
Int_t GetRawBlockValue(size_t blocknum) const
Int_t GetRawSoftwareSum() const
static void PrintErrorCounterTail()
Int_t fSequenceNo_Prev
Keep the sequence number of the last event.
void MultiplyBy(const VQwHardwareChannel *valueptr) override
void SetHardwareSum(Double_t hwsum, UInt_t sequencenumber=0)
Bool_t MatchSequenceNumber(size_t seqnum)
VQwHardwareChannel * Clone(VQwDataElement::EDataToSave datatosave) const override
void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override
UInt_t fSequenceNumber
Event sequence number for this channel.
Double_t fSaturationABSLimit
absolute value of the VQWK saturation volt
static const Bool_t kDEBUG
Bool_t bDevice_Error_Code
VQwHardwareChannel & operator/=(const VQwHardwareChannel &input) override
void ConstructHistograms(TDirectory *folder, TString &prefix) override
Construct the histograms for this data element.
void DeaccumulateRunningSum(const QwVQWK_Channel &value, Int_t ErrorMask=0xFFFFFFF)
Double_t GetHardwareSumWidth() const
void SetDefaultSampleSize(size_t num_samples_map)
void ClearEventData() override
Clear the event data in this element.
void CalculateRunningAverage() override
QwVQWK_Channel & operator*=(const QwVQWK_Channel &value)
const QwVQWK_Channel operator-(const QwVQWK_Channel &value) const
void LoadChannelParameters(QwParameterFile &paramfile) override
Int_t fHardwareBlockSum_raw
Module-based sum of the four sub-blocks as read from the module.
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
void AddChannelOffset(Double_t Offset)
void PrintErrorCounters() const override
report number of events failed due to HW and event cut failure
void FillHistograms() override
Fill the histograms for this data element.
void ForceMapfileSampleSize()
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t index=0) override
Decode the event data from a CODA buffer.
void AccumulateRunningSum(const QwVQWK_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
void ConstructBranch(TTree *tree, TString &prefix) override
Int_t fSoftwareBlockSum_raw
Sum of the data in the four sub-blocks raw.
Double_t GetHardwareSumError() const
Double_t GetValue(size_t element) const override
static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex)
Int_t fErrorCount_ZeroHW
check to see ADC returning zero
Int_t fErrorCount_sample
for sample size check
QwVQWK_Channel & operator+=(const QwVQWK_Channel &value)
void RandomizeEventData(int helicity=0.0, double time=0.0) override
Internally generate random event data.
static const Bool_t bDEBUG
debugging display purposes
Double_t GetBlockErrorValue(size_t blocknum) const
Double_t GetHardwareSum() const
UInt_t fPreviousSequenceNumber
Previous event sequence number for this channel.
size_t GetSequenceNumber() const
Double_t fHardwareBlockSumError
Uncertainty on the hardware sum.
void ProcessEvent() override
Process the event data according to pedestal and calibration factor.
void SmearByResolution(double resolution) override
QwVQWK_Channel(TString name, TString datatosave="raw")
static const Int_t kWordsPerChannel
Int_t fErrorCount_Sequence
sequence number check
void FillTreeVector(QwRootTreeBranchVector &values) const override
void AssignValueFrom(const VQwDataElement *valueptr) override
UInt_t fNumberOfSamples_map
Number of samples in the expected to read through the module. This value is set in the QwBeamline map...
QwVQWK_Channel & operator-=(const QwVQWK_Channel &value)
void Ratio(const QwVQWK_Channel &numer, const QwVQWK_Channel &denom)
void Difference(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
void AddValueFrom(const VQwHardwareChannel *valueptr) override
const QwVQWK_Channel operator+(const QwVQWK_Channel &value) const
void PrintValue() const override
Print single line of value and error of this data element.
Int_t ApplyHWChecks() override
Double_t fPrev_HardwareBlockSum
Previous Module-based sum of the four sub-blocks.
void PrintInfo() const override
Print multiple lines of information about this data element.
Double_t GetAverageVolts() const
void Blind(const QwBlinder *blinder)
Blind this channel as an asymmetry.
static const Double_t kTimePerSample
static const Double_t kVQWK_VoltsPerBit
void CopyFrom(const QwVQWK_Channel &value)
QwVQWK_Channel(const QwVQWK_Channel &value, VQwDataElement::EDataToSave datatosave)
void DivideBy(const VQwHardwareChannel *valueptr) override
void Product(const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
Int_t fErrorCount_SameHW
check to see ADC returning same HW value
Double_t fBlock[4]
Array of the sub-block data.
Double_t GetValueWidth() const
void SetEventData(Double_t *block, UInt_t sequencenumber=0)
static const Int_t kMaxChannels
Bool_t ApplySingleEventCuts() override
void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
const QwVQWK_Channel operator*(const QwVQWK_Channel &value) const
void AssignScaledValue(const QwVQWK_Channel &value, Double_t scale)
friend std::ostream & operator<<(std::ostream &stream, const QwVQWK_Channel &channel)
Double_t fBlockM2[4]
Second moment of the sub-block.
Double_t GetValue() const
void IncrementErrorCounters() override
Int_t fBlock_raw[4]
Array of the sub-block data as read from the module.
Double_t GetValueM2(size_t element) const override
void Scale(Double_t Offset) override
Double_t fHardwareBlockSum
Module-based sum of the four sub-blocks.
void EncodeEventData(std::vector< UInt_t > &buffer) override
Encode the event data into a CODA buffer.
void ArcTan(const QwVQWK_Channel &value)
Int_t fErrorCount_HWSat
check to see ADC channel is saturated
Double_t GetBlockValue(size_t blocknum) const
Int_t fErrorCount_SW_HW
HW_sum==SW_sum check.
Int_t GetRawValue() const
Int_t fSequenceNo_Counter
Internal counter to keep track of the sequence number.
Int_t GetRawHardwareSum() const
Double_t fHardwareBlockSumM2
Second moment of the hardware sum.
Int_t GetRawValue(size_t element) const override
size_t GetNumberOfSamples() const
static void PrintErrorCounterHead()
Int_t fNumEvtsWithEventCutsRejected
Counts the Event cut rejected events.
QwVQWK_Channel(const QwVQWK_Channel &value)
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)
void CopyFrom(const VQwHardwareChannel &value)
virtual VQwHardwareChannel * Clone() const
virtual void AddErrEntriesToList(std::vector< QwErrDBInterface > &)
Double_t GetValueWidth() const
void RangeCheck(size_t element) const
Checks that the requested element is in range, to be used in accesses to subelements similar to std::...
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