JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
VQwHardwareChannel.h
Go to the documentation of this file.
1/**********************************************************\
2* File: VQwHardwareChannel.h *
3* *
4* Author: P. King *
5* Date: Tue Mar 29 13:08:12 EDT 2011 *
6\**********************************************************/
7
8#pragma once
9
10// System headers
11#include <cmath>
12#include <vector>
13#include <stdexcept>
14
15// Qweak headers
16#include "VQwDataElement.h"
17
18// ROOT headers
19#ifdef HAS_RNTUPLE_SUPPORT
20#include <ROOT/RNTupleModel.hxx>
21#endif // HAS_RNTUPLE_SUPPORT
22
23// ROOT forward declarations
24class TTree;
25
26// Qweak forward declarations
27class QwDBInterface;
28class QwParameterFile;
31
32/**
33 * \class VQwHardwareChannel
34 * \ingroup QwAnalysis
35 * \brief Abstract base for concrete hardware channels implementing dual-operator pattern
36 *
37 * This class extends VQwDataElement to provide common services for hardware
38 * channel implementations that represent single physical readouts (ADC channels,
39 * scalers, etc.). It enforces the dual-operator architectural pattern at the
40 * channel level and provides infrastructure for calibration, cuts, and statistics.
41 *
42 * \par Dual-Operator Pattern Implementation:
43 * VQwHardwareChannel inherits the dual-operator requirement from VQwDataElement
44 * and adds channel-specific enforcement. Derived classes must implement:
45 *
46 * **Required Operator Pairs:**
47 * - `QwVQWK_Channel& operator+=(const QwVQWK_Channel&)` (type-specific)
48 * - `VQwHardwareChannel& operator+=(const VQwHardwareChannel&)` (polymorphic)
49 *
50 * **Polymorphic Delegation Pattern:**
51 * \code
52 * VQwHardwareChannel& QwVQWK_Channel::operator+=(const VQwHardwareChannel& source) {
53 * const QwVQWK_Channel* tmpptr = dynamic_cast<const QwVQWK_Channel*>(&source);
54 * if (tmpptr != NULL) {
55 * *this += *tmpptr; // Calls type-specific version
56 * } else {
57 * throw std::invalid_argument("Type mismatch in operator+=");
58 * }
59 * return *this;
60 * }
61 * \endcode
62 *
63 * \par Assignment + Operators Pattern:
64 * Sum/Difference methods follow the canonical pattern:
65 * \code
66 * void Sum(const QwVQWK_Channel& value1, const QwVQWK_Channel& value2) {
67 * *this = value1; // Uses derived class assignment operator
68 * *this += value2; // Uses type-specific operator+=
69 * }
70 * \endcode
71 *
72 * \par Channel Infrastructure:
73 * - **Calibration**: Pedestal subtraction and gain application
74 * - **Event Cuts**: Single-event limits with error flag propagation
75 * - **Statistics**: Running sums with error mask support
76 * - **Hardware Checks**: Burp detection and error counting
77 * - **Subelements**: Support for multi-element channels
78 *
79 * \par Representative Example:
80 * QwVQWK_Channel provides the canonical implementation demonstrating:
81 * - Complete dual-operator pattern with proper delegation
82 * - Six-word VQWK data processing (blocks 0-5)
83 * - Pedestal/calibration application
84 * - Single-event cuts and error propagation
85 * - Histogram and tree branch construction
86 *
87 * \par Error Handling Strategy:
88 * - Type mismatches in operators throw std::invalid_argument
89 * - Hardware errors set device-specific error codes
90 * - Event cuts use configurable upper/lower limits
91 * - Burp detection compares against reference channels
92 */
94/****************************************************************//**
95 * Class: VQwHardwareChannel
96 * Virtual base class to support common functions for all
97 * hardware channel data elements.
98 * Only the data element classes which contain raw data
99 * from one physical channel (such as QwVQWK_Channel,
100 * QwScaler_Channel, etc.) should inherit from this class.
101 ******************************************************************/
102public:
106 ~VQwHardwareChannel() override { };
107
108 void CopyFrom(const VQwHardwareChannel& value);
109
111
112 virtual VQwHardwareChannel* Clone() const{
113 return Clone(this->fDataToSave);
114 };
116
118
119 /*! \brief Get the number of data words in this data element */
121
122 /*! \brief Get the number of subelements in this data element */
124
125 Int_t GetRawValue() const {return this->GetRawValue(0);};
126 Double_t GetValue() const {return this->GetValue(0);};
127 Double_t GetValueM2() const {return this->GetValueM2(0);};
128 Double_t GetValueError() const {return this->GetValueError(0);};
129 Double_t GetValueWidth() const {return this->GetValueWidth(0);};
130 virtual Int_t GetRawValue(size_t element) const = 0;
131 virtual Double_t GetValue(size_t element) const = 0;
132 virtual Double_t GetValueM2(size_t element) const = 0;
133 virtual Double_t GetValueError(size_t element) const = 0;
134 Double_t GetValueWidth(size_t element) const {
135 RangeCheck(element);
136 Double_t width;
137 if (fGoodEventCount>0){
138 width = (GetValueError(element)*std::sqrt(Double_t(fGoodEventCount)));
139 } else {
140 width = 0.0;
141 }
142 return width;
143 };
144
145 void ClearEventData() override{
147 };
148
149 /* virtual void AddChannelOffset(Double_t Offset) = 0; */
150 virtual void Scale(Double_t Offset) = 0;
151
152
153 /// \brief Initialize the fields in this object
154 void InitializeChannel(TString name){InitializeChannel(name, "raw");};
155 virtual void InitializeChannel(TString name, TString datatosave) = 0;
156 virtual void InitializeChannel(TString subsystem, TString instrumenttype, TString name, TString datatosave) = 0;
157
158 //Check for hardware errors in the devices. This will return the device error code.
159 virtual Int_t ApplyHWChecks() = 0;
160
161 void SetEventCutMode(Int_t bcuts){bEVENTCUTMODE=bcuts;};
162
163 virtual Bool_t ApplySingleEventCuts() = 0;//check values read from modules are at desired level
164
165 virtual Bool_t CheckForBurpFail(const VQwHardwareChannel *event){
166 Bool_t foundburp = kFALSE;
167 if (fBurpThreshold>0){
168 Double_t diff = this->GetValue() - event->GetValue();
169 if (fabs(diff)>fBurpThreshold){
170 foundburp = kTRUE;
172 } else if (fBurpCountdown>0) {
173 foundburp = kTRUE;
175 }
176 }
177 if (foundburp){
179 }
180
181 return foundburp;
182 }
183
184 /*! \brief Set the upper and lower limits (fULimit and fLLimit)
185 * for this channel */
186 void SetSingleEventCuts(Double_t min, Double_t max);
187 /*! \brief Inherited from VQwDataElement to set the upper and lower
188 * limits (fULimit and fLLimit), stability % and the
189 * error flag on this channel */
190 void SetSingleEventCuts(UInt_t errorflag,Double_t min, Double_t max, Double_t stability=-1.0, Double_t BurpLevel=-1.0);
191
192 Double_t GetEventCutUpperLimit() const { return fULimit; };
193 Double_t GetEventCutLowerLimit() const { return fLLimit; };
194
195 Double_t GetStabilityLimit() const { return fStability;};
196
197 UInt_t UpdateErrorFlag() override {return GetEventcutErrorFlag();};
199 virtual UInt_t GetErrorCode() const {return (fErrorFlag);};
200
201 virtual void IncrementErrorCounters()=0;
202 virtual void ProcessEvent()=0;
203
204
205 virtual void CalculateRunningAverage() = 0;
206// virtual void AccumulateRunningSum(const VQwHardwareChannel *value) = 0;
207
208 /// Arithmetic assignment operator: Should only copy event-based data
210 if (this != &value) {
212 }
213 return *this;
214 }
215 void AssignScaledValue(const VQwHardwareChannel &value, Double_t scale){
216 AssignValueFrom(&value);
217 Scale(scale);
218 };
219 virtual void Ratio(const VQwHardwareChannel* numer, const VQwHardwareChannel* denom){
220 if (!IsNameEmpty()){
221 this->AssignValueFrom(numer);
222 this->operator/=(*denom);
223
224 // Remaining variables
226 fErrorFlag = (numer->fErrorFlag|denom->fErrorFlag);//error code is ORed.
227 }
228 }
229
230 void AssignValueFrom(const VQwDataElement* valueptr) override = 0;
235
236
237 virtual void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) = 0;
238
239 void SetPedestal(Double_t ped) { fPedestal = ped; kFoundPedestal = 1; };
240 Double_t GetPedestal() const { return fPedestal; };
241 void SetCalibrationFactor(Double_t factor) { fCalibrationFactor = factor; kFoundGain = 1; };
242 Double_t GetCalibrationFactor() const { return fCalibrationFactor; };
243
244 void AddEntriesToList(std::vector<QwDBInterface> &row_list);
245 virtual void AddErrEntriesToList(std::vector<QwErrDBInterface> & /*row_list*/) {};
246
247
248 virtual void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF){
249 if(count==0){
250 count = value->fGoodEventCount;
251 }
252 if(ErrorMask == kPreserveError){QwError << "VQwHardwareChannel count=" << count << QwLog::endl;}
253 AccumulateRunningSum(value, count, ErrorMask);
254 };
255 virtual void DeaccumulateRunningSum(const VQwHardwareChannel *value, Int_t ErrorMask=0xFFFFFFF){
256 AccumulateRunningSum(value, -1, ErrorMask);
257 };
258
259 virtual void AddValueFrom(const VQwHardwareChannel* valueptr) = 0;
260 virtual void SubtractValueFrom(const VQwHardwareChannel* valueptr) = 0;
261 virtual void MultiplyBy(const VQwHardwareChannel* valueptr) = 0;
262 virtual void DivideBy(const VQwHardwareChannel* valueptr) = 0;
263
264 virtual void ConstructBranchAndVector(TTree *tree, TString& prefix, QwRootTreeBranchVector& values) = 0;
265 virtual void ConstructBranch(TTree *tree, TString &prefix) = 0;
266 void ConstructBranch(TTree *tree, TString &prefix, QwParameterFile& modulelist);
267 virtual void FillTreeVector(QwRootTreeBranchVector& values) const = 0;
268#ifdef HAS_RNTUPLE_SUPPORT
269 virtual void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) = 0;
270 virtual void FillNTupleVector(std::vector<Double_t>& values) const = 0;
271#endif // HAS_RNTUPLE_SUPPORT
272
273 virtual void CopyParameters(const VQwHardwareChannel* /*valueptr*/){};
274
275 protected:
276 /*! \brief Set the number of data words in this data element */
277 void SetNumberOfDataWords(const UInt_t &numwords) {fNumberOfDataWords = numwords;}
278 /*! \brief Set the number of data words in this data element */
279 void SetNumberOfSubElements(const size_t elements) {fNumberOfSubElements = elements;};
280
281 /*! \brief Set the flag indicating if raw or derived values are
282 * in this data element */
283 void SetDataToSave(TString datatosave) {
284 if (datatosave == "raw")
286 else if (datatosave == "derived")
288 else
289 fDataToSave = kRaw; // wdc, added default fall-through
290 }
291 /*! \brief Set the flag indicating if raw or derived values are
292 * in this data element */
294 fDataToSave = datatosave;
295 }
296 /*! \brief Set the flag indicating if raw or derived values are
297 * in this data element based on prefix */
298 void SetDataToSaveByPrefix(const TString& prefix) {
299 if (prefix.Contains("asym_")
300 || prefix.Contains("diff_")
301 || prefix.Contains("yield_"))
303 if (prefix.Contains("stat"))
304 fDataToSave = kMoments; // stat has priority
305 }
306
307 /*! \brief Checks that the requested element is in range, to be
308 * used in accesses to subelements similar to
309 * std::vector::at(). */
310 void RangeCheck(size_t element) const {
311 if (element >= fNumberOfSubElements){
312 TString loc="VQwDataElement::RangeCheck for "
313 +this->GetElementName()+" failed for subelement "+Form("%zu",element);
314 throw std::out_of_range(loc.Data());
315
316 }
317 };
318
319protected:
320 UInt_t fNumberOfDataWords; ///< Number of raw data words in this data element
321 UInt_t fNumberOfSubElements; ///< Number of subelements in this data element
322
324
325 /* Ntuple array indices */
328
329 /*! \name Channel calibration */
330 // @{
331 Double_t fPedestal; /*!< Pedestal of the hardware sum signal,
332 we assume the pedestal level is constant over time
333 and can be divided by four for use with each block,
334 units: [counts / number of samples] */
338 //@}
339
340 /*! \name Single event cuts and errors */
341 // @{
342 Int_t bEVENTCUTMODE;/*!<If this set to kFALSE then Event cuts are OFF*/
343 Double_t fULimit, fLLimit;/*!<this sets the upper and lower limits*/
344 Double_t fStability;/*!<how much deviaton from the stable reading is allowed*/
345
348 //@}
349
350
351public:
352 /*! \name Global event cuts */
353 static void SetBurpHoldoff(Int_t holdoff) {
354 fBurpHoldoff = holdoff;
355 }
356
357protected:
358 // @{
359 static Int_t fBurpHoldoff;
360 // @}
361
362}; // class VQwHardwareChannel
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
Definition of the pure virtual base class of all data elements.
static const UInt_t kErrorFlag_BurpCut
Definition QwTypes.h:179
static const UInt_t kPreserveError
Definition QwTypes.h:186
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
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
UInt_t fGoodEventCount
Number of good events accumulated in this element.
VQwDataElement()
Default constructor.
virtual void ClearEventData()
Clear the event data in this element.
virtual UInt_t UpdateErrorFlag()
Update the error flag based on the error flags of internally contained objects Return parameter is th...
virtual UInt_t GetEventcutErrorFlag()
return the error flag on this channel/device
virtual const TString & GetElementName() const
Get the name of this element.
UInt_t fErrorFlag
This the standard error code generated for the channel that contains the global/local/stability flags...
Bool_t IsNameEmpty() const
Is the name of this element empty?
VQwDataElement & operator=(const VQwDataElement &value)
Arithmetic assignment operator: Should only copy event-based data.
Abstract base for concrete hardware channels implementing dual-operator pattern.
void AddEntriesToList(std::vector< QwDBInterface > &row_list)
virtual void CalculateRunningAverage()=0
void SetSingleEventCuts(Double_t min, Double_t max)
Set the upper and lower limits (fULimit and fLLimit) for this channel.
virtual VQwHardwareChannel & operator/=(const VQwHardwareChannel &input)=0
void SetDataToSave(TString datatosave)
Set the flag indicating if raw or derived values are in this data element.
Double_t GetValueError() const
virtual void MultiplyBy(const VQwHardwareChannel *valueptr)=0
virtual void ScaledAdd(Double_t scale, const VQwHardwareChannel *value)=0
void UpdateErrorFlag(const VQwHardwareChannel &elem)
size_t GetNumberOfSubelements()
Get the number of subelements in this data element.
virtual VQwHardwareChannel & operator+=(const VQwHardwareChannel &input)=0
Double_t GetCalibrationFactor() const
virtual void FillTreeVector(QwRootTreeBranchVector &values) const =0
Double_t GetValueM2() const
size_t GetNumberOfDataWords()
Get the number of data words in this data element.
static void SetBurpHoldoff(Int_t holdoff)
void SetCalibrationFactor(Double_t factor)
UInt_t fNumberOfSubElements
Number of subelements in this data element.
virtual Double_t GetValue(size_t element) const =0
virtual Double_t GetValueM2(size_t element) const =0
virtual UInt_t GetErrorCode() const
void SetDataToSave(VQwDataElement::EDataToSave datatosave)
Set the flag indicating if raw or derived values are in this data element.
Double_t GetEventCutUpperLimit() const
virtual void CopyParameters(const VQwHardwareChannel *)
void AssignScaledValue(const VQwHardwareChannel &value, Double_t scale)
virtual void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)=0
Double_t GetValueWidth(size_t element) const
virtual Int_t ApplyHWChecks()=0
void CopyFrom(const VQwHardwareChannel &value)
virtual void ConstructBranch(TTree *tree, TString &prefix)=0
UInt_t UpdateErrorFlag() override
Update the error flag based on the error flags of internally contained objects Return parameter is th...
virtual void InitializeChannel(TString subsystem, TString instrumenttype, TString name, TString datatosave)=0
virtual Bool_t CheckForBurpFail(const VQwHardwareChannel *event)
void SetNumberOfSubElements(const size_t elements)
Set the number of data words in this data element.
virtual void AddValueFrom(const VQwHardwareChannel *valueptr)=0
void SetPedestal(Double_t ped)
virtual Bool_t ApplySingleEventCuts()=0
virtual void SubtractValueFrom(const VQwHardwareChannel *valueptr)=0
virtual VQwHardwareChannel * Clone() const
VQwHardwareChannel & operator=(const VQwHardwareChannel &value)
Arithmetic assignment operator: Should only copy event-based data.
void SetNumberOfDataWords(const UInt_t &numwords)
Set the number of data words in this data element.
virtual void InitializeChannel(TString name, TString datatosave)=0
void InitializeChannel(TString name)
Initialize the fields in this object.
virtual VQwHardwareChannel & operator*=(const VQwHardwareChannel &input)=0
Double_t GetPedestal() const
virtual void Scale(Double_t Offset)=0
void ClearEventData() override
Clear the event data in this element.
virtual void AddErrEntriesToList(std::vector< QwErrDBInterface > &)
UInt_t fNumberOfDataWords
Number of raw data words in this data element.
void SetDataToSaveByPrefix(const TString &prefix)
Set the flag indicating if raw or derived values are in this data element based on prefix.
Double_t GetValueWidth() const
virtual void IncrementErrorCounters()=0
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 VQwHardwareChannel & operator-=(const VQwHardwareChannel &input)=0
virtual void DivideBy(const VQwHardwareChannel *valueptr)=0
virtual void ProcessEvent()=0
virtual Int_t GetRawValue(size_t element) const =0
virtual VQwHardwareChannel * Clone(VQwDataElement::EDataToSave datatosave) const =0
virtual void Ratio(const VQwHardwareChannel *numer, const VQwHardwareChannel *denom)
void AssignValueFrom(const VQwDataElement *valueptr) override=0
virtual Double_t GetValueError(size_t element) const =0
virtual void DeaccumulateRunningSum(const VQwHardwareChannel *value, Int_t ErrorMask=0xFFFFFFF)
Double_t GetStabilityLimit() const
void SetEventCutMode(Int_t bcuts)
Double_t GetEventCutLowerLimit() const
virtual void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)