JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwScaler_Channel.h
Go to the documentation of this file.
1/*!
2 * \file QwScaler_Channel.h
3 * \brief Base and derived classes for scaler channel data handling
4 *
5 * \author J. Pan
6 * \date Thu Sep 16 18:08:33 CDT 2009
7 */
8
9#pragma once
10
11// System headers
12#include <vector>
13
14// ROOT headers
15#include "TTree.h"
16
17// Qweak headers
18#include "VQwHardwareChannel.h"
19#include "MQwMockable.h"
20
21// Forward declarations
22class QwParameterFile;
23
24/**
25 * \class VQwScaler_Channel
26 * \ingroup QwAnalysis_ADC
27 * \brief Abstract base class for scaler channels
28 *
29 * Provides common infrastructure for scaler-based hardware channels,
30 * including differential scaler support, external clock normalization,
31 * and mock data generation. Concrete implementations handle the
32 * hardware-specific details of scaler decoding.
33 */
35
36public:
37 static Int_t GetBufferOffset(Int_t scalerindex, Int_t wordindex, UInt_t header = 1);
38 static void PrintErrorCounterHead();
39 static void PrintErrorCounterTail();
40
42
48
51
52public:
56
57 VQwScaler_Channel(TString name, TString datatosave = "raw"): MQwMockable() {
58 InitializeChannel(name,datatosave);
59 };
61 : VQwHardwareChannel(source),MQwMockable(source),
63 fValue_Raw(source.fValue_Raw),
64 fValue(source.fValue),
65 fValueM2(source.fValueM2),
67 // TODO: Don't copy the pointer; we need to regenerate it somehow.
68 //fNormChannelPtr(source.fNormChannelPtr);
73 { }
75 : VQwHardwareChannel(source,datatosave),MQwMockable(source),
77 fValue_Raw(source.fValue_Raw),
78 fValue(source.fValue),
79 fValueM2(source.fValueM2),
81 // TODO: Don't copy the pointer; we need to regenerate it somehow.
82 //fNormChannelPtr(source.fNormChannelPtr);
87 { }
88 ~VQwScaler_Channel() override { };
89
102
103
104 /// \brief Initialize the fields in this object
105 void InitializeChannel(TString name, TString datatosave = "raw") override;
106
107 /// \brief Initialize the fields in this object
108 void InitializeChannel(TString subsystem, TString instrumenttype, TString name, TString datatosave) override;
109
110 void SetDefaultSampleSize(size_t /*NumberOfSamples_map*/) {
111 //std::cerr << "QwScaler_Channel SetDefaultSampleSize does nothing!"
112 // << std::endl;
113 }
114
115 void LoadChannelParameters(QwParameterFile &paramfile) override;
116
117 void ClearEventData() override;
118
119 void RandomizeEventData(int helicity = 0, double time = 0.0) override;
120 void SetEventData(Double_t value);
121 void SetRawEventData() override{
122 //fValue = fCalibrationFactor * (Double_t(fValue_Raw) - Double_t(fValue_Raw_Old) - fPedestal);
123
127 else
128 fValue_Raw_Old = 0;
129
130};
131
132
133 void EncodeEventData(std::vector<UInt_t> &buffer) override = 0;
134 Int_t ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t index = 0) override = 0;
135
136//-----------------------------------------------------------------------------------------------
137 void SmearByResolution(double resolution) override;
138//-----------------------------------------------------------------------------------------------
139
140 void ProcessEvent() override;
141
142 Int_t GetRawValue(size_t /*element*/) const override { return fValue_Raw; };
143 Double_t GetValue(size_t /*element*/) const override { return fValue; };
144 Double_t GetValueM2(size_t /*element*/) const override { return fValueM2; };
145 Double_t GetValueError(size_t /*element*/) const override { return fValueError; };
146
148 void AssignScaledValue(const VQwScaler_Channel &value, Double_t scale);
149 void AssignValueFrom(const VQwDataElement* valueptr) override;
150 void AddValueFrom(const VQwHardwareChannel* valueptr) override;
151 void SubtractValueFrom(const VQwHardwareChannel* valueptr) override;
152 void MultiplyBy(const VQwHardwareChannel* valueptr) override;
153 void DivideBy(const VQwHardwareChannel* valueptr) override;
154
158
159 VQwHardwareChannel& operator+=(const VQwHardwareChannel& input) override;
160 VQwHardwareChannel& operator-=(const VQwHardwareChannel& input) override;
161 VQwHardwareChannel& operator*=(const VQwHardwareChannel& input) override;
162 VQwHardwareChannel& operator/=(const VQwHardwareChannel& input) override;
163
164 void Sum(VQwScaler_Channel &value1, VQwScaler_Channel &value2);
165 void Difference(VQwScaler_Channel &value1, VQwScaler_Channel &value2);
166 void Ratio(const VQwScaler_Channel &numer, const VQwScaler_Channel &denom);
167 void Product(VQwScaler_Channel &numer, VQwScaler_Channel &denom);
168 void AddChannelOffset(Double_t Offset);
169 void Scale(Double_t Offset) override;
170 void DivideBy(const VQwScaler_Channel &denom);
171
172
173 Int_t ApplyHWChecks() override; //Check for hardware errors in the devices. This will return the device error code.
174
175 Bool_t ApplySingleEventCuts() override;//check values read from modules are at desired level
176
177 Bool_t CheckForBurpFail(const VQwDataElement * /*ev_error*/) {return kFALSE;};
178
179 void IncrementErrorCounters() override;
180
181 /// report number of events failed due to HW and event cut failure
182 void PrintErrorCounters() const override;
183
184// UInt_t GetDeviceErrorCode(){//return the device error code
185// return fDeviceErrorCode;
186// };
187
188 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
189 void FillHistograms() override;
190
191 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override = 0;
192 void FillTreeVector(QwRootTreeBranchVector &values) const override = 0;
193 void ConstructBranch(TTree *tree, TString &prefix) override;
194#ifdef HAS_RNTUPLE_SUPPORT
195 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) override = 0;
196 void FillNTupleVector(std::vector<Double_t>& values) const override = 0;
197#endif // HAS_RNTUPLE_SUPPORT
198
199
200 void AccumulateRunningSum(const VQwScaler_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF);
201 void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override{
202 const VQwScaler_Channel *tmp_ptr = dynamic_cast<const VQwScaler_Channel*>(value);
203 if (tmp_ptr != NULL) {
204 AccumulateRunningSum(*tmp_ptr, count, ErrorMask);
205 } else {
206 throw std::invalid_argument("Standard exception from VQwScaler_Channel::AccumulateRunningSum: incompatible hardware channel type");
207 }
208 };
209 inline void DeaccumulateRunningSum(const VQwScaler_Channel& value, Int_t ErrorMask){
210 AccumulateRunningSum(value, -1, ErrorMask);
211 };
212
213 void PrintValue() const override;
214 void PrintInfo() const override;
215 void CalculateRunningAverage() override;
216
217 // These are related to those hardware channels that need to normalize
218 // to an external clock
219 Bool_t NeedsExternalClock() override { return fNeedsExternalClock; };
220 void SetNeedsExternalClock(Bool_t needed) override { fNeedsExternalClock = needed; };
221 std::string GetExternalClockName() override { return fNormChannelName; };
222 void SetExternalClockPtr( const VQwHardwareChannel* clock) override { fNormChannelPtr = clock; };
223 void SetExternalClockName( const std::string name) override { fNormChannelName = name; };
224
225 // Differential scalers automatically subtract the previous value
226 virtual Bool_t IsDifferentialScaler() { return fIsDifferentialScaler; };
227 virtual void SetDifferentialScaler(Bool_t diff) { fIsDifferentialScaler = diff; };
228
229 void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override;
230
231protected:
233
234protected:
235 static const Bool_t kDEBUG;
236
237 UInt_t fHeader;
240 Double_t fValue;
241 Double_t fValueM2;
242 Double_t fValueError;
245 std::string fNormChannelName;
246
249
250
251 Int_t fNumEvtsWithHWErrors;//counts the HW failed events
252 Int_t fNumEvtsWithEventCutsRejected;////counts the Event cut rejected events
253};
254
255
256// Derived templated class
257/**
258 * \class QwScaler_Channel
259 * \ingroup QwAnalysis_ADC
260 * \brief Templated concrete scaler channel with configurable data masking
261 *
262 * Template specialization of VQwScaler_Channel that applies data_mask
263 * and data_shift to the raw scaler values during processing. Commonly
264 * used to handle different scaler module types with varying data formats.
265 */
266template <UInt_t data_mask=0xffffffff, UInt_t data_shift=0 >
268{
269 public:
270
271 // Define the constructors (cascade)
274 : VQwScaler_Channel(source) { };
275 QwScaler_Channel(TString name, TString datatosave = "raw")
276 : VQwScaler_Channel(name,datatosave) { };
278 : VQwScaler_Channel(source,datatosave) { };
279
283
284 public:
285
286 // Implement the templated methods
287 void EncodeEventData(std::vector<UInt_t> &buffer) override;
288 Int_t ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t index = 0) override;
289
290 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override;
291 void FillTreeVector(QwRootTreeBranchVector &values) const override;
292#ifdef HAS_RNTUPLE_SUPPORT
293 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) override;
294 void FillNTupleVector(std::vector<Double_t>& values) const override;
295#endif // HAS_RNTUPLE_SUPPORT
296
297
298};
299
300// These typedef's should be the last things in the file.
301// Class template instationation must be made in the source
302// file for anything defined here.
303typedef class QwScaler_Channel<0x00ffffff,0> QwSIS3801D24_Channel;
304typedef class QwScaler_Channel<0xffffffff,0> QwSIS3801D32_Channel;
305typedef class QwScaler_Channel<0xffffffff,0> QwSIS3801_Channel;
306typedef class QwScaler_Channel<0xffffffff,0> QwSTR7200_Channel;
class QwScaler_Channel< 0x00ffffff, 0 > QwSIS3801D24_Channel
class QwScaler_Channel< 0xffffffff, 0 > QwSIS3801_Channel
class QwScaler_Channel< 0xffffffff, 0 > QwSTR7200_Channel
class QwScaler_Channel< 0xffffffff, 0 > QwSIS3801D32_Channel
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
Abstract base class for scaler channels.
void AccumulateRunningSum(const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
static const Bool_t kDEBUG
void SetDefaultSampleSize(size_t)
void Difference(VQwScaler_Channel &value1, VQwScaler_Channel &value2)
void InitializeChannel(TString name, TString datatosave="raw") override
Initialize the fields in this object.
void Product(VQwScaler_Channel &numer, VQwScaler_Channel &denom)
std::string GetExternalClockName() override
void ProcessEvent() override
static void PrintErrorCounterTail()
void ClearEventData() override
Clear the event data in this element.
Double_t GetValue(size_t) const override
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t index=0) override=0
Process the CODA event buffer for this element.
void Scale(Double_t Offset) override
Bool_t NeedsExternalClock() override
void AddValueFrom(const VQwHardwareChannel *valueptr) override
void LoadChannelParameters(QwParameterFile &paramfile) override
void SubtractValueFrom(const VQwHardwareChannel *valueptr) override
void FillTreeVector(QwRootTreeBranchVector &values) const override=0
VQwScaler_Channel(const VQwScaler_Channel &source, VQwDataElement::EDataToSave datatosave)
void ConstructHistograms(TDirectory *folder, TString &prefix) override
Construct the histograms for this data element.
Double_t GetValueError(size_t) const override
Bool_t ApplySingleEventCuts() override
void RandomizeEventData(int helicity=0, double time=0.0) override
Internally generate random event data.
void ConstructBranch(TTree *tree, TString &prefix) override
void SetNeedsExternalClock(Bool_t needed) override
VQwScaler_Channel & operator=(const VQwScaler_Channel &value)
VQwScaler_Channel(const VQwScaler_Channel &source)
void Sum(VQwScaler_Channel &value1, VQwScaler_Channel &value2)
void SetExternalClockName(const std::string name) override
void FillHistograms() override
Fill the histograms for this data element.
VQwHardwareChannel & operator/=(const VQwHardwareChannel &input) override
void PrintInfo() const override
Print multiple lines of information about this data element.
virtual void SetDifferentialScaler(Bool_t diff)
void DeaccumulateRunningSum(const VQwScaler_Channel &value, Int_t ErrorMask)
static Int_t GetBufferOffset(Int_t scalerindex, Int_t wordindex, UInt_t header=1)
void PrintErrorCounters() const override
report number of events failed due to HW and event cut failure
void ScaledAdd(Double_t scale, const VQwHardwareChannel *value) override
std::string fNormChannelName
void SmearByResolution(double resolution) override
void SetRawEventData() override
void MultiplyBy(const VQwHardwareChannel *valueptr) override
Bool_t CheckForBurpFail(const VQwDataElement *)
Int_t ApplyHWChecks() override
Int_t GetRawValue(size_t) const override
VQwScaler_Channel & operator+=(const VQwScaler_Channel &value)
void SetExternalClockPtr(const VQwHardwareChannel *clock) override
virtual Bool_t IsDifferentialScaler()
const VQwHardwareChannel * fNormChannelPtr
static void PrintErrorCounterHead()
void PrintValue() const override
Print single line of value and error of this data element.
void SetEventData(Double_t value)
VQwScaler_Channel & operator*=(const VQwScaler_Channel &value)
void CalculateRunningAverage() override
void Ratio(const VQwScaler_Channel &numer, const VQwScaler_Channel &denom)
void DivideBy(const VQwHardwareChannel *valueptr) override
VQwScaler_Channel & operator-=(const VQwScaler_Channel &value)
Double_t GetValueM2(size_t) const override
~VQwScaler_Channel() override
void AddChannelOffset(Double_t Offset)
void IncrementErrorCounters() override
void AssignValueFrom(const VQwDataElement *valueptr) override
void AssignScaledValue(const VQwScaler_Channel &value, Double_t scale)
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override=0
void AccumulateRunningSum(const VQwScaler_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
void CopyFrom(const VQwScaler_Channel &value)
VQwScaler_Channel(TString name, TString datatosave="raw")
void EncodeEventData(std::vector< UInt_t > &buffer) override=0
Encode the event data into a CODA buffer.
Templated concrete scaler channel with configurable data masking.
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
QwScaler_Channel(const QwScaler_Channel &source, VQwDataElement::EDataToSave datatosave)
QwScaler_Channel(const QwScaler_Channel &source)
QwScaler_Channel(TString name, TString datatosave="raw")
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.
VQwHardwareChannel * Clone(VQwDataElement::EDataToSave datatosave) const override
void FillTreeVector(QwRootTreeBranchVector &values) const override
The pure virtual base class of all data elements.
VQwDataElement()
Default constructor.
Abstract base for concrete hardware channels implementing dual-operator pattern.
Double_t GetValueError() const
Double_t GetValueM2() const
void CopyFrom(const VQwHardwareChannel &value)
virtual VQwHardwareChannel * Clone() const
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)