JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches

Concrete hardware channel for VQWK ADC modules (6x32-bit words) More...

#include <QwVQWK_Channel.h>

+ Inheritance diagram for QwVQWK_Channel:
+ Collaboration diagram for QwVQWK_Channel:

Public Member Functions

 QwVQWK_Channel ()
 
 QwVQWK_Channel (TString name, TString datatosave="raw")
 
 QwVQWK_Channel (const QwVQWK_Channel &value)
 
 QwVQWK_Channel (const QwVQWK_Channel &value, VQwDataElement::EDataToSave datatosave)
 
 ~QwVQWK_Channel () override
 
void CopyFrom (const QwVQWK_Channel &value)
 
VQwHardwareChannelClone (VQwDataElement::EDataToSave datatosave) const override
 
void InitializeChannel (TString name, TString datatosave) override
 Initialize the fields in this object.
 
void InitializeChannel (TString subsystem, TString instrumenttype, TString name, TString datatosave) override
 Initialize the fields in this object.
 
void LoadChannelParameters (QwParameterFile &paramfile) override
 
void SetDefaultSampleSize (size_t num_samples_map)
 
void ClearEventData () override
 Clear the event data in this element.
 
void RandomizeEventData (int helicity=0.0, double time=0.0) override
 Internally generate random event data.
 
void ForceMapfileSampleSize ()
 
void SmearByResolution (double resolution) override
 
void SetHardwareSum (Double_t hwsum, UInt_t sequencenumber=0)
 
void SetEventData (Double_t *block, UInt_t sequencenumber=0)
 
void SetRawEventData () override
 
void EncodeEventData (std::vector< UInt_t > &buffer) override
 Encode the event data into a CODA buffer.
 
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 ProcessEvent () override
 Process the event data according to pedestal and calibration factor.
 
QwVQWK_Channeloperator= (const QwVQWK_Channel &value)
 
void AssignScaledValue (const QwVQWK_Channel &value, Double_t scale)
 
void AssignValueFrom (const VQwDataElement *valueptr) override
 
void AddValueFrom (const VQwHardwareChannel *valueptr) override
 
void SubtractValueFrom (const VQwHardwareChannel *valueptr) override
 
void MultiplyBy (const VQwHardwareChannel *valueptr) override
 
void DivideBy (const VQwHardwareChannel *valueptr) override
 
void ArcTan (const QwVQWK_Channel &value)
 
QwVQWK_Channeloperator+= (const QwVQWK_Channel &value)
 
QwVQWK_Channeloperator-= (const QwVQWK_Channel &value)
 
QwVQWK_Channeloperator*= (const QwVQWK_Channel &value)
 
VQwHardwareChanneloperator+= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator-= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator*= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator/= (const VQwHardwareChannel &input) override
 
const QwVQWK_Channel operator+ (const QwVQWK_Channel &value) const
 
const QwVQWK_Channel operator- (const QwVQWK_Channel &value) const
 
const QwVQWK_Channel operator* (const QwVQWK_Channel &value) const
 
void Sum (const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
 
void Difference (const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
 
void Ratio (const QwVQWK_Channel &numer, const QwVQWK_Channel &denom)
 
void Product (const QwVQWK_Channel &value1, const QwVQWK_Channel &value2)
 
void DivideBy (const QwVQWK_Channel &denom)
 
void AddChannelOffset (Double_t Offset)
 
void Scale (Double_t Offset) override
 
void AccumulateRunningSum (const QwVQWK_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
 
void AccumulateRunningSum (const VQwHardwareChannel *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
 
void DeaccumulateRunningSum (const QwVQWK_Channel &value, Int_t ErrorMask=0xFFFFFFF)
 
void CalculateRunningAverage () override
 
Bool_t MatchSequenceNumber (size_t seqnum)
 
Bool_t MatchNumberOfSamples (size_t numsamp)
 
Bool_t ApplySingleEventCuts (Double_t LL, Double_t UL)
 
Bool_t ApplySingleEventCuts () override
 
void PrintErrorCounters () const override
 report number of events failed due to HW and event cut failure
 
void SetVQWKSaturationLimt (Double_t sat_volts=8.5)
 
Double_t GetVQWKSaturationLimt ()
 
Int_t ApplyHWChecks () override
 
void IncrementErrorCounters () override
 
void ConstructHistograms (TDirectory *folder, TString &prefix) override
 Construct the histograms for this data element.
 
void FillHistograms () override
 Fill the histograms for this data element.
 
void ConstructBranchAndVector (TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
 
void ConstructBranch (TTree *tree, TString &prefix) override
 
void FillTreeVector (QwRootTreeBranchVector &values) const override
 
Int_t GetRawValue (size_t element) const override
 
Double_t GetValue (size_t element) const override
 
Double_t GetValueM2 (size_t element) const override
 
Double_t GetValueError (size_t element) const override
 
Double_t GetAverageVolts () const
 
size_t GetSequenceNumber () const
 
size_t GetNumberOfSamples () const
 
void SetCalibrationToVolts ()
 
void PrintValue () const override
 Print single line of value and error of this data element.
 
void PrintInfo () const override
 Print multiple lines of information about this data element.
 
void Blind (const QwBlinder *blinder)
 Blind this channel as an asymmetry.
 
void Blind (const QwBlinder *blinder, const QwVQWK_Channel &yield)
 Blind this channel as a difference.
 
void ScaledAdd (Double_t scale, const VQwHardwareChannel *value) override
 
virtual void LoadMockDataParameters (QwParameterFile &paramfile)
 Load the mock data parameters from the current line in the param file.
 
Int_t GetRawValue () const
 
Double_t GetValue () const
 
Double_t GetValueM2 () const
 
Double_t GetValueError () const
 
Double_t GetValueWidth () const
 
Double_t GetValueWidth (size_t element) const
 
virtual void DeaccumulateRunningSum (const VQwHardwareChannel *value, Int_t ErrorMask=0xFFFFFFF)
 
virtual VQwHardwareChannelClone () const
 
- Public Member Functions inherited from VQwHardwareChannel
 VQwHardwareChannel ()
 
 VQwHardwareChannel (const VQwHardwareChannel &value)
 
 VQwHardwareChannel (const VQwHardwareChannel &value, VQwDataElement::EDataToSave datatosave)
 
 ~VQwHardwareChannel () override
 
void CopyFrom (const VQwHardwareChannel &value)
 
void ProcessOptions ()
 
size_t GetNumberOfDataWords ()
 Get the number of data words in this data element.
 
size_t GetNumberOfSubelements ()
 Get the number of subelements in this data element.
 
Int_t GetRawValue () const
 
Double_t GetValue () const
 
Double_t GetValueM2 () const
 
Double_t GetValueError () const
 
Double_t GetValueWidth () const
 
Double_t GetValueWidth (size_t element) const
 
void ClearEventData () override
 Clear the event data in this element.
 
void InitializeChannel (TString name)
 Initialize the fields in this object.
 
void SetEventCutMode (Int_t bcuts)
 
Bool_t CheckForBurpFail (const VQwHardwareChannel *event)
 
void SetSingleEventCuts (Double_t min, Double_t max)
 Set the upper and lower limits (fULimit and fLLimit) for this channel.
 
void SetSingleEventCuts (UInt_t errorflag, Double_t min, Double_t max, Double_t stability=-1.0, Double_t BurpLevel=-1.0)
 Inherited from VQwDataElement to set the upper and lower limits (fULimit and fLLimit), stability % and the error flag on this channel.
 
Double_t GetEventCutUpperLimit () const
 
Double_t GetEventCutLowerLimit () const
 
Double_t GetStabilityLimit () const
 
UInt_t UpdateErrorFlag () override
 Update the error flag based on the error flags of internally contained objects Return parameter is the "Eventcut Error Flag".
 
void UpdateErrorFlag (const VQwHardwareChannel &elem)
 
virtual UInt_t GetErrorCode () const
 
VQwHardwareChanneloperator= (const VQwHardwareChannel &value)
 Arithmetic assignment operator: Should only copy event-based data.
 
void AssignScaledValue (const VQwHardwareChannel &value, Double_t scale)
 
void Ratio (const VQwHardwareChannel *numer, const VQwHardwareChannel *denom)
 
void SetPedestal (Double_t ped)
 
Double_t GetPedestal () const
 
void SetCalibrationFactor (Double_t factor)
 
Double_t GetCalibrationFactor () const
 
void AddEntriesToList (std::vector< QwDBInterface > &row_list)
 
virtual void AddErrEntriesToList (std::vector< QwErrDBInterface > &)
 
void ConstructBranch (TTree *tree, TString &prefix, QwParameterFile &modulelist)
 
virtual void CopyParameters (const VQwHardwareChannel *)
 
void UpdateErrorFlag (const UInt_t &error)
 
- Public Member Functions inherited from VQwDataElement
 VQwDataElement ()
 Default constructor.
 
 VQwDataElement (const VQwDataElement &value)
 Copy constructor.
 
 ~VQwDataElement () override
 Virtual destructor.
 
void CopyFrom (const VQwDataElement &value)
 
Bool_t IsNameEmpty () const
 Is the name of this element empty?
 
void SetElementName (const TString &name)
 Set the name of this element.
 
virtual const TString & GetElementName () const
 Get the name of this element.
 
size_t GetNumberOfDataWords ()
 Get the number of data words in this data element.
 
UInt_t GetGoodEventCount () const
 
VQwDataElementoperator+= (const VQwDataElement &)
 Addition-assignment operator.
 
VQwDataElementoperator-= (const VQwDataElement &)
 Subtraction-assignment operator.
 
void Sum (const VQwDataElement &, const VQwDataElement &)
 Sum operator (base class fallback throws runtime error)
 
void Difference (const VQwDataElement &, const VQwDataElement &)
 Difference operator (base class fallback throws runtime error)
 
void Ratio (const VQwDataElement &, const VQwDataElement &)
 Ratio operator (base class fallback throws runtime error)
 
void SetSingleEventCuts (UInt_t, Double_t, Double_t, Double_t)
 set the upper and lower limits (fULimit and fLLimit), stability % and the error flag on this channel
 
Bool_t CheckForBurpFail (const VQwDataElement *)
 
virtual UInt_t GetEventcutErrorFlag ()
 return the error flag on this channel/device
 
virtual void SetNeedsExternalClock (Bool_t)
 
virtual Bool_t NeedsExternalClock ()
 
virtual std::string GetExternalClockName ()
 
virtual void SetExternalClockPtr (const VQwHardwareChannel *)
 
virtual void SetExternalClockName (const std::string)
 
virtual Double_t GetNormClockValue ()
 
TString GetSubsystemName () const
 Return the name of the inheriting subsystem name.
 
void SetSubsystemName (TString sysname)
 Set the name of the inheriting subsystem name.
 
TString GetModuleType () const
 Return the type of the beam instrument.
 
void SetModuleType (TString ModuleType)
 set the type of the beam instrument
 
- Public Member Functions inherited from MQwHistograms
void ShareHistograms (const MQwHistograms *source)
 Share histogram pointers between objects.
 
- Public Member Functions inherited from MQwMockable
 MQwMockable ()
 
virtual ~MQwMockable ()
 
void SetRandomEventDriftParameters (Double_t amplitude, Double_t phase, Double_t frequency)
 Set a single set of harmonic drift parameters.
 
void AddRandomEventDriftParameters (Double_t amplitude, Double_t phase, Double_t frequency)
 Add drift parameters to the internal set.
 
void SetRandomEventParameters (Double_t mean, Double_t sigma)
 Set the normal random event parameters.
 
void SetRandomEventAsymmetry (Double_t asymmetry)
 Set the helicity asymmetry.
 
Double_t GetRandomValue ()
 
void UseExternalRandomVariable ()
 Set the flag to use an externally provided random variable.
 
void SetExternalRandomVariable (Double_t random_variable)
 Set the externally provided random variable.
 
void SetMockDataAsDiff ()
 

Static Public Member Functions

static Int_t GetBufferOffset (Int_t moduleindex, Int_t channelindex)
 
static void PrintErrorCounterHead ()
 
static void PrintErrorCounterTail ()
 
static void SetBurpHoldoff (Int_t holdoff)
 

Static Public Attributes

static const Double_t kTimePerSample = (2.0/30.0) * Qw::us
 

Protected Member Functions

QwVQWK_Channeloperator/= (const QwVQWK_Channel &value)
 
- Protected Member Functions inherited from VQwHardwareChannel
void SetNumberOfDataWords (const UInt_t &numwords)
 Set the number of data words in this data element.
 
void SetNumberOfSubElements (const size_t elements)
 Set the number of data words in this data element.
 
void SetDataToSave (TString datatosave)
 Set the flag indicating if raw or derived values are in this data element.
 
void SetDataToSave (VQwDataElement::EDataToSave datatosave)
 Set the flag indicating if raw or derived values are 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.
 
void RangeCheck (size_t element) const
 Checks that the requested element is in range, to be used in accesses to subelements similar to std::vector::at().
 
- Protected Member Functions inherited from VQwDataElement
void SetNumberOfDataWords (const UInt_t &numwords)
 Set the number of data words in this data element.
 
VQwDataElementoperator= (const VQwDataElement &value)
 Arithmetic assignment operator: Should only copy event-based data.
 
void UpdateErrorFlag (const UInt_t &error)
 
- Protected Member Functions inherited from MQwHistograms
 MQwHistograms ()
 Default constructor.
 
 MQwHistograms (const MQwHistograms &source)
 Copy constructor.
 
virtual ~MQwHistograms ()
 Virtual destructor.
 
MQwHistogramsoperator= (const MQwHistograms &value)
 
void Fill_Pointer (TH1_ptr hist_ptr, Double_t value)
 
void AddHistogram (TH1 *h)
 Register a histogram.
 

Private Member Functions

Double_t GetBlockValue (size_t blocknum) const
 
Double_t GetBlockErrorValue (size_t blocknum) const
 
Double_t GetHardwareSum () const
 
Double_t GetHardwareSumM2 () const
 
Double_t GetHardwareSumWidth () const
 
Double_t GetHardwareSumError () const
 
Int_t GetRawBlockValue (size_t blocknum) const
 
Int_t GetRawHardwareSum () const
 
Int_t GetRawSoftwareSum () const
 

Private Attributes

Channel configuration data members
Short_t fBlocksPerEvent
 
Event data members—Raw values
Int_t fBlock_raw [4]
 Array of the sub-block data as read from the module.
 
Int_t fHardwareBlockSum_raw
 Module-based sum of the four sub-blocks as read from the module.
 
Int_t fSoftwareBlockSum_raw
 Sum of the data in the four sub-blocks raw.
 
Event data members—Potentially calibrated values
Double_t fBlock [4]
 Array of the sub-block data.
 
Double_t fHardwareBlockSum
 Module-based sum of the four sub-blocks.
 

Static Private Attributes

static const Bool_t kDEBUG = kFALSE
 
static const Int_t kWordsPerChannel = 6
 
static const Int_t kMaxChannels = 8
 
ADC Calibration
static const Double_t kVQWK_VoltsPerBit = (20./(1<<18))
 

Friends

std::ostream & operator<< (std::ostream &stream, const QwVQWK_Channel &channel)
 

Calculation of the statistical moments

static const Bool_t bDEBUG =kFALSE
 debugging display purposes
 
Double_t fBlockM2 [4]
 Second moment of the sub-block.
 
Double_t fBlockError [4]
 Uncertainty on the sub-block.
 
Double_t fHardwareBlockSumM2
 Second moment of the hardware sum.
 
Double_t fHardwareBlockSumError
 Uncertainty on the hardware sum.
 
UInt_t fSequenceNumber
 Event sequence number for this channel.
 
UInt_t fPreviousSequenceNumber
 Previous event sequence number for this channel.
 
UInt_t fNumberOfSamples
 Number of samples read through the module.
 
UInt_t fNumberOfSamples_map
 Number of samples in the expected to read through the module. This value is set in the QwBeamline map file.
 
Int_t fErrorCount_HWSat
 check to see ADC channel is saturated
 
Int_t fErrorCount_sample
 for sample size check
 
Int_t fErrorCount_SW_HW
 HW_sum==SW_sum check.
 
Int_t fErrorCount_Sequence
 sequence number check
 
Int_t fErrorCount_SameHW
 check to see ADC returning same HW value
 
Int_t fErrorCount_ZeroHW
 check to see ADC returning zero
 
Int_t fNumEvtsWithEventCutsRejected
 Counts the Event cut rejected events.
 
Int_t fADC_Same_NumEvt
 Keep track of how many events with same ADC value returned.
 
Int_t fSequenceNo_Prev
 Keep the sequence number of the last event.
 
Int_t fSequenceNo_Counter
 Internal counter to keep track of the sequence number.
 
Double_t fPrev_HardwareBlockSum
 Previous Module-based sum of the four sub-blocks.
 
Double_t fSaturationABSLimit
 absolute value of the VQWK saturation volt
 
Bool_t bHw_sum
 
Bool_t bHw_sum_raw
 
Bool_t bBlock
 
Bool_t bBlock_raw
 
Bool_t bNum_samples
 
Bool_t bDevice_Error_Code
 
Bool_t bSequence_number
 

Additional Inherited Members

- Public Types inherited from VQwDataElement
enum  EDataToSave { kRaw = 0 , kDerived , kMoments }
 
- Protected Attributes inherited from VQwHardwareChannel
UInt_t fNumberOfDataWords
 Number of raw data words in this data element.
 
UInt_t fNumberOfSubElements
 Number of subelements in this data element.
 
EDataToSave fDataToSave
 
size_t fTreeArrayIndex
 
size_t fTreeArrayNumEntries
 
Double_t fPedestal
 
Double_t fCalibrationFactor
 
Bool_t kFoundPedestal
 
Bool_t kFoundGain
 
Int_t bEVENTCUTMODE
 
Double_t fULimit
 
Double_t fLLimit
 
Double_t fStability
 
Double_t fBurpThreshold
 
Int_t fBurpCountdown
 
- Protected Attributes inherited from VQwDataElement
TString fElementName
 Name of this data element.
 
UInt_t fNumberOfDataWords
 Number of raw data words in this data element.
 
UInt_t fGoodEventCount
 Number of good events accumulated in this element.
 
TString fSubsystemName
 
TString fModuleType
 
UInt_t fErrorFlag
 This the standard error code generated for the channel that contains the global/local/stability flags and the Device error code (Unique error code for HW failures)
 
UInt_t fErrorConfigFlag
 contains the global/local/stability flags
 
- Protected Attributes inherited from MQwHistograms
std::vector< TH1_ptrfHistograms
 Histograms associated with this data element.
 
bool fUseExternalRandomVariable
 Flag to use an externally provided normal random variable.
 
double fExternalRandomVariable
 Externally provided normal random variable.
 
bool fCalcMockDataAsDiff
 
Double_t fMockAsymmetry
 Helicity asymmetry.
 
Double_t fMockGaussianMean
 Mean of normal distribution.
 
Double_t fMockGaussianSigma
 Sigma of normal distribution.
 
std::vector< Double_t > fMockDriftAmplitude
 Harmonic drift amplitude.
 
std::vector< Double_t > fMockDriftFrequency
 Harmonic drift frequency.
 
std::vector< Double_t > fMockDriftPhase
 Harmonic drift phase.
 
static Int_t fBurpHoldoff = 10
 
static std::mt19937 fRandomnessGenerator
 Internal randomness generator.
 
static std::normal_distribution< double > fNormalDistribution
 Internal normal probability distribution.
 
static std::function< double()> fNormalRandomVariable = []() -> double { return MQwMockable::fNormalDistribution(MQwMockable::fRandomnessGenerator); }
 Internal normal random variable.
 

Detailed Description

Concrete hardware channel for VQWK ADC modules (6x32-bit words)

Decodes and processes the data for a single VQWK channel, providing access to sub-blocks and hardware sums, statistical moments, single-event cuts, and running statistics. Implements the dual-operator pattern by offering efficient type-specific operators and delegating polymorphic operators from the VQwHardwareChannel interface using dynamic_cast.

Definition at line 48 of file QwVQWK_Channel.h.

Constructor & Destructor Documentation

◆ QwVQWK_Channel() [1/4]

QwVQWK_Channel::QwVQWK_Channel ( )
inline

Definition at line 76 of file QwVQWK_Channel.h.

76 : MQwMockable() {
77 InitializeChannel("","");
78 SetVQWKSaturationLimt(8.5);//set the default saturation limit
79 };
void InitializeChannel(TString name, TString datatosave) override
Initialize the fields in this object.
void SetVQWKSaturationLimt(Double_t sat_volts=8.5)

References InitializeChannel(), MQwMockable::MQwMockable(), and SetVQWKSaturationLimt().

Referenced by AccumulateRunningSum(), AccumulateRunningSum(), AddValueFrom(), ArcTan(), AssignScaledValue(), AssignValueFrom(), Blind(), Clone(), CopyFrom(), DeaccumulateRunningSum(), Difference(), DivideBy(), DivideBy(), MultiplyBy(), operator*(), operator*=(), operator*=(), operator+(), operator+=(), operator+=(), operator-(), operator-=(), operator-=(), operator/=(), operator/=(), operator<<, operator=(), Product(), QwVQWK_Channel(), QwVQWK_Channel(), Ratio(), ScaledAdd(), SubtractValueFrom(), and Sum().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ QwVQWK_Channel() [2/4]

QwVQWK_Channel::QwVQWK_Channel ( TString name,
TString datatosave = "raw" )
inline

Definition at line 80 of file QwVQWK_Channel.h.

80 : MQwMockable() {
81 InitializeChannel(name, datatosave);
82 SetVQWKSaturationLimt(8.5);//set the default saturation limit
83 };

References InitializeChannel(), MQwMockable::MQwMockable(), and SetVQWKSaturationLimt().

+ Here is the call graph for this function:

◆ QwVQWK_Channel() [3/4]

QwVQWK_Channel::QwVQWK_Channel ( const QwVQWK_Channel & value)
inline

Definition at line 84 of file QwVQWK_Channel.h.

84 :
85 VQwHardwareChannel(value), MQwMockable(value),
89 {
90 *this = value;
91 };
Short_t fBlocksPerEvent
Double_t fSaturationABSLimit
absolute value of the VQWK saturation volt
UInt_t fNumberOfSamples_map
Number of samples in the expected to read through the module. This value is set in the QwBeamline map...

References fBlocksPerEvent, fNumberOfSamples_map, fSaturationABSLimit, MQwMockable::MQwMockable(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ QwVQWK_Channel() [4/4]

QwVQWK_Channel::QwVQWK_Channel ( const QwVQWK_Channel & value,
VQwDataElement::EDataToSave datatosave )
inline

Definition at line 92 of file QwVQWK_Channel.h.

92 :
93 VQwHardwareChannel(value,datatosave), MQwMockable(value),
97 {
98 *this = value;
99 };

References fBlocksPerEvent, fNumberOfSamples_map, fSaturationABSLimit, MQwMockable::MQwMockable(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ ~QwVQWK_Channel()

QwVQWK_Channel::~QwVQWK_Channel ( )
inlineoverride

Definition at line 100 of file QwVQWK_Channel.h.

100{ };

Member Function Documentation

◆ AccumulateRunningSum() [1/2]

void QwVQWK_Channel::AccumulateRunningSum ( const QwVQWK_Channel & value,
Int_t count = 0,
Int_t ErrorMask = 0xFFFFFFF )

Accumulate event values into the running sum with optional scaling.

Parameters
valueSource channel to accumulate from.
countEvent count scaling (0 means use value.fGoodEventCount).
ErrorMaskBit mask of error flags to exclude when accumulating.

Moments and uncertainty calculation on the running sums and averages The calculation of the first and second moments of the running sum is not completely straightforward due to numerical instabilities associated with small variances and large average values. The naive computation taking the difference of the square of the average and the average of the squares leads to the subtraction of two very large numbers to get a small number.

Alternative algorithms (including for higher order moments) are supplied in Pebay, Philippe (2008), "Formulas for Robust, One-Pass Parallel Computation of Covariances and Arbitrary-Order Statistical Moments", Technical Report SAND2008-6212, Sandia National Laboratories. http://infoserve.sandia.gov/sand_doc/2008/086212.pdf

In the following formulas the moments \( M^1 \) and \( M^2 \) are defined

\begin{eqnarray} M^1 & = & \frac{1}{n} \sum^n y \\ M^2 & = & \sum^n (y - \mu) \end{eqnarray}

Recurrence relations for the addition of a single event:

\begin{eqnarray} M^1_n & = & M^1_{n-1} + \frac{y - M^1_{n-1}}{n} \\ M^2_n & = & M^2_{n-1} + (y - M^1_{n-1})(y - M^1_n) \end{eqnarray}

For the addition of an already accumulated sum:

\begin{eqnarray} M^1 & = & M^1_1 + n_2 \frac{M^1_2 - M^1_1}{n} \\ M^2 & = & M^2_1 + M^2_2 + n_1 n_2 \frac{(M^1_2 - M^1_1)^2}{n} \end{eqnarray}

In these recursive formulas we start from \( M^1 = y \) and \( M^2 = 0 \).

To calculate the mean and standard deviation we use

\begin{eqnarray} \mu & = & M^1 \\ \sigma^2 & = & \frac{1}{n} M^2 \end{eqnarray}

The standard deviation is a biased estimator, but this is what ROOT uses. Better would be to divide by \( (n-1) \).

We use the formulas provided there for the calculation of the first and second moments (i.e. average and variance).

Definition at line 1570 of file QwVQWK_Channel.cc.

1571{
1572 /*
1573 note:
1574 The AccumulateRunningSum is called on a dedicated subsystem array object and
1575 for the standard running avg computations we only need value.fErrorFlag==0
1576 events to be included in the running avg. So the "berror" conditions is only
1577 used for the stability check purposes.
1578
1579 The need for this check below came due to fact that when routine
1580 DeaccumulateRunningSum is called the errorflag is updated with
1581 the kBeamStabilityError flag (+ configuration flags for global errors) and
1582 need to make sure we remove this flag and any configuration flags before
1583 checking the (fErrorFlag != 0) condition
1584
1585 See how the stability check is implemented in the QwEventRing class
1586
1587 Rakitha
1588 */
1589
1590 if(count==0){
1591 count = value.fGoodEventCount;
1592 }
1593
1594 Int_t n1 = fGoodEventCount;
1595 Int_t n2 = count;
1596
1597 // If there are no good events, check the error flag
1598 if (n2 == 0 && (value.fErrorFlag == 0)) {
1599 n2 = 1;
1600 }
1601
1602 // If a single event is removed from the sum, check all but stability fail flags
1603 if (n2 == -1) {
1604 if ((value.fErrorFlag & ErrorMask) == 0) {
1605 n2 = -1;
1606 } else {
1607 n2 = 0;
1608 }
1609 }
1610
1611 if (ErrorMask == kPreserveError){
1612 //n = 1;
1613 if (n2 == 0) {
1614 n2 = 1;
1615 }
1616 if (count == -1) {
1617 n2 = -1;
1618 }
1619 }
1620
1621 // New total number of good events
1622 Int_t n = n1 + n2;
1623
1624 // Set up variables
1625 Double_t M11 = fHardwareBlockSum;
1626 Double_t M12 = value.fHardwareBlockSum;
1627 Double_t M22 = value.fHardwareBlockSumM2;
1628
1629 //if(this->GetElementName() == "bcm_an_ds3" && ErrorMask == kPreserveError){QwError << "count=" << fGoodEventCount << " n=" << n << QwLog::endl; }
1630 if (n2 == 0) {
1631 // no good events for addition
1632 return;
1633 } else if (n2 == -1) {
1634 // simple version for removal of single event from the sum
1636 if (n > 1) {
1637 fHardwareBlockSum -= (M12 - M11) / n;
1638 fHardwareBlockSumM2 -= (M12 - M11)
1639 * (M12 - fHardwareBlockSum); // note: using updated mean
1640 // and for individual blocks
1641 for (Int_t i = 0; i < 4; i++) {
1642 M11 = fBlock[i];
1643 M12 = value.fBlock[i];
1644 M22 = value.fBlockM2[i];
1645 fBlock[i] -= (M12 - M11) / n;
1646 fBlockM2[i] -= (M12 - M11) * (M12 - fBlock[i]); // note: using updated mean
1647 }
1648 } else if (n == 1) {
1649 fHardwareBlockSum -= (M12 - M11) / n;
1650 fHardwareBlockSumM2 -= (M12 - M11)
1651 * (M12 - fHardwareBlockSum); // note: using updated mean
1652 if (fabs(fHardwareBlockSumM2) < 10.*std::numeric_limits<double>::epsilon())
1653 fHardwareBlockSumM2 = 0; // rounding
1654 // and for individual blocks
1655 for (Int_t i = 0; i < 4; i++) {
1656 M11 = fBlock[i];
1657 M12 = value.fBlock[i];
1658 M22 = value.fBlockM2[i];
1659 fBlock[i] -= (M12 - M11) / n;
1660 fBlockM2[i] -= (M12 - M11) * (M12 - fBlock[i]); // note: using updated mean
1661 if (fabs(fBlockM2[i]) < 10.*std::numeric_limits<double>::epsilon())
1662 fBlockM2[i] = 0; // rounding
1663 }
1664 } else if (n == 0) {
1665 fHardwareBlockSum -= M12;
1666 fHardwareBlockSumM2 -= M22;
1667 if (fabs(fHardwareBlockSum) < 10.*std::numeric_limits<double>::epsilon())
1668 fHardwareBlockSum = 0; // rounding
1669 if (fabs(fHardwareBlockSumM2) < 10.*std::numeric_limits<double>::epsilon())
1670 fHardwareBlockSumM2 = 0; // rounding
1671 // and for individual blocks
1672 for (Int_t i = 0; i < 4; i++) {
1673 M11 = fBlock[i];
1674 M12 = value.fBlock[i];
1675 M22 = value.fBlockM2[i];
1676 fBlock[i] -= M12;
1677 fBlockM2[i] -= M22;
1678 if (fabs(fBlock[i]) < 10.*std::numeric_limits<double>::epsilon())
1679 fBlock[i] = 0; // rounding
1680 if (fabs(fBlockM2[i]) < 10.*std::numeric_limits<double>::epsilon())
1681 fBlockM2[i] = 0; // rounding
1682 }
1683 } else {
1684 QwWarning << "Running sum has deaccumulated to negative good events." << QwLog::endl;
1685 }
1686 } else if (n2 == 1) {
1687 // simple version for addition of single event
1689 fHardwareBlockSum += (M12 - M11) / n;
1690 fHardwareBlockSumM2 += (M12 - M11)
1691 * (M12 - fHardwareBlockSum); // note: using updated mean
1692 // and for individual blocks
1693 for (Int_t i = 0; i < 4; i++) {
1694 M11 = fBlock[i];
1695 M12 = value.fBlock[i];
1696 M22 = value.fBlockM2[i];
1697 fBlock[i] += (M12 - M11) / n;
1698 fBlockM2[i] += (M12 - M11) * (M12 - fBlock[i]); // note: using updated mean
1699 }
1700 } else if (n2 > 1) {
1701 // general version for addition of multi-event sets
1702 fGoodEventCount += n2;
1703 fHardwareBlockSum += n2 * (M12 - M11) / n;
1704 fHardwareBlockSumM2 += M22 + n1 * n2 * (M12 - M11) * (M12 - M11) / n;
1705 // and for individual blocks
1706 for (Int_t i = 0; i < 4; i++) {
1707 M11 = fBlock[i];
1708 M12 = value.fBlock[i];
1709 M22 = value.fBlockM2[i];
1710 fBlock[i] += n2 * (M12 - M11) / n;
1711 fBlockM2[i] += M22 + n1 * n2 * (M12 - M11) * (M12 - M11) / n;
1712 }
1713 }
1714
1715 // Nanny
1717 QwWarning << "Angry Nanny: NaN detected in " << GetElementName() << QwLog::endl;
1718}
#define QwWarning
Predefined log drain for warnings.
Definition QwLog.h:44
static const UInt_t kPreserveError
Definition QwTypes.h:186
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
Double_t fBlock[4]
Array of the sub-block data.
Double_t fBlockM2[4]
Second moment of the sub-block.
Double_t fHardwareBlockSum
Module-based sum of the four sub-blocks.
Double_t fHardwareBlockSumM2
Second moment of the hardware sum.
UInt_t fGoodEventCount
Number of good events accumulated in this element.
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...

References QwLog::endl(), fBlock, fBlockM2, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fHardwareBlockSum, fHardwareBlockSumM2, VQwDataElement::GetElementName(), kPreserveError, QwVQWK_Channel(), and QwWarning.

Referenced by AccumulateRunningSum(), and DeaccumulateRunningSum().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ AccumulateRunningSum() [2/2]

void QwVQWK_Channel::AccumulateRunningSum ( const VQwHardwareChannel * value,
Int_t count = 0,
Int_t ErrorMask = 0xFFFFFFF )
inlineoverridevirtual

Reimplemented from VQwHardwareChannel.

Definition at line 194 of file QwVQWK_Channel.h.

194 {
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 };
void AccumulateRunningSum(const QwVQWK_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)

References AccumulateRunningSum(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ AddChannelOffset()

void QwVQWK_Channel::AddChannelOffset ( Double_t offset)

This function will add a offset to the hw_sum and add the same offset for blocks.

Definition at line 1491 of file QwVQWK_Channel.cc.

1492{
1493 if (!IsNameEmpty()){
1494 fHardwareBlockSum += offset;
1495 for (Int_t i=0; i<fBlocksPerEvent; i++)
1496 fBlock[i] += offset;
1497 }
1498 return;
1499}
Bool_t IsNameEmpty() const
Is the name of this element empty?

References fBlock, fBlocksPerEvent, fHardwareBlockSum, and VQwDataElement::IsNameEmpty().

+ Here is the call graph for this function:

◆ AddValueFrom()

void QwVQWK_Channel::AddValueFrom ( const VQwHardwareChannel * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1185 of file QwVQWK_Channel.cc.

1186{
1187 const QwVQWK_Channel* tmpptr;
1188 tmpptr = dynamic_cast<const QwVQWK_Channel*>(valueptr);
1189 if (tmpptr!=NULL){
1190 *this += *tmpptr;
1191 } else {
1192 TString loc="Standard exception from QwVQWK_Channel::AddValueFrom = "
1193 +valueptr->GetElementName()+" is an incompatible type.";
1194 throw std::invalid_argument(loc.Data());
1195 }
1196}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ ApplyHWChecks()

Int_t QwVQWK_Channel::ApplyHWChecks ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 65 of file QwVQWK_Channel.cc.

66{
67 Bool_t fEventIsGood=kTRUE;
68 Bool_t bStatus;
69 if (bEVENTCUTMODE>0){//Global switch to ON/OFF event cuts set at the event cut file
70
71 if (bDEBUG)
72 QwWarning<<" QwQWVK_Channel "<<GetElementName()<<" "<<GetNumberOfSamples()<<QwLog::endl;
73
74 // Sample size check
75 bStatus = MatchNumberOfSamples(fNumberOfSamples_map);//compare the default sample size with no.of samples read by the module
76 if (!bStatus) {
78 }
79
80 // Check SW and HW return the same sum
81 bStatus = (GetRawHardwareSum() == GetRawSoftwareSum());
82 //fEventIsGood = bStatus;
83 if (!bStatus) {
85 }
86
87
88
89 //check sequence number
91 if (fSequenceNo_Counter==0 || GetSequenceNumber()==0){//starting the data run
93 }
94
95 if (!MatchSequenceNumber(fSequenceNo_Prev)){//we have a sequence number error
96 fEventIsGood=kFALSE;
98 if (bDEBUG) QwWarning<<" QwQWVK_Channel "<<GetElementName()<<" Sequence number previous value = "<<fSequenceNo_Prev<<" Current value= "<< GetSequenceNumber()<<QwLog::endl;
99 }
100
102
103 //Checking for HW_sum is returning same value.
105 //std::cout<<" BCM hardware sum is different "<<std::endl;
108 }else
109 fADC_Same_NumEvt++;//hw_sum is same increment the counter
110
111 //check for the hw_sum is giving the same value
112 if (fADC_Same_NumEvt>0){//we have ADC stuck with same value
113 if (bDEBUG) QwWarning<<" BCM hardware sum is same for more than "<<fADC_Same_NumEvt<<" time consecutively "<<QwLog::endl;
115 }
116
117 //check for the hw_sum is zero
118 if (GetRawHardwareSum()==0){
120 }
121 if (!fEventIsGood)
122 fSequenceNo_Counter=0;//resetting the counter after ApplyHWChecks() a failure
123
125 if (bDEBUG)
126 QwWarning << this->GetElementName()<<" "<<GetRawHardwareSum() << "Saturating VQWK invoked! " <<TMath::Abs(GetRawHardwareSum())*kVQWK_VoltsPerBit/fNumberOfSamples<<" Limit "<<GetVQWKSaturationLimt() << QwLog::endl;
128 }
129
130 }
131 else {
132 fGoodEventCount = 1;
133 fErrorFlag = 0;
134 }
135
136 return fErrorFlag;
137}
static const UInt_t kErrorFlag_ZeroHW
Definition QwTypes.h:164
static const UInt_t kErrorFlag_SW_HW
Definition QwTypes.h:161
static const UInt_t kErrorFlag_sample
Definition QwTypes.h:160
static const UInt_t kErrorFlag_VQWK_Sat
Definition QwTypes.h:159
static const UInt_t kErrorFlag_SameHW
Definition QwTypes.h:163
static const UInt_t kErrorFlag_Sequence
Definition QwTypes.h:162
Int_t fADC_Same_NumEvt
Keep track of how many events with same ADC value returned.
Bool_t MatchNumberOfSamples(size_t numsamp)
Double_t GetVQWKSaturationLimt()
UInt_t fNumberOfSamples
Number of samples read through the module.
Int_t GetRawSoftwareSum() const
Int_t fSequenceNo_Prev
Keep the sequence number of the last event.
Bool_t MatchSequenceNumber(size_t seqnum)
static const Bool_t bDEBUG
debugging display purposes
size_t GetSequenceNumber() const
Double_t fPrev_HardwareBlockSum
Previous Module-based sum of the four sub-blocks.
static const Double_t kVQWK_VoltsPerBit
Int_t fSequenceNo_Counter
Internal counter to keep track of the sequence number.
Int_t GetRawHardwareSum() const
size_t GetNumberOfSamples() const

References bDEBUG, VQwHardwareChannel::bEVENTCUTMODE, QwLog::endl(), fADC_Same_NumEvt, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fNumberOfSamples, fNumberOfSamples_map, fPrev_HardwareBlockSum, fSequenceNo_Counter, fSequenceNo_Prev, VQwDataElement::GetElementName(), GetNumberOfSamples(), GetRawHardwareSum(), GetRawSoftwareSum(), GetSequenceNumber(), GetVQWKSaturationLimt(), kErrorFlag_SameHW, kErrorFlag_sample, kErrorFlag_Sequence, kErrorFlag_SW_HW, kErrorFlag_VQWK_Sat, kErrorFlag_ZeroHW, kVQWK_VoltsPerBit, MatchNumberOfSamples(), MatchSequenceNumber(), and QwWarning.

+ Here is the call graph for this function:

◆ ApplySingleEventCuts() [1/2]

Bool_t QwVQWK_Channel::ApplySingleEventCuts ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1876 of file QwVQWK_Channel.cc.

1877{
1878 Bool_t status;
1879
1880 if (bEVENTCUTMODE>=2){//Global switch to ON/OFF event cuts set at the event cut file
1881
1882 if (fULimit < fLLimit){
1883 status=kTRUE;
1884 } else if (GetHardwareSum()<=fULimit && GetHardwareSum()>=fLLimit){
1885 if ((fErrorFlag)==0)
1886 status=kTRUE;
1887 else
1888 status=kFALSE;//If the device HW is failed
1889 }
1890 else{
1891 if (GetHardwareSum()> fULimit)
1893 else
1895 status=kFALSE;
1896 }
1897
1898 if (bEVENTCUTMODE==3){
1899 status=kTRUE; //Update the event cut fail flag but pass the event.
1900 }
1901
1902
1903 }
1904 else{
1905 status=kTRUE;
1906 //fErrorFlag=0;//we need to keep the device error codes
1907 }
1908
1909 return status;
1910}
static const UInt_t kErrorFlag_EventCut_L
Definition QwTypes.h:165
static const UInt_t kErrorFlag_EventCut_U
Definition QwTypes.h:166
Double_t GetHardwareSum() const

References VQwHardwareChannel::bEVENTCUTMODE, VQwDataElement::fErrorFlag, VQwHardwareChannel::fLLimit, VQwHardwareChannel::fULimit, GetHardwareSum(), kErrorFlag_EventCut_L, and kErrorFlag_EventCut_U.

+ Here is the call graph for this function:

◆ ApplySingleEventCuts() [2/2]

Bool_t QwVQWK_Channel::ApplySingleEventCuts ( Double_t LL,
Double_t UL )

Definition at line 1860 of file QwVQWK_Channel.cc.

1861{
1862 Bool_t status = kFALSE;
1863
1864 if (UL < LL){
1865 status=kTRUE;
1866 } else if (GetHardwareSum()<=UL && GetHardwareSum()>=LL){
1867 if ((fErrorFlag & kPreserveError)!=0)
1868 status=kTRUE;
1869 else
1870 status=kFALSE;//If the device HW is failed
1871 }
1872 std::cout<<(this->fErrorFlag & kPreserveError)<<std::endl;
1873 return status;
1874}

References VQwDataElement::fErrorFlag, GetHardwareSum(), and kPreserveError.

+ Here is the call graph for this function:

◆ ArcTan()

void QwVQWK_Channel::ArcTan ( const QwVQWK_Channel & value)

Definition at line 1453 of file QwVQWK_Channel.cc.

1454{
1455 if (!IsNameEmpty()) {
1456 this->fHardwareBlockSum = 0.0;
1457 for (Int_t i=0; i<fBlocksPerEvent; i++) {
1458 this->fBlock[i] = atan(value.fBlock[i]);
1459 this->fHardwareBlockSum += this->fBlock[i];
1460 }
1462 }
1463
1464 return;
1465}

References fBlock, fBlocksPerEvent, fHardwareBlockSum, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ AssignScaledValue()

void QwVQWK_Channel::AssignScaledValue ( const QwVQWK_Channel & value,
Double_t scale )

Definition at line 1153 of file QwVQWK_Channel.cc.

1155{
1156 if(this == &value) return;
1157
1158 if (!IsNameEmpty()) {
1159 for (Int_t i=0; i<fBlocksPerEvent; i++){
1160 this->fBlock[i] = value.fBlock[i] * scale;
1161 this->fBlockM2[i] = value.fBlockM2[i] * scale * scale;
1162 }
1163 this->fHardwareBlockSum = value.fHardwareBlockSum * scale;
1164 this->fHardwareBlockSumM2 = value.fHardwareBlockSumM2 * scale * scale;
1165 this->fHardwareBlockSumError = value.fHardwareBlockSumError; // Keep this?
1166 this->fGoodEventCount = value.fGoodEventCount;
1167 this->fNumberOfSamples = value.fNumberOfSamples;
1168 this->fSequenceNumber = value.fSequenceNumber;
1169 this->fErrorFlag = value.fErrorFlag;
1170 }
1171}
UInt_t fSequenceNumber
Event sequence number for this channel.
Double_t fHardwareBlockSumError
Uncertainty on the hardware sum.

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fHardwareBlockSum, fHardwareBlockSumError, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ AssignValueFrom()

void QwVQWK_Channel::AssignValueFrom ( const VQwDataElement * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1173 of file QwVQWK_Channel.cc.

1174{
1175 const QwVQWK_Channel* tmpptr;
1176 tmpptr = dynamic_cast<const QwVQWK_Channel*>(valueptr);
1177 if (tmpptr!=NULL){
1178 *this = *tmpptr;
1179 } else {
1180 TString loc="Standard exception from QwVQWK_Channel::AssignValueFrom = "
1181 +valueptr->GetElementName()+" is an incompatible type.";
1182 throw std::invalid_argument(loc.Data());
1183 }
1184}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwDataElement::VQwDataElement().

+ Here is the call graph for this function:

◆ Blind() [1/2]

void QwVQWK_Channel::Blind ( const QwBlinder * blinder)

Blind this channel as an asymmetry.

Blind this channel as an asymmetry

Parameters
blinderBlinder

Definition at line 1793 of file QwVQWK_Channel.cc.

1794{
1795 if (!IsNameEmpty()) {
1796 if (blinder->IsBlinderOkay() && ((fErrorFlag)==0) ){
1797 for (Int_t i = 0; i < fBlocksPerEvent; i++)
1798 blinder->BlindValue(fBlock[i]);
1799 blinder->BlindValue(fHardwareBlockSum);
1800 } else {
1802 for (Int_t i = 0; i < fBlocksPerEvent; i++)
1805 }
1806 }
1807 return;
1808}
const Bool_t & IsBlinderOkay() const
Definition QwBlinder.h:206
void ModifyThisErrorCode(UInt_t &errorcode) const
Definition QwBlinder.h:119
void BlindValue(Double_t &value) const
Asymmetry blinding.
Definition QwBlinder.h:124
static constexpr const Double_t kValue_BlinderFail
Definition QwBlinder.h:79

References QwBlinder::BlindValue(), fBlock, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, QwBlinder::IsBlinderOkay(), VQwDataElement::IsNameEmpty(), QwBlinder::kValue_BlinderFail, and QwBlinder::ModifyThisErrorCode().

+ Here is the call graph for this function:

◆ Blind() [2/2]

void QwVQWK_Channel::Blind ( const QwBlinder * blinder,
const QwVQWK_Channel & yield )

Blind this channel as a difference.

Blind this channel as a difference with specified yield

Parameters
blinderBlinder
yieldCorresponding yield

Definition at line 1815 of file QwVQWK_Channel.cc.

1816{
1817 if (!IsNameEmpty()) {
1818 if (blinder->IsBlinderOkay() && ((fErrorFlag) ==0) ){
1819 for (Int_t i = 0; i < fBlocksPerEvent; i++)
1820 blinder->BlindValue(fBlock[i], yield.fBlock[i]);
1822 } else {
1823 blinder->ModifyThisErrorCode(fErrorFlag);//update the HW error code
1824 for (Int_t i = 0; i < fBlocksPerEvent; i++)
1827 }
1828 }
1829 return;
1830}

References QwBlinder::BlindValue(), fBlock, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, QwBlinder::IsBlinderOkay(), VQwDataElement::IsNameEmpty(), QwBlinder::kValue_BlinderFail, QwBlinder::ModifyThisErrorCode(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ CalculateRunningAverage()

void QwVQWK_Channel::CalculateRunningAverage ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1721 of file QwVQWK_Channel.cc.

1722{
1723 if (fGoodEventCount <= 0)
1724 {
1725 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
1726 fBlockError[i] = 0.0;
1727 }
1729 }
1730 else
1731 {
1732 // We use a biased estimator by dividing by n. Use (n - 1) to get the
1733 // unbiased estimator for the standard deviation.
1734 //
1735 // Note we want to calculate the error here, not sigma:
1736 // sigma = sqrt(M2 / n);
1737 // error = sigma / sqrt (n) = sqrt(M2) / n;
1738 for (Int_t i = 0; i < fBlocksPerEvent; i++)
1739 fBlockError[i] = sqrt(fBlockM2[i]) / fGoodEventCount;
1741
1742 // Stability check 83951872
1744 // check to see the channel has stability cut activated in the event cut file
1745 if (GetValueWidth() > fStability){
1746 // if the width is greater than the stability required flag the event
1748 } else
1749 fErrorFlag = 0;
1750 }
1751 }
1752}
static const UInt_t kBeamStabilityError
Definition QwTypes.h:180
static const UInt_t kStabilityCut
Definition QwTypes.h:184
Double_t fBlockError[4]
Uncertainty on the sub-block.
Double_t GetValueWidth() const
UInt_t fErrorConfigFlag
contains the global/local/stability flags

References fBlockError, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorConfigFlag, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fHardwareBlockSumError, fHardwareBlockSumM2, VQwHardwareChannel::fStability, GetValueWidth(), kBeamStabilityError, and kStabilityCut.

+ Here is the call graph for this function:

◆ ClearEventData()

void QwVQWK_Channel::ClearEventData ( )
overridevirtual

Clear the event data in this element.

Reimplemented from VQwDataElement.

Definition at line 250 of file QwVQWK_Channel.cc.

251{
252 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
253 fBlock_raw[i] = 0;
254 fBlock[i] = 0.0;
255 fBlockM2[i] = 0.0;
256 fBlockError[i] = 0.0;
257 }
260 fHardwareBlockSum = 0.0;
263 fSequenceNumber = 0;
265 fGoodEventCount = 0;
266 fErrorFlag=0;
267 return;
268}
Int_t fHardwareBlockSum_raw
Module-based sum of the four sub-blocks as read from the module.
Int_t fSoftwareBlockSum_raw
Sum of the data in the four sub-blocks raw.
Int_t fBlock_raw[4]
Array of the sub-block data as read from the module.

References fBlock, fBlock_raw, fBlockError, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fHardwareBlockSum, fHardwareBlockSum_raw, fHardwareBlockSumError, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, and fSoftwareBlockSum_raw.

Referenced by InitializeChannel(), QwLinearDiodeArray::ProcessEvent(), and QwQPD::ProcessEvent().

+ Here is the caller graph for this function:

◆ Clone() [1/2]

virtual VQwHardwareChannel * VQwHardwareChannel::Clone ( ) const
inlinevirtual

Reimplemented from VQwHardwareChannel.

Definition at line 112 of file VQwHardwareChannel.h.

112 {
113 return Clone(this->fDataToSave);
114 };
virtual VQwHardwareChannel * Clone() const

◆ Clone() [2/2]

VQwHardwareChannel * QwVQWK_Channel::Clone ( VQwDataElement::EDataToSave datatosave) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 112 of file QwVQWK_Channel.h.

112 {
113 return new QwVQWK_Channel(*this,datatosave);
114 };

References QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ ConstructBranch()

void QwVQWK_Channel::ConstructBranch ( TTree * tree,
TString & prefix )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 827 of file QwVQWK_Channel.cc.

828{
829 // This channel is not used, so skip setting up the tree.
830 if (IsNameEmpty()) return;
831
832 TString basename = prefix + GetElementName();
833 tree->Branch(basename,&fHardwareBlockSum,basename+"/D");
834 if (kDEBUG){
835 std::cerr << "QwVQWK_Channel::ConstructBranchAndVector: fTreeArrayIndex==" << fTreeArrayIndex
836 << "; fTreeArrayNumEntries==" << fTreeArrayNumEntries
837 << std::endl;
838 }
839}
static const Bool_t kDEBUG

References fHardwareBlockSum, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), and kDEBUG.

+ Here is the call graph for this function:

◆ ConstructBranchAndVector()

void QwVQWK_Channel::ConstructBranchAndVector ( TTree * tree,
TString & prefix,
QwRootTreeBranchVector & values )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 743 of file QwVQWK_Channel.cc.

744{
745 // This channel is not used, so skip setting up the tree.
746 if (IsNameEmpty()) return;
747
748 // Decide what to store based on prefix
749 SetDataToSaveByPrefix(prefix);
750
751 TString basename = prefix(0, (prefix.First("|") >= 0)? prefix.First("|"): prefix.Length()) + GetElementName();
752 fTreeArrayIndex = values.size();
753
759 bDevice_Error_Code = gQwHists.MatchVQWKElementFromList(GetSubsystemName().Data(), GetModuleType().Data(), "Device_Error_Code");
761
762 if (bHw_sum) {
763 values.push_back("hw_sum", 'I');
764 if (fDataToSave == kMoments) {
765 values.push_back("hw_sum_m2", 'D');
766 values.push_back("hw_sum_err", 'D');
767 }
768 }
769
770 if (bBlock) {
771 values.push_back("block0", 'D');
772 values.push_back("block1", 'D');
773 values.push_back("block2", 'D');
774 values.push_back("block3", 'D');
775 }
776
777 if (bNum_samples) {
778 values.push_back("num_samples", 'I');
779 }
780
781 if (bDevice_Error_Code) {
782 values.push_back("Device_Error_Code", 'i');
783 }
784
785 if (fDataToSave == kRaw) {
786 if (bHw_sum_raw) {
787 values.push_back("hw_sum_raw", 'I');
788 }
789 if (bBlock_raw) {
790 values.push_back("block0_raw", 'I');
791 values.push_back("block1_raw", 'I');
792 values.push_back("block2_raw", 'I');
793 values.push_back("block3_raw", 'I');
794 }
795 if (bSequence_number) {
796 values.push_back("sequence_number", 'i');
797 }
798 }
799
801
802 std::string leaf_list = values.LeafList();
803
804 if (gQwHists.MatchDeviceParamsFromList(basename.Data())
807
808 // This is for the RT mode
809 if (leaf_list == "hw_sum/D")
810 leaf_list = basename+"/D";
811
812 if (kDEBUG)
813 QwMessage << "base name " << basename << " List " << leaf_list << QwLog::endl;
814
815 tree->Branch(basename, &(values[fTreeArrayIndex]), leaf_list.c_str());
816 }
817
818 if (kDEBUG) {
819 std::cerr << "QwVQWK_Channel::ConstructBranchAndVector: fTreeArrayIndex==" << fTreeArrayIndex
820 << "; fTreeArrayNumEntries==" << fTreeArrayNumEntries
821 << "; values.size()==" << values.size()
822 << "; list==" << leaf_list
823 << std::endl;
824 }
825}
QwHistogramHelper gQwHists
Globally defined instance of the QwHistogramHelper class.
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
Bool_t MatchVQWKElementFromList(const std::string &subsystemname, const std::string &moduletype, const std::string &devicename)
Bool_t MatchDeviceParamsFromList(const std::string &devicename)
size_type size() const noexcept
Definition QwRootFile.h:83
std::string LeafList(size_type start_index=0) const
Definition QwRootFile.h:230
void push_back(const std::string &name, const char type='D')
Definition QwRootFile.h:197
Bool_t bDevice_Error_Code
TString GetSubsystemName() const
Return the name of the inheriting subsystem name.
TString GetModuleType() const
Return the type of the beam instrument.
void SetDataToSaveByPrefix(const TString &prefix)
Set the flag indicating if raw or derived values are in this data element based on prefix.

References bBlock, bBlock_raw, bDevice_Error_Code, bHw_sum, bHw_sum_raw, bNum_samples, bSequence_number, QwLog::endl(), VQwHardwareChannel::fDataToSave, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, VQwDataElement::GetElementName(), VQwDataElement::GetModuleType(), VQwDataElement::GetSubsystemName(), gQwHists, VQwDataElement::IsNameEmpty(), kDEBUG, VQwDataElement::kMoments, VQwDataElement::kRaw, QwRootTreeBranchVector::LeafList(), QwRootTreeBranchVector::push_back(), QwMessage, VQwHardwareChannel::SetDataToSaveByPrefix(), and QwRootTreeBranchVector::size().

+ Here is the call graph for this function:

◆ ConstructHistograms()

void QwVQWK_Channel::ConstructHistograms ( TDirectory * folder,
TString & prefix )
overridevirtual

Construct the histograms for this data element.

Implements VQwDataElement.

Definition at line 637 of file QwVQWK_Channel.cc.

638{
639 // If we have defined a subdirectory in the ROOT file, then change into it.
640 if (folder != NULL) folder->cd();
641
642 if (IsNameEmpty()){
643 // This channel is not used, so skip filling the histograms.
644 } else {
645 // Now create the histograms.
646 SetDataToSaveByPrefix(prefix);
647
648 TString basename = prefix + GetElementName();
649
650 if(fDataToSave==kRaw)
651 {
652 fHistograms.resize(8+2+1, NULL);
653 size_t index=0;
654 for (Int_t i=0; i<fBlocksPerEvent; i++){
655 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_block%d_raw",i));
656 fHistograms[index+1] = gQwHists.Construct1DHist(basename+Form("_block%d",i));
657 index += 2;
658 }
659 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_hw_raw"));
660 fHistograms[index+1] = gQwHists.Construct1DHist(basename+Form("_hw"));
661 index += 2;
662 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_sw-hw_raw"));
663 }
664 else if(fDataToSave==kDerived)
665 {
666 fHistograms.resize(4+1+1, NULL);
667 Int_t index=0;
668 for (Int_t i=0; i<fBlocksPerEvent; i++){
669 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_block%d",i));
670 index += 1;
671 }
672 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_hw"));
673 index += 1;
674 fHistograms[index] = gQwHists.Construct1DHist(basename+Form("_dev_err"));
675 index += 1;
676 }
677 else
678 {
679 // this is not recognized
680 }
681
682 }
683}
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
TH1F * Construct1DHist(const TString &inputfile, const TString &name_title)

References fBlocksPerEvent, VQwHardwareChannel::fDataToSave, MQwHistograms::fHistograms, VQwDataElement::GetElementName(), gQwHists, VQwDataElement::IsNameEmpty(), VQwDataElement::kDerived, VQwDataElement::kRaw, and VQwHardwareChannel::SetDataToSaveByPrefix().

+ Here is the call graph for this function:

◆ CopyFrom()

void QwVQWK_Channel::CopyFrom ( const QwVQWK_Channel & value)
inline

Definition at line 102 of file QwVQWK_Channel.h.

102 {
107 *this = value;
108 };
void CopyFrom(const VQwHardwareChannel &value)

References VQwHardwareChannel::CopyFrom(), fBlocksPerEvent, fNumberOfSamples_map, fSaturationABSLimit, and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ DeaccumulateRunningSum() [1/2]

void QwVQWK_Channel::DeaccumulateRunningSum ( const QwVQWK_Channel & value,
Int_t ErrorMask = 0xFFFFFFF )
inline

Definition at line 203 of file QwVQWK_Channel.h.

203 {
204 AccumulateRunningSum(value, -1, ErrorMask);
205 };

References AccumulateRunningSum(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ DeaccumulateRunningSum() [2/2]

virtual void VQwHardwareChannel::DeaccumulateRunningSum ( const VQwHardwareChannel * value,
Int_t ErrorMask = 0xFFFFFFF )
inlinevirtual

Reimplemented from VQwHardwareChannel.

Definition at line 248 of file VQwHardwareChannel.h.

248 {
249 AccumulateRunningSum(value, -1, ErrorMask);
250 };

◆ Difference()

void QwVQWK_Channel::Difference ( const QwVQWK_Channel & value1,
const QwVQWK_Channel & value2 )

Definition at line 1373 of file QwVQWK_Channel.cc.

1374{
1375 *this = value1;
1376 *this -= value2;
1377}

References QwVQWK_Channel().

Referenced by QwQPD::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ DivideBy() [1/2]

void QwVQWK_Channel::DivideBy ( const QwVQWK_Channel & denom)

Definition at line 1514 of file QwVQWK_Channel.cc.

1515{
1516 *this /= denom;
1517}

References QwVQWK_Channel().

+ Here is the call graph for this function:

◆ DivideBy() [2/2]

void QwVQWK_Channel::DivideBy ( const VQwHardwareChannel * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1221 of file QwVQWK_Channel.cc.

1222{
1223 const QwVQWK_Channel* tmpptr;
1224 tmpptr = dynamic_cast<const QwVQWK_Channel*>(valueptr);
1225 if (tmpptr!=NULL){
1226 *this /= *tmpptr;
1227 } else {
1228 TString loc="Standard exception from QwVQWK_Channel::DivideBy = "
1229 +valueptr->GetElementName()+" is an incompatible type.";
1230 throw std::invalid_argument(loc.Data());
1231 }
1232}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ EncodeEventData()

void QwVQWK_Channel::EncodeEventData ( std::vector< UInt_t > & buffer)
overridevirtual

Encode the event data into a CODA buffer.

Implements MQwMockable.

Definition at line 443 of file QwVQWK_Channel.cc.

444{
445 Long_t localbuf[6] = {0};
446
447 if (IsNameEmpty()) {
448 // This channel is not used, but is present in the data stream.
449 // Skip over this data.
450 } else {
451 // localbuf[4] = 0;
452 for (Int_t i = 0; i < 4; i++) {
453 localbuf[i] = fBlock_raw[i];
454 // localbuf[4] += localbuf[i]; // fHardwareBlockSum_raw
455 }
456 // The following causes many rounding errors and skips due to the check
457 // that fHardwareBlockSum_raw == fSoftwareBlockSum_raw in IsGoodEvent().
458 localbuf[4] = fHardwareBlockSum_raw;
459 localbuf[5] = (fNumberOfSamples << 16 & 0xFFFF0000)
460 | (fSequenceNumber << 8 & 0x0000FF00);
461
462 for (Int_t i = 0; i < 6; i++){
463 buffer.push_back(localbuf[i]);
464 }
465 }
466}

References fBlock_raw, fHardwareBlockSum_raw, fNumberOfSamples, fSequenceNumber, and VQwDataElement::IsNameEmpty().

+ Here is the call graph for this function:

◆ FillHistograms()

void QwVQWK_Channel::FillHistograms ( )
overridevirtual

Fill the histograms for this data element.

Implements VQwDataElement.

Definition at line 685 of file QwVQWK_Channel.cc.

686{
687 Int_t index=0;
688
689 if (IsNameEmpty())
690 {
691 // This channel is not used, so skip creating the histograms.
692 } else
693 {
694 if(fDataToSave==kRaw)
695 {
696 for (Int_t i=0; i<fBlocksPerEvent; i++)
697 {
698 if (fHistograms[index] != NULL && (fErrorFlag)==0)
699 fHistograms[index]->Fill(this->GetRawBlockValue(i));
700 if (fHistograms[index+1] != NULL && (fErrorFlag)==0)
701 fHistograms[index+1]->Fill(this->GetBlockValue(i));
702 index+=2;
703 }
704 if (fHistograms[index] != NULL && (fErrorFlag)==0)
705 fHistograms[index]->Fill(this->GetRawHardwareSum());
706 if (fHistograms[index+1] != NULL && (fErrorFlag)==0)
707 fHistograms[index+1]->Fill(this->GetHardwareSum());
708 index+=2;
709 if (fHistograms[index] != NULL && (fErrorFlag)==0)
710 fHistograms[index]->Fill(this->GetRawSoftwareSum()-this->GetRawHardwareSum());
711 }
712 else if(fDataToSave==kDerived)
713 {
714 for (Int_t i=0; i<fBlocksPerEvent; i++)
715 {
716 if (fHistograms[index] != NULL && (fErrorFlag)==0)
717 fHistograms[index]->Fill(this->GetBlockValue(i));
718 index+=1;
719 }
720 if (fHistograms[index] != NULL && (fErrorFlag)==0)
721 fHistograms[index]->Fill(this->GetHardwareSum());
722 index+=1;
723 if (fHistograms[index] != NULL){
725 fHistograms[index]->Fill(kErrorFlag_sample);
727 fHistograms[index]->Fill(kErrorFlag_SW_HW);
729 fHistograms[index]->Fill(kErrorFlag_Sequence);
731 fHistograms[index]->Fill(kErrorFlag_ZeroHW);
733 fHistograms[index]->Fill(kErrorFlag_VQWK_Sat);
735 fHistograms[index]->Fill(kErrorFlag_SameHW);
736 }
737
738 }
739
740 }
741}
Int_t GetRawBlockValue(size_t blocknum) const
Double_t GetBlockValue(size_t blocknum) const

References fBlocksPerEvent, VQwHardwareChannel::fDataToSave, VQwDataElement::fErrorFlag, MQwHistograms::fHistograms, GetBlockValue(), GetHardwareSum(), GetRawBlockValue(), GetRawHardwareSum(), GetRawSoftwareSum(), VQwDataElement::IsNameEmpty(), VQwDataElement::kDerived, kErrorFlag_SameHW, kErrorFlag_sample, kErrorFlag_Sequence, kErrorFlag_SW_HW, kErrorFlag_VQWK_Sat, kErrorFlag_ZeroHW, and VQwDataElement::kRaw.

+ Here is the call graph for this function:

◆ FillTreeVector()

void QwVQWK_Channel::FillTreeVector ( QwRootTreeBranchVector & values) const
overridevirtual

Implements VQwHardwareChannel.

Definition at line 842 of file QwVQWK_Channel.cc.

843{
844 if (IsNameEmpty()) {
845 // This channel is not used, so skip filling the tree vector.
846 } else if (fTreeArrayNumEntries <= 0) {
847 if (bDEBUG) std::cerr << "QwVQWK_Channel::FillTreeVector: fTreeArrayNumEntries=="
848 << fTreeArrayNumEntries << std::endl;
849 } else if (values.size() < fTreeArrayIndex+fTreeArrayNumEntries){
850 if (bDEBUG) std::cerr << "QwVQWK_Channel::FillTreeVector: values.size()=="
851 << values.size()
852 << "; fTreeArrayIndex+fTreeArrayNumEntries=="
854 << std::endl;
855 } else {
856
857 UInt_t index = fTreeArrayIndex;
858
859 // hw_sum
860 if (bHw_sum) {
861 //values.SetValue(fTreeArrayIndex, "hw_sum"_h32, this->GetHardwareSum());
862 values.SetValue(index++, this->GetHardwareSum());
863 if (fDataToSave == kMoments) {
864 values.SetValue(index++, this->GetHardwareSumM2());
865 values.SetValue(index++, this->GetHardwareSumError());
866 }
867 }
868
869 if (bBlock) {
870 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
871 // blocki
872 values.SetValue(index++, this->GetBlockValue(i));
873 }
874 }
875
876 // num_samples
877 if (bNum_samples)
878 values.SetValue(index++, (fDataToSave == kMoments)? this->fGoodEventCount: this->fNumberOfSamples);
879
880 // Device_Error_Code
882 values.SetValue(index++, this->fErrorFlag);
883
884 if (fDataToSave == kRaw)
885 {
886 // hw_sum_raw
887 if (bHw_sum_raw)
888 values.SetValue(index++, this->GetRawHardwareSum());
889
890 if (bBlock_raw) {
891 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
892 // blocki_raw
893 values.SetValue(index++, this->GetRawBlockValue(i));
894 }
895 }
896
897 // sequence_number
899 values.SetValue(index++, this->fSequenceNumber);
900 }
901 }
902}
void SetValue(size_type index, Double_t val)
Definition QwRootFile.h:110
Double_t GetHardwareSumM2() const
Double_t GetHardwareSumError() const

References bBlock, bBlock_raw, bDEBUG, bDevice_Error_Code, bHw_sum, bHw_sum_raw, bNum_samples, bSequence_number, fBlocksPerEvent, VQwHardwareChannel::fDataToSave, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fNumberOfSamples, fSequenceNumber, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, GetBlockValue(), GetHardwareSum(), GetHardwareSumError(), GetHardwareSumM2(), GetRawBlockValue(), GetRawHardwareSum(), VQwDataElement::IsNameEmpty(), VQwDataElement::kMoments, VQwDataElement::kRaw, QwRootTreeBranchVector::SetValue(), and QwRootTreeBranchVector::size().

+ Here is the call graph for this function:

◆ ForceMapfileSampleSize()

void QwVQWK_Channel::ForceMapfileSampleSize ( )
inline

Forces the event "number of samples" variable to be what was expected from the mapfile. NOTE: this should only be used in mock data generation!

Definition at line 137 of file QwVQWK_Channel.h.

References fNumberOfSamples, and fNumberOfSamples_map.

◆ GetAverageVolts()

Double_t QwVQWK_Channel::GetAverageVolts ( ) const

Definition at line 598 of file QwVQWK_Channel.cc.

599{
600 //Double_t avgVolts = (fBlock[0]+fBlock[1]+fBlock[2]+fBlock[3])*kVQWK_VoltsPerBit/fNumberOfSamples;
602 //std::cout<<"QwVQWK_Channel::GetAverageVolts() = "<<avgVolts<<std::endl;
603 return avgVolts;
604
605}

References fHardwareBlockSum, fNumberOfSamples, and kVQWK_VoltsPerBit.

◆ GetBlockErrorValue()

Double_t QwVQWK_Channel::GetBlockErrorValue ( size_t blocknum) const
inlineprivate

Definition at line 304 of file QwVQWK_Channel.h.

304{ return GetValueError(blocknum+1);};
Double_t GetValueError() const

References GetValueError().

Referenced by PrintValue().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetBlockValue()

Double_t QwVQWK_Channel::GetBlockValue ( size_t blocknum) const
inlineprivate

Definition at line 303 of file QwVQWK_Channel.h.

303{ return GetValue(blocknum+1);};
Double_t GetValue() const

References GetValue().

Referenced by FillHistograms(), FillTreeVector(), and PrintValue().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetBufferOffset()

Int_t QwVQWK_Channel::GetBufferOffset ( Int_t moduleindex,
Int_t channelindex )
static

Class: QwVQWK_Channel Base class containing decoding functions for the VQWK_Channel 6 32-bit datawords. The functions in this class will decode a single channel worth of VQWK_Channel data and provide the components through member functions.

Static member function to return the word offset within a data buffer given the module number index and the channel number index.

Parameters
moduleindexModule index within this buffer; counts from zero
channelindexChannel index within this module; counts from zero
Returns
The number of words offset to the beginning of this channel's data from the beginning of the VQWK buffer.

Definition at line 44 of file QwVQWK_Channel.cc.

44 {
45 Int_t offset = -1;
46 if (moduleindex<0 ){
47 QwError << "QwVQWK_Channel::GetBufferOffset: Invalid module index,"
48 << moduleindex
49 << ". Must be zero or greater."
50 << QwLog::endl;
51 } else if (channelindex<0 || channelindex>kMaxChannels){
52 QwError << "QwVQWK_Channel::GetBufferOffset: Invalid channel index,"
53 << channelindex
54 << ". Must be in range [0," << kMaxChannels << "]."
55 << QwLog::endl;
56 } else {
57 offset = ( (moduleindex * kMaxChannels) + channelindex )
59 }
60 return offset;
61 }
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
static const Int_t kWordsPerChannel
static const Int_t kMaxChannels

References QwLog::endl(), kMaxChannels, kWordsPerChannel, and QwError.

Referenced by QwBeamDetectorID::QwBeamDetectorID(), and QwModChannelID::QwModChannelID().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetHardwareSum()

Double_t QwVQWK_Channel::GetHardwareSum ( ) const
inlineprivate

Definition at line 306 of file QwVQWK_Channel.h.

306{ return GetValue(0);};

References GetValue().

Referenced by ApplySingleEventCuts(), ApplySingleEventCuts(), FillHistograms(), FillTreeVector(), operator<<, and PrintValue().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetHardwareSumError()

Double_t QwVQWK_Channel::GetHardwareSumError ( ) const
inlineprivate

Definition at line 309 of file QwVQWK_Channel.h.

309{ return GetValueError(0); };

References GetValueError().

Referenced by FillTreeVector(), and PrintValue().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetHardwareSumM2()

Double_t QwVQWK_Channel::GetHardwareSumM2 ( ) const
inlineprivate

Definition at line 307 of file QwVQWK_Channel.h.

307{ return GetValueM2(0); };
Double_t GetValueM2() const

References GetValueM2().

Referenced by FillTreeVector().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetHardwareSumWidth()

Double_t QwVQWK_Channel::GetHardwareSumWidth ( ) const
inlineprivate

Definition at line 308 of file QwVQWK_Channel.h.

308{ return GetValueWidth(0); };

References GetValueWidth().

Referenced by PrintValue().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetNumberOfSamples()

size_t QwVQWK_Channel::GetNumberOfSamples ( ) const
inline

Definition at line 276 of file QwVQWK_Channel.h.

276{return (fNumberOfSamples);};

References fNumberOfSamples.

Referenced by ApplyHWChecks().

+ Here is the caller graph for this function:

◆ GetRawBlockValue()

Int_t QwVQWK_Channel::GetRawBlockValue ( size_t blocknum) const
inlineprivate

Definition at line 312 of file QwVQWK_Channel.h.

312{return GetRawValue(blocknum+1);};
Int_t GetRawValue() const

References GetRawValue().

Referenced by FillHistograms(), and FillTreeVector().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetRawHardwareSum()

Int_t QwVQWK_Channel::GetRawHardwareSum ( ) const
inlineprivate

Definition at line 313 of file QwVQWK_Channel.h.

313{ return GetRawValue(0);};

References GetRawValue().

Referenced by ApplyHWChecks(), FillHistograms(), and FillTreeVector().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetRawSoftwareSum()

Int_t QwVQWK_Channel::GetRawSoftwareSum ( ) const
inlineprivate

Definition at line 314 of file QwVQWK_Channel.h.

314{return fSoftwareBlockSum_raw;};

References fSoftwareBlockSum_raw.

Referenced by ApplyHWChecks(), and FillHistograms().

+ Here is the caller graph for this function:

◆ GetRawValue() [1/2]

Int_t VQwHardwareChannel::GetRawValue ( ) const
inline

Definition at line 125 of file VQwHardwareChannel.h.

125{return this->GetRawValue(0);};

Referenced by GetRawBlockValue(), and GetRawHardwareSum().

+ Here is the caller graph for this function:

◆ GetRawValue() [2/2]

Int_t QwVQWK_Channel::GetRawValue ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 251 of file QwVQWK_Channel.h.

251 {
252 RangeCheck(element);
253 if (element==0) return fHardwareBlockSum_raw;
254 return fBlock_raw[element-1];
255 }
void RangeCheck(size_t element) const
Checks that the requested element is in range, to be used in accesses to subelements similar to std::...

References fBlock_raw, fHardwareBlockSum_raw, and VQwHardwareChannel::RangeCheck().

+ Here is the call graph for this function:

◆ GetSequenceNumber()

size_t QwVQWK_Channel::GetSequenceNumber ( ) const
inline

Definition at line 275 of file QwVQWK_Channel.h.

275{return (fSequenceNumber);};

References fSequenceNumber.

Referenced by ApplyHWChecks().

+ Here is the caller graph for this function:

◆ GetValue() [1/2]

Double_t VQwHardwareChannel::GetValue ( ) const
inline

Definition at line 126 of file VQwHardwareChannel.h.

126{return this->GetValue(0);};

Referenced by GetBlockValue(), and GetHardwareSum().

+ Here is the caller graph for this function:

◆ GetValue() [2/2]

Double_t QwVQWK_Channel::GetValue ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 256 of file QwVQWK_Channel.h.

256 {
257 RangeCheck(element);
258 if (element==0) return fHardwareBlockSum;
259 return fBlock[element-1];
260 }

References fBlock, fHardwareBlockSum, and VQwHardwareChannel::RangeCheck().

Referenced by QwQPD::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ GetValueError() [1/2]

Double_t VQwHardwareChannel::GetValueError ( ) const
inline

Definition at line 128 of file VQwHardwareChannel.h.

128{return this->GetValueError(0);};

Referenced by GetBlockErrorValue(), and GetHardwareSumError().

+ Here is the caller graph for this function:

◆ GetValueError() [2/2]

Double_t QwVQWK_Channel::GetValueError ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 266 of file QwVQWK_Channel.h.

266 {
267 RangeCheck(element);
268 if (element==0) return fHardwareBlockSumError;
269 return fBlockError[element-1];
270 }

References fBlockError, fHardwareBlockSumError, and VQwHardwareChannel::RangeCheck().

+ Here is the call graph for this function:

◆ GetValueM2() [1/2]

Double_t VQwHardwareChannel::GetValueM2 ( ) const
inline

Definition at line 127 of file VQwHardwareChannel.h.

127{return this->GetValueM2(0);};

Referenced by GetHardwareSumM2().

+ Here is the caller graph for this function:

◆ GetValueM2() [2/2]

Double_t QwVQWK_Channel::GetValueM2 ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 261 of file QwVQWK_Channel.h.

261 {
262 RangeCheck(element);
263 if (element==0) return fHardwareBlockSumM2;
264 return fBlockM2[element-1];
265 }

References fBlockM2, fHardwareBlockSumM2, and VQwHardwareChannel::RangeCheck().

+ Here is the call graph for this function:

◆ GetValueWidth() [1/2]

Double_t VQwHardwareChannel::GetValueWidth ( ) const
inline

Definition at line 129 of file VQwHardwareChannel.h.

129{return this->GetValueWidth(0);};

Referenced by CalculateRunningAverage(), and GetHardwareSumWidth().

+ Here is the caller graph for this function:

◆ GetValueWidth() [2/2]

Double_t VQwHardwareChannel::GetValueWidth ( size_t element) const
inline

Definition at line 134 of file VQwHardwareChannel.h.

134 {
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 };

◆ GetVQWKSaturationLimt()

Double_t QwVQWK_Channel::GetVQWKSaturationLimt ( )
inline

Definition at line 227 of file QwVQWK_Channel.h.

227 {//Get the absolute staturation limit in volts.
228 return fSaturationABSLimit;
229 }

References fSaturationABSLimit.

Referenced by ApplyHWChecks().

+ Here is the caller graph for this function:

◆ IncrementErrorCounters()

void QwVQWK_Channel::IncrementErrorCounters ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 141 of file QwVQWK_Channel.cc.

141 {
143 fErrorCount_sample++; //increment the hw error counter
145 fErrorCount_SW_HW++; //increment the hw error counter
147 fErrorCount_Sequence++; //increment the hw error counter
149 fErrorCount_SameHW++; //increment the hw error counter
151 fErrorCount_ZeroHW++; //increment the hw error counter
153 fErrorCount_HWSat++; //increment the hw saturation error counter
156 fNumEvtsWithEventCutsRejected++; //increment the event cut error counter
157 }
158}
Int_t fErrorCount_ZeroHW
check to see ADC returning zero
Int_t fErrorCount_sample
for sample size check
Int_t fErrorCount_Sequence
sequence number check
Int_t fErrorCount_SameHW
check to see ADC returning same HW value
Int_t fErrorCount_HWSat
check to see ADC channel is saturated
Int_t fErrorCount_SW_HW
HW_sum==SW_sum check.
Int_t fNumEvtsWithEventCutsRejected
Counts the Event cut rejected events.

References fErrorCount_HWSat, fErrorCount_SameHW, fErrorCount_sample, fErrorCount_Sequence, fErrorCount_SW_HW, fErrorCount_ZeroHW, VQwDataElement::fErrorFlag, fNumEvtsWithEventCutsRejected, kErrorFlag_EventCut_L, kErrorFlag_EventCut_U, kErrorFlag_SameHW, kErrorFlag_sample, kErrorFlag_Sequence, kErrorFlag_SW_HW, kErrorFlag_VQWK_Sat, and kErrorFlag_ZeroHW.

◆ InitializeChannel() [1/2]

void QwVQWK_Channel::InitializeChannel ( TString name,
TString datatosave )
overridevirtual

Initialize the fields in this object.

Implements VQwHardwareChannel.

Definition at line 162 of file QwVQWK_Channel.cc.

163{
164 SetElementName(name);
165 SetDataToSave(datatosave);
168
169 kFoundPedestal = 0;
170 kFoundGain = 0;
171
172 fPedestal = 0.0;
173 fCalibrationFactor = 1.0;
174
175 fBlocksPerEvent = 4;
176
177 fTreeArrayIndex = 0;
179
181
185
186 // Use internal random variable by default
188
189 // Mock drifts
190 fMockDriftAmplitude.clear();
191 fMockDriftFrequency.clear();
192 fMockDriftPhase.clear();
193
194 // Mock asymmetries
195 fMockAsymmetry = 0.0;
196 fMockGaussianMean = 0.0;
197 fMockGaussianSigma = 0.0;
198
199 // Event cuts
200 fULimit=-1;
201 fLLimit=1;
203
204 fErrorFlag=0; //Initialize the error flag
205 fErrorConfigFlag=0; //Initialize the error config. flag
206
207 //init error counters//
214
219
220 fGoodEventCount = 0;
221
222 bEVENTCUTMODE = 0;
223
224 //std::cout<< "name = "<<name<<" error count same _HW = "<<fErrorCount_SameHW <<std::endl;
225 return;
226}
std::vector< Double_t > fMockDriftAmplitude
Harmonic drift amplitude.
std::vector< Double_t > fMockDriftFrequency
Harmonic drift frequency.
Double_t fMockAsymmetry
Helicity asymmetry.
bool fUseExternalRandomVariable
Flag to use an externally provided normal random variable.
Double_t fMockGaussianSigma
Sigma of normal distribution.
Double_t fMockGaussianMean
Mean of normal distribution.
std::vector< Double_t > fMockDriftPhase
Harmonic drift phase.
void ClearEventData() override
Clear the event data in this element.
UInt_t fPreviousSequenceNumber
Previous event sequence number for this channel.
void SetElementName(const TString &name)
Set the name of this element.
void SetDataToSave(TString datatosave)
Set the flag indicating if raw or derived values are in this data element.
void SetNumberOfSubElements(const size_t elements)
Set the number of data words in this data element.
void SetNumberOfDataWords(const UInt_t &numwords)
Set the number of data words in this data element.

References VQwHardwareChannel::bEVENTCUTMODE, ClearEventData(), fADC_Same_NumEvt, fBlocksPerEvent, VQwHardwareChannel::fCalibrationFactor, VQwDataElement::fErrorConfigFlag, fErrorCount_HWSat, fErrorCount_SameHW, fErrorCount_sample, fErrorCount_Sequence, fErrorCount_SW_HW, fErrorCount_ZeroHW, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, VQwHardwareChannel::fLLimit, MQwMockable::fMockAsymmetry, MQwMockable::fMockDriftAmplitude, MQwMockable::fMockDriftFrequency, MQwMockable::fMockDriftPhase, MQwMockable::fMockGaussianMean, MQwMockable::fMockGaussianSigma, fNumberOfSamples, fNumberOfSamples_map, fNumEvtsWithEventCutsRejected, VQwHardwareChannel::fPedestal, fPrev_HardwareBlockSum, fPreviousSequenceNumber, fSequenceNo_Counter, fSequenceNo_Prev, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, VQwHardwareChannel::fULimit, MQwMockable::fUseExternalRandomVariable, VQwHardwareChannel::kFoundGain, VQwHardwareChannel::kFoundPedestal, VQwHardwareChannel::SetDataToSave(), VQwDataElement::SetElementName(), VQwHardwareChannel::SetNumberOfDataWords(), and VQwHardwareChannel::SetNumberOfSubElements().

Referenced by InitializeChannel(), QwLinearDiodeArray::ProcessEvent(), QwQPD::ProcessEvent(), QwVQWK_Channel(), and QwVQWK_Channel().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ InitializeChannel() [2/2]

void QwVQWK_Channel::InitializeChannel ( TString subsystem,
TString instrumenttype,
TString name,
TString datatosave )
overridevirtual

Initialize the fields in this object.

Implements VQwHardwareChannel.

Definition at line 230 of file QwVQWK_Channel.cc.

230 {
231 InitializeChannel(name,datatosave);
232 SetSubsystemName(subsystem);
233 SetModuleType(instrumenttype);
234 //PrintInfo();
235}
void SetSubsystemName(TString sysname)
Set the name of the inheriting subsystem name.
void SetModuleType(TString ModuleType)
set the type of the beam instrument

References InitializeChannel(), VQwDataElement::SetModuleType(), and VQwDataElement::SetSubsystemName().

+ Here is the call graph for this function:

◆ LoadChannelParameters()

void QwVQWK_Channel::LoadChannelParameters ( QwParameterFile & paramfile)
overridevirtual

Reimplemented from VQwDataElement.

Definition at line 237 of file QwVQWK_Channel.cc.

237 {
238 UInt_t value = 0;
239 if (paramfile.ReturnValue("sample_size",value)){
241 } else {
242 QwWarning << "VQWK Channel "
243 << GetElementName()
244 << " cannot set the default sample size."
245 << QwLog::endl;
246 }
247};
Bool_t ReturnValue(const std::string keyname, T &retvalue)
void SetDefaultSampleSize(size_t num_samples_map)

References QwLog::endl(), VQwDataElement::GetElementName(), QwWarning, QwParameterFile::ReturnValue(), and SetDefaultSampleSize().

+ Here is the call graph for this function:

◆ LoadMockDataParameters()

void MQwMockable::LoadMockDataParameters ( QwParameterFile & paramfile)
virtual

Load the mock data parameters from the current line in the param file.

Reimplemented from VQwDataElement.

Definition at line 60 of file MQwMockable.cc.

19 {
20 Bool_t ldebug=kFALSE;
21 Double_t asym=0.0, mean=0.0, sigma=0.0;
22 Double_t amplitude=0.0, phase=0.0, frequency=0.0;
23
24 //Check to see if this line contains "drift"
25 if (paramfile.GetLine().find("drift")!=std::string::npos){
26 // "drift" appears somewhere in the line. Assume it is the next token and move on...
27 paramfile.GetNextToken(); //Throw away this token. Now the next three should be the drift parameters.
28 //read 3 parameters
29 amplitude = paramfile.GetTypedNextToken<Double_t>();
30 phase = paramfile.GetTypedNextToken<Double_t>();
31 frequency = paramfile.GetTypedNextToken<Double_t>(); // The frequency we read in should be in Hz.
32 this->AddRandomEventDriftParameters(amplitude, phase, frequency*Qw::Hz);
33 // std::cout << "In MQwMockable::LoadMockDataParameters: amp = " << amplitude << "\t phase = " << phase << "\t freq = " << frequency << std::endl;
34 }
35 else {
36 asym = paramfile.GetTypedNextToken<Double_t>();
37 mean = paramfile.GetTypedNextToken<Double_t>();
38 sigma = paramfile.GetTypedNextToken<Double_t>();
39 if (ldebug==1) {
40 std::cout << "#################### \n";
41 std::cout << "asym, mean, sigma \n" << std::endl;
42 std::cout << asym << " / "
43 << mean << " / "
44 << sigma << " / "
45 << std::endl;
46 }
47 this->SetRandomEventParameters(mean, sigma);
48 this->SetRandomEventAsymmetry(asym);
49 }
50}
static const double Hz
Frequency units: base unit is kHz.
Definition QwUnits.h:88
void AddRandomEventDriftParameters(Double_t amplitude, Double_t phase, Double_t frequency)
Add drift parameters to the internal set.
void SetRandomEventAsymmetry(Double_t asymmetry)
Set the helicity asymmetry.
void SetRandomEventParameters(Double_t mean, Double_t sigma)
Set the normal random event parameters.
T GetTypedNextToken()
Get next token into specific type.
std::string GetLine()
std::string GetNextToken(const std::string &separatorchars)
Get next token as a string.

◆ MatchNumberOfSamples()

Bool_t QwVQWK_Channel::MatchNumberOfSamples ( size_t numsamp)

Definition at line 1842 of file QwVQWK_Channel.cc.

1843{
1844 Bool_t status = kTRUE;
1845 if (!IsNameEmpty()){
1846 status = (fNumberOfSamples==numsamp);
1847 if (! status){
1848 if (bDEBUG)
1849 std::cerr << "QwVQWK_Channel::MatchNumberOfSamples: Channel "
1850 << GetElementName()
1851 << " had fNumberOfSamples==" << fNumberOfSamples
1852 << " and was supposed to have " << numsamp
1853 << std::endl;
1854 // PrintChannel();
1855 }
1856 }
1857 return status;
1858}

References bDEBUG, fNumberOfSamples, VQwDataElement::GetElementName(), and VQwDataElement::IsNameEmpty().

Referenced by ApplyHWChecks().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ MatchSequenceNumber()

Bool_t QwVQWK_Channel::MatchSequenceNumber ( size_t seqnum)

Definition at line 1832 of file QwVQWK_Channel.cc.

1833{
1834
1835 Bool_t status = kTRUE;
1836 if (!IsNameEmpty()){
1837 status = (fSequenceNumber==seqnum);
1838 }
1839 return status;
1840}

References fSequenceNumber, and VQwDataElement::IsNameEmpty().

Referenced by ApplyHWChecks().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ MultiplyBy()

void QwVQWK_Channel::MultiplyBy ( const VQwHardwareChannel * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1209 of file QwVQWK_Channel.cc.

1210{
1211 const QwVQWK_Channel* tmpptr;
1212 tmpptr = dynamic_cast<const QwVQWK_Channel*>(valueptr);
1213 if (tmpptr!=NULL){
1214 *this *= *tmpptr;
1215 } else {
1216 TString loc="Standard exception from QwVQWK_Channel::MultiplyBy = "
1217 +valueptr->GetElementName()+" is an incompatible type.";
1218 throw std::invalid_argument(loc.Data());
1219 }
1220}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ operator*()

const QwVQWK_Channel QwVQWK_Channel::operator* ( const QwVQWK_Channel & value) const

Definition at line 1285 of file QwVQWK_Channel.cc.

1286{
1287 QwVQWK_Channel result = *this;
1288 result *= value;
1289 return result;
1290}

References QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator*=() [1/2]

QwVQWK_Channel & QwVQWK_Channel::operator*= ( const QwVQWK_Channel & value)

Definition at line 1292 of file QwVQWK_Channel.cc.

1293{
1294 if (!IsNameEmpty()){
1295 for (Int_t i=0; i<fBlocksPerEvent; i++){
1296 this->fBlock[i] *= value.fBlock[i];
1297 this->fBlockM2[i] = 0.0;
1298 }
1299 this->fHardwareBlockSum *= value.fHardwareBlockSum;
1300 this->fHardwareBlockSumM2 = 0.0;
1301 this->fNumberOfSamples *= value.fNumberOfSamples;
1302 this->fSequenceNumber = 0;
1303 this->fErrorFlag |= (value.fErrorFlag);
1304 }
1305
1306 return *this;
1307}

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator*=() [2/2]

VQwHardwareChannel & QwVQWK_Channel::operator*= ( const VQwHardwareChannel & input)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1337 of file QwVQWK_Channel.cc.

1338{
1339 const QwVQWK_Channel* tmpptr;
1340 tmpptr = dynamic_cast<const QwVQWK_Channel*>(&source);
1341 if (tmpptr!=NULL){
1342 *this *= *tmpptr;
1343 } else {
1344 TString loc="Standard exception from QwVQWK_Channel::operator*= "
1345 +source.GetElementName()+" "
1346 +this->GetElementName()+" are not of the same type";
1347 throw(std::invalid_argument(loc.Data()));
1348 }
1349 return *this;
1350}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ operator+()

const QwVQWK_Channel QwVQWK_Channel::operator+ ( const QwVQWK_Channel & value) const

Definition at line 1235 of file QwVQWK_Channel.cc.

1236{
1237 QwVQWK_Channel result = *this;
1238 result += value;
1239 return result;
1240}

References QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator+=() [1/2]

QwVQWK_Channel & QwVQWK_Channel::operator+= ( const QwVQWK_Channel & value)

Definition at line 1242 of file QwVQWK_Channel.cc.

1243{
1244
1245 if (!IsNameEmpty()) {
1246 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
1247 this->fBlock[i] += value.fBlock[i];
1248 this->fBlockM2[i] = 0.0;
1249 }
1250 this->fHardwareBlockSum += value.fHardwareBlockSum;
1251 this->fHardwareBlockSumM2 = 0.0;
1252 this->fNumberOfSamples += value.fNumberOfSamples;
1253 this->fSequenceNumber = 0;
1254 this->fErrorFlag |= (value.fErrorFlag);
1255
1256 }
1257
1258 return *this;
1259}

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator+=() [2/2]

VQwHardwareChannel & QwVQWK_Channel::operator+= ( const VQwHardwareChannel & input)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1309 of file QwVQWK_Channel.cc.

1310{
1311 const QwVQWK_Channel* tmpptr;
1312 tmpptr = dynamic_cast<const QwVQWK_Channel*>(&source);
1313 if (tmpptr!=NULL){
1314 *this += *tmpptr;
1315 } else {
1316 TString loc="Standard exception from QwVQWK_Channel::operator+= "
1317 +source.GetElementName()+" "
1318 +this->GetElementName()+" are not of the same type";
1319 throw(std::invalid_argument(loc.Data()));
1320 }
1321 return *this;
1322}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ operator-()

const QwVQWK_Channel QwVQWK_Channel::operator- ( const QwVQWK_Channel & value) const

Definition at line 1261 of file QwVQWK_Channel.cc.

1262{
1263 QwVQWK_Channel result = *this;
1264 result -= value;
1265 return result;
1266}

References QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator-=() [1/2]

QwVQWK_Channel & QwVQWK_Channel::operator-= ( const QwVQWK_Channel & value)

Definition at line 1268 of file QwVQWK_Channel.cc.

1269{
1270 if (!IsNameEmpty()){
1271 for (Int_t i=0; i<fBlocksPerEvent; i++){
1272 this->fBlock[i] -= value.fBlock[i];
1273 this->fBlockM2[i] = 0.0;
1274 }
1275 this->fHardwareBlockSum -= value.fHardwareBlockSum;
1276 this->fHardwareBlockSumM2 = 0.0;
1277 this->fNumberOfSamples += value.fNumberOfSamples;
1278 this->fSequenceNumber = 0;
1279 this->fErrorFlag |= (value.fErrorFlag);
1280 }
1281
1282 return *this;
1283}

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ operator-=() [2/2]

VQwHardwareChannel & QwVQWK_Channel::operator-= ( const VQwHardwareChannel & input)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1323 of file QwVQWK_Channel.cc.

1324{
1325 const QwVQWK_Channel* tmpptr;
1326 tmpptr = dynamic_cast<const QwVQWK_Channel*>(&source);
1327 if (tmpptr!=NULL){
1328 *this -= *tmpptr;
1329 } else {
1330 TString loc="Standard exception from QwVQWK_Channel::operator-= "
1331 +source.GetElementName()+" "
1332 +this->GetElementName()+" are not of the same type";
1333 throw(std::invalid_argument(loc.Data()));
1334 }
1335 return *this;
1336}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ operator/=() [1/2]

QwVQWK_Channel & QwVQWK_Channel::operator/= ( const QwVQWK_Channel & value)
protected

Definition at line 1393 of file QwVQWK_Channel.cc.

1394{
1395 // In this function, leave the "raw" variables untouched.
1396 //
1397 Double_t ratio;
1398 Double_t variance;
1399 if (!IsNameEmpty()) {
1400 // The variances are calculated using the following formula:
1401 // Var[ratio] = ratio^2 (Var[numer] / numer^2 + Var[denom] / denom^2)
1402 //
1403 // This requires that both the numerator and denominator are non-zero!
1404 //
1405 for (Int_t i = 0; i < 4; i++) {
1406 if (this->fBlock[i] != 0.0 && denom.fBlock[i] != 0.0){
1407 ratio = (this->fBlock[i]) / (denom.fBlock[i]);
1408 variance = ratio * ratio *
1409 (this->fBlockM2[i] / this->fBlock[i] / this->fBlock[i]
1410 + denom.fBlockM2[i] / denom.fBlock[i] / denom.fBlock[i]);
1411 fBlock[i] = ratio;
1412 fBlockM2[i] = variance;
1413 } else if (this->fBlock[i] == 0.0) {
1414 this->fBlock[i] = 0.0;
1415 this->fBlockM2[i] = 0.0;
1416 } else {
1417 QwVerbose << "Attempting to divide by zero block in "
1419 fBlock[i] = 0.0;
1420 fBlockM2[i] = 0.0;
1421 }
1422 }
1423 if (this->fHardwareBlockSum != 0.0 && denom.fHardwareBlockSum != 0.0){
1424 ratio = (this->fHardwareBlockSum) / (denom.fHardwareBlockSum);
1425 variance = ratio * ratio *
1427 + denom.fHardwareBlockSumM2 / denom.fHardwareBlockSum / denom.fHardwareBlockSum);
1428 fHardwareBlockSum = ratio;
1429 fHardwareBlockSumM2 = variance;
1430 } else if (this->fHardwareBlockSum == 0.0) {
1431 fHardwareBlockSum = 0.0;
1432 fHardwareBlockSumM2 = 0.0;
1433 } else {
1434 QwVerbose << "Attempting to divide by zero sum in "
1436 fHardwareBlockSumM2 = 0.0;
1437 }
1438 // Remaining variables
1439 // Don't change fNumberOfSamples, fSequenceNumber, fGoodEventCount,
1440 // 'OR' the HW error codes in the fErrorFlag values together.
1441 fErrorFlag |= (denom.fErrorFlag);//mix only the hardware error codes
1442 }
1443
1444 // Nanny
1446 QwWarning << "Angry Nanny: NaN detected in " << GetElementName() << QwLog::endl;
1447
1448 return *this;
1449}
#define QwVerbose
Predefined log drain for verbose messages.
Definition QwLog.h:54

References QwLog::endl(), fBlock, fBlockM2, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), QwVerbose, QwVQWK_Channel(), and QwWarning.

+ Here is the call graph for this function:

◆ operator/=() [2/2]

VQwHardwareChannel & QwVQWK_Channel::operator/= ( const VQwHardwareChannel & input)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1351 of file QwVQWK_Channel.cc.

1352{
1353 const QwVQWK_Channel* tmpptr;
1354 tmpptr = dynamic_cast<const QwVQWK_Channel*>(&source);
1355 if (tmpptr!=NULL){
1356 *this /= *tmpptr;
1357 } else {
1358 TString loc="Standard exception from QwVQWK_Channel::operator/= "
1359 +source.GetElementName()+" "
1360 +this->GetElementName()+" are not of the same type";
1361 throw(std::invalid_argument(loc.Data()));
1362 }
1363 return *this;
1364}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ operator=()

QwVQWK_Channel & QwVQWK_Channel::operator= ( const QwVQWK_Channel & value)

Definition at line 1126 of file QwVQWK_Channel.cc.

1127{
1128 if(this ==&value) return *this;
1129
1130 if (!IsNameEmpty()) {
1132 for (Int_t i=0; i<fBlocksPerEvent; i++){
1133 this->fBlock[i] = value.fBlock[i];
1134 this->fBlockM2[i] = value.fBlockM2[i];
1135 }
1139 this->fNumberOfSamples = value.fNumberOfSamples;
1140 this->fSequenceNumber = value.fSequenceNumber;
1141
1142 if (this->fDataToSave == kRaw){
1143 for (Int_t i=0; i<fBlocksPerEvent; i++){
1144 this->fBlock_raw[i] = value.fBlock_raw[i];
1145 }
1148 }
1149 }
1150 return *this;
1151}
VQwHardwareChannel & operator=(const VQwHardwareChannel &value)
Arithmetic assignment operator: Should only copy event-based data.

References fBlock, fBlock_raw, fBlockM2, fBlocksPerEvent, VQwHardwareChannel::fDataToSave, fHardwareBlockSum, fHardwareBlockSum_raw, fHardwareBlockSumError, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, fSoftwareBlockSum_raw, VQwDataElement::IsNameEmpty(), VQwDataElement::kRaw, VQwHardwareChannel::operator=(), and QwVQWK_Channel().

+ Here is the call graph for this function:

◆ PrintErrorCounterHead()

void QwVQWK_Channel::PrintErrorCounterHead ( )
static

Definition at line 1912 of file QwVQWK_Channel.cc.

1913{
1914 TString message;
1915 message = Form("%30s","Device name");
1916 message += Form("%9s", "HW Sat");
1917 message += Form("%9s", "Sample");
1918 message += Form("%9s", "SW_HW");
1919 message += Form("%9s", "Sequence");
1920 message += Form("%9s", "SameHW");
1921 message += Form("%9s", "ZeroHW");
1922 message += Form("%9s", "EventCut");
1923 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1924 QwMessage << message << QwLog::endl;
1925 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1926 return;
1927}

References QwLog::endl(), and QwMessage.

Referenced by QwBeamLine::PrintErrorCounters(), and QwBeamMod::PrintErrorCounters().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ PrintErrorCounters()

void QwVQWK_Channel::PrintErrorCounters ( ) const
overridevirtual

report number of events failed due to HW and event cut failure

Reimplemented from VQwDataElement.

Definition at line 1935 of file QwVQWK_Channel.cc.

1936{
1937 TString message;
1939 message = Form("%30s", GetElementName().Data());
1940 message += Form("%9d", fErrorCount_HWSat);
1941 message += Form("%9d", fErrorCount_sample);
1942 message += Form("%9d", fErrorCount_SW_HW);
1943 message += Form("%9d", fErrorCount_Sequence);
1944 message += Form("%9d", fErrorCount_SameHW);
1945 message += Form("%9d", fErrorCount_ZeroHW);
1946 message += Form("%9d", fNumEvtsWithEventCutsRejected);
1947
1948 if((fDataToSave == kRaw) && (!kFoundPedestal||!kFoundGain)){
1949 message += " >>>>> No Pedestal or Gain in map file";
1950 }
1951
1952 QwMessage << message << QwLog::endl;
1953 }
1954 return;
1955}

References QwLog::endl(), VQwHardwareChannel::fDataToSave, fErrorCount_HWSat, fErrorCount_SameHW, fErrorCount_sample, fErrorCount_Sequence, fErrorCount_SW_HW, fErrorCount_ZeroHW, fNumEvtsWithEventCutsRejected, VQwDataElement::GetElementName(), VQwHardwareChannel::kFoundGain, VQwHardwareChannel::kFoundPedestal, VQwDataElement::kRaw, and QwMessage.

+ Here is the call graph for this function:

◆ PrintErrorCounterTail()

void QwVQWK_Channel::PrintErrorCounterTail ( )
static

Definition at line 1929 of file QwVQWK_Channel.cc.

1930{
1931 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1932 return;
1933}

References QwLog::endl(), and QwMessage.

Referenced by QwBeamLine::PrintErrorCounters(), and QwBeamMod::PrintErrorCounters().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ PrintInfo()

void QwVQWK_Channel::PrintInfo ( ) const
overridevirtual

Print multiple lines of information about this data element.

Reimplemented from VQwDataElement.

Definition at line 607 of file QwVQWK_Channel.cc.

608{
609 std::cout<<"***************************************"<<"\n";
610 std::cout<<"Subsystem "<<GetSubsystemName()<<"\n"<<"\n";
611 std::cout<<"Beam Instrument Type: "<<GetModuleType()<<"\n"<<"\n";
612 std::cout<<"QwVQWK channel: "<<GetElementName()<<"\n"<<"\n";
613 std::cout<<"fPedestal= "<< fPedestal<<"\n";
614 std::cout<<"fCalibrationFactor= "<<fCalibrationFactor<<"\n";
615 std::cout<<"fBlocksPerEvent= "<<fBlocksPerEvent<<"\n"<<"\n";
616 std::cout<<"fSequenceNumber= "<<fSequenceNumber<<"\n";
617 std::cout<<"fNumberOfSamples= "<<fNumberOfSamples<<"\n";
618 std::cout<<"fBlock_raw ";
619
620 for (Int_t i = 0; i < fBlocksPerEvent; i++)
621 std::cout << " : " << fBlock_raw[i];
622 std::cout<<"\n";
623 std::cout<<"fHardwareBlockSum_raw= "<<fHardwareBlockSum_raw<<"\n";
624 std::cout<<"fSoftwareBlockSum_raw= "<<fSoftwareBlockSum_raw<<"\n";
625 std::cout<<"fBlock ";
626 for (Int_t i = 0; i < fBlocksPerEvent; i++)
627 std::cout << " : " <<std::setprecision(8) << fBlock[i];
628 std::cout << std::endl;
629
630 std::cout << "fHardwareBlockSum = "<<std::setprecision(8) <<fHardwareBlockSum << std::endl;
631 std::cout << "fHardwareBlockSumM2 = "<<fHardwareBlockSumM2 << std::endl;
632 std::cout << "fHardwareBlockSumError = "<<fHardwareBlockSumError << std::endl;
633
634 return;
635}

References fBlock, fBlock_raw, fBlocksPerEvent, VQwHardwareChannel::fCalibrationFactor, fHardwareBlockSum, fHardwareBlockSum_raw, fHardwareBlockSumError, fHardwareBlockSumM2, fNumberOfSamples, VQwHardwareChannel::fPedestal, fSequenceNumber, fSoftwareBlockSum_raw, VQwDataElement::GetElementName(), VQwDataElement::GetModuleType(), and VQwDataElement::GetSubsystemName().

+ Here is the call graph for this function:

◆ PrintValue()

void QwVQWK_Channel::PrintValue ( ) const
overridevirtual

Print single line of value and error of this data element.

Reimplemented from VQwDataElement.

Definition at line 1755 of file QwVQWK_Channel.cc.

1756{
1757 QwMessage << std::setprecision(8)
1758 << std::setw(18) << std::left << GetSubsystemName() << " "
1759 << std::setw(18) << std::left << GetModuleType() << " "
1760 << std::setw(18) << std::left << GetElementName() << " "
1761 << std::setw(12) << std::left << GetHardwareSum() << " +/- "
1762 << std::setw(12) << std::left << GetHardwareSumError() << " sig "
1763 << std::setw(12) << std::left << GetHardwareSumWidth() << " "
1764 << std::setw(10) << std::left << GetGoodEventCount() << " "
1765 << std::setw(12) << std::left << GetBlockValue(0) << " +/- "
1766 << std::setw(12) << std::left << GetBlockErrorValue(0) << " "
1767 << std::setw(12) << std::left << GetBlockValue(1) << " +/- "
1768 << std::setw(12) << std::left << GetBlockErrorValue(1) << " "
1769 << std::setw(12) << std::left << GetBlockValue(2) << " +/- "
1770 << std::setw(12) << std::left << GetBlockErrorValue(2) << " "
1771 << std::setw(12) << std::left << GetBlockValue(3) << " +/- "
1772 << std::setw(12) << std::left << GetBlockErrorValue(3) << " "
1773 << std::setw(12) << std::left << fGoodEventCount << " "
1774 << QwLog::endl;
1775 /*
1776 //for Debudding
1777 << std::setw(12) << std::left << fErrorFlag << " err "
1778 << std::setw(12) << std::left << fErrorConfigFlag << " c-err "
1779
1780 */
1781}
Double_t GetHardwareSumWidth() const
Double_t GetBlockErrorValue(size_t blocknum) const
UInt_t GetGoodEventCount() const

References QwLog::endl(), VQwDataElement::fGoodEventCount, GetBlockErrorValue(), GetBlockValue(), VQwDataElement::GetElementName(), VQwDataElement::GetGoodEventCount(), GetHardwareSum(), GetHardwareSumError(), GetHardwareSumWidth(), VQwDataElement::GetModuleType(), VQwDataElement::GetSubsystemName(), and QwMessage.

+ Here is the call graph for this function:

◆ ProcessEvBuffer()

Int_t QwVQWK_Channel::ProcessEvBuffer ( UInt_t * buffer,
UInt_t num_words_left,
UInt_t index = 0 )
overridevirtual

Decode the event data from a CODA buffer.

Process raw event buffer data for a VQWK ADC channel.

Parameters
bufferPointer to raw data buffer from DAQ system.
num_words_leftNumber of words remaining in the buffer.
indexChannel index within the ADC module (0-7).
Returns
Number of words consumed from the buffer.

This is a critical data processing function that decodes the 6-word VQWK ADC data format used throughout the Qweak/MOLLER analysis framework:

VQWK Data Format (6 words per channel):

  • Words 0-3: Individual block sums for 4 integration periods
  • Word 4: Hardware-calculated sum of all 4 blocks
  • Word 5: Combined sequence number (bits 8-15) and sample count (bits 16-31)

Data Processing Steps:

  1. Validates sufficient buffer space (6 words minimum)
  2. Copies raw data to local buffer with sign conversion (UInt_t -> Int_t)
  3. Extracts individual block sums and hardware sum
  4. Decodes sequence number for event ordering verification
  5. Extracts sample count for integration time normalization
  6. Calculates software block sum for hardware validation

Channel State Handling:

  • Empty channel names are skipped but consume buffer space
  • Insufficient buffer words trigger error messages
  • Raw data is stored for subsequent ProcessEvent() calibration

Error Detection:

  • Hardware vs software sum comparison (done in ProcessEvent())
  • Sequence number continuity checking
  • Sample count validation for proper integration

Buffer Management:

  • Always consumes exactly kWordsPerChannel (6) words when successful
  • Returns 0 on buffer underrun to indicate processing failure
  • Thread-safe local buffer prevents data corruption
Note
This function only processes raw data extraction. Calibration, pedestal subtraction, and physics calculations are performed in ProcessEvent().
Warning
Buffer underrun conditions will print error messages but may not halt processing, potentially causing downstream data corruption.

Implements VQwDataElement.

Definition at line 514 of file QwVQWK_Channel.cc.

515{
516 UInt_t words_read = 0;
517 UInt_t localbuf[kWordsPerChannel] = {0};
518 // The conversion from UInt_t to Double_t discards the sign, so we need an intermediate
519 // static_cast from UInt_t to Int_t.
520 Int_t localbuf_signed[kWordsPerChannel] = {0};
521
522 if (IsNameEmpty()){
523 // This channel is not used, but is present in the data stream.
524 // Skip over this data.
525 words_read = fNumberOfDataWords;
526 } else if (num_words_left >= fNumberOfDataWords)
527 {
528 for (Int_t i=0; i<kWordsPerChannel; i++){
529 localbuf[i] = buffer[i];
530 localbuf_signed[i] = static_cast<Int_t>(localbuf[i]);
531 }
532
534 for (Int_t i=0; i<fBlocksPerEvent; i++){
535 fBlock_raw[i] = localbuf_signed[i];
537 }
538 fHardwareBlockSum_raw = localbuf_signed[4];
539
540 /* Permanent change in the structure of the 6th word of the ADC readout.
541 * The upper 16 bits are the number of samples, and the upper 8 of the
542 * lower 16 are the sequence number. This matches the structure of
543 * the ADC readout in block read mode, and now also in register read mode.
544 * P.King, 2007sep04.
545 */
546 fSequenceNumber = (localbuf[5]>>8) & 0xFF;
547 fNumberOfSamples = (localbuf[5]>>16) & 0xFFFF;
548
549 words_read = fNumberOfDataWords;
550
551 } else
552 {
553 std::cerr << "QwVQWK_Channel::ProcessEvBuffer: Not enough words!"
554 << std::endl;
555 }
556 return words_read;
557}
UInt_t fNumberOfDataWords
Number of raw data words in this data element.

References fBlock_raw, fBlocksPerEvent, fHardwareBlockSum_raw, VQwHardwareChannel::fNumberOfDataWords, fNumberOfSamples, fSequenceNumber, fSoftwareBlockSum_raw, VQwDataElement::IsNameEmpty(), and kWordsPerChannel.

+ Here is the call graph for this function:

◆ ProcessEvent()

void QwVQWK_Channel::ProcessEvent ( )
overridevirtual

Process the event data according to pedestal and calibration factor.

Implements VQwHardwareChannel.

Definition at line 561 of file QwVQWK_Channel.cc.

562{
563 if (fNumberOfSamples == 0 && fHardwareBlockSum_raw == 0) {
564 // There isn't valid data for this channel. Just flag it and
565 // move on.
566 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
567 fBlock[i] = 0.0;
568 fBlockM2[i] = 0.0;
569 }
570 fHardwareBlockSum = 0.0;
573 } else if (fNumberOfSamples == 0) {
574 // This is probably a more serious problem.
575 QwWarning << "QwVQWK_Channel::ProcessEvent: Channel "
576 << this->GetElementName().Data()
577 << " has fNumberOfSamples==0 but has valid data in the hardware sum. "
578 << "Flag this as an error."
579 << QwLog::endl;
580 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
581 fBlock[i] = 0.0;
582 fBlockM2[i] = 0.0;
583 }
584 fHardwareBlockSum = 0.0;
587 } else {
588 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
590 fBlockM2[i] = 0.0; // second moment is zero for single events
591 }
593 fHardwareBlockSumM2 = 0.0; // second moment is zero for single events
594 }
595 return;
596}

References QwLog::endl(), fBlock, fBlock_raw, fBlockM2, fBlocksPerEvent, VQwHardwareChannel::fCalibrationFactor, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSum_raw, fHardwareBlockSumM2, fNumberOfSamples, VQwHardwareChannel::fPedestal, VQwDataElement::GetElementName(), kErrorFlag_sample, and QwWarning.

+ Here is the call graph for this function:

◆ Product()

void QwVQWK_Channel::Product ( const QwVQWK_Channel & value1,
const QwVQWK_Channel & value2 )

Definition at line 1468 of file QwVQWK_Channel.cc.

1469{
1470 if (!IsNameEmpty()){
1471 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
1472 this->fBlock[i] = (value1.fBlock[i]) * (value2.fBlock[i]);
1473 // For a single event the second moment is still zero
1474 this->fBlockM2[i] = 0.0;
1475 }
1476
1477 // For a single event the second moment is still zero
1478 this->fHardwareBlockSumM2 = 0.0;
1479
1481 this->fNumberOfSamples = value1.fNumberOfSamples;
1482 this->fSequenceNumber = 0;
1483 this->fErrorFlag = (value1.fErrorFlag|value2.fErrorFlag);
1484 }
1485 return;
1486}

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

Referenced by QwLinearDiodeArray::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ RandomizeEventData()

void QwVQWK_Channel::RandomizeEventData ( int helicity = 0.0,
double time = 0.0 )
overridevirtual

Internally generate random event data.

Implements MQwMockable.

Definition at line 270 of file QwVQWK_Channel.cc.

271{
272 //std::cout << "In channel: " << GetElementName() << std::endl;
273 // The blocks are assumed to be independent measurements
274 // Double_t* block = new Double_t[fBlocksPerEvent];
275 // Double_t sqrt_fBlocksPerEvent = 0.0;
276 // sqrt_fBlocksPerEvent = sqrt(fBlocksPerEvent);
277
278 // Calculate drift (if time is not specified, it stays constant at zero)
279 //Double_t drift = 0.0;
280 //for (UInt_t i = 0; i < fMockDriftFrequency.size(); i++) {
281 // drift += fMockDriftAmplitude[i] * sin(2.0 * Qw::pi * fMockDriftFrequency[i] * time + fMockDriftPhase[i]);
282 //}
283
284 // updated to calculate the drift for each block individually
285 Double_t drift = 0.0;
286 for (Int_t i = 0; i < fBlocksPerEvent; i++){
287 drift = 0.0;
288 if (i >= 1){
290 }
291 for (UInt_t i = 0; i < fMockDriftFrequency.size(); i++) {
292 drift += fMockDriftAmplitude[i] * sin(2.0 * Qw::pi * fMockDriftFrequency[i] * time + fMockDriftPhase[i]);
293 //std::cout << "Drift: " << drift << std::endl;
294 }
295 }
296
297 // Calculate signal
298 fHardwareBlockSum = 0.0;
299 fHardwareBlockSumM2 = 0.0; // second moment is zero for single events
300
301 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
302 double tmpvar = GetRandomValue();
303 //std::cout << "tmpvar: " << tmpvar << std::endl;
304 //std::cout << "->fMockSigma: " << fMockGaussianSigma << std::endl;
305
306 fBlock[i] = fMockGaussianMean + drift;
307 //std::cout << "(Start of loop) " << this->GetElementName() << "-> "<< "fBlock[" << i << "]: " << fBlock[i] << ", Drift: " << drift <<", Mean: " << fMockGaussianMean<< std::endl;
309 fBlock[i] += helicity*fMockAsymmetry;
310 } else {
311 fBlock[i] *= 1.0 + helicity*fMockAsymmetry;
312 }
313 fBlock[i] += fMockGaussianSigma*tmpvar*sqrt(fBlocksPerEvent);
314 //std::cout << "(End of loop) " << this->GetElementName() << "-> "<< "fBlock[" << i << "]: " << fBlock[i] << ", Drift: " << drift <<", Mean: " << fMockGaussianMean<< std::endl;
315
316
317/*
318 fBlock[i] = //GetRandomValue();
319 fMockGaussianMean * (1 + helicity * fMockAsymmetry)
320 + fMockGaussianSigma*sqrt(fBlocksPerEvent) * tmpvar
321 + drift; */
322
323 fBlockM2[i] = 0.0; // second moment is zero for single events
325
326/* std::cout << "In RandomizeEventData: " << tmpvar << " " << fMockGaussianMean << " "<< (1 + helicity * fMockAsymmetry) << " "
327 << fMockGaussianSigma << " " << fMockGaussianSigma*tmpvar << " "
328 << drift << " " << block[i] << std::endl; */
329 }
331 fSequenceNumber = 0;
333 // SetEventData(block);
334 // delete block;
335 return;
336}
static const double pi
Angles: base unit is radian.
Definition QwUnits.h:107
Double_t GetRandomValue()
bool fCalcMockDataAsDiff
static const Double_t kTimePerSample

References fBlock, fBlockM2, fBlocksPerEvent, MQwMockable::fCalcMockDataAsDiff, fHardwareBlockSum, fHardwareBlockSumM2, MQwMockable::fMockAsymmetry, MQwMockable::fMockDriftAmplitude, MQwMockable::fMockDriftFrequency, MQwMockable::fMockDriftPhase, MQwMockable::fMockGaussianMean, MQwMockable::fMockGaussianSigma, fNumberOfSamples, fNumberOfSamples_map, fSequenceNumber, MQwMockable::GetRandomValue(), kTimePerSample, and Qw::pi.

+ Here is the call graph for this function:

◆ Ratio()

void QwVQWK_Channel::Ratio ( const QwVQWK_Channel & numer,
const QwVQWK_Channel & denom )

Definition at line 1379 of file QwVQWK_Channel.cc.

1380{
1381 if (!IsNameEmpty()) {
1382 *this = numer;
1383 *this /= denom;
1384
1385 // Remaining variables
1387 fSequenceNumber = 0;
1389 fErrorFlag = (numer.fErrorFlag|denom.fErrorFlag);
1390 }
1391}

References VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fNumberOfSamples, fSequenceNumber, VQwDataElement::IsNameEmpty(), and QwVQWK_Channel().

Referenced by QwLinearDiodeArray::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ Scale()

void QwVQWK_Channel::Scale ( Double_t Offset)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1501 of file QwVQWK_Channel.cc.

1502{
1503 if (!IsNameEmpty()){
1504 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
1505 fBlock[i] *= scale;
1506 fBlockM2[i] *= scale * scale;
1507 }
1508 fHardwareBlockSum *= scale;
1509 fHardwareBlockSumM2 *= scale * scale;
1510 }
1511}

References fBlock, fBlockM2, fBlocksPerEvent, fHardwareBlockSum, fHardwareBlockSumM2, and VQwDataElement::IsNameEmpty().

Referenced by QwLinearDiodeArray::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ ScaledAdd()

void QwVQWK_Channel::ScaledAdd ( Double_t scale,
const VQwHardwareChannel * value )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1957 of file QwVQWK_Channel.cc.

1958{
1959 const QwVQWK_Channel* input = dynamic_cast<const QwVQWK_Channel*>(value);
1960
1961 // follows same steps as += but w/ scaling factor
1962 if(input!=NULL && !IsNameEmpty()){
1963 // QwWarning << "Adding " << input->GetElementName()
1964 // << " to " << GetElementName()
1965 // << " with scale factor " << scale
1966 // << QwLog::endl;
1967 // PrintValue();
1968 // input->PrintValue();
1969 for(Int_t i = 0; i < fBlocksPerEvent; i++){
1970 this -> fBlock[i] += scale * input->fBlock[i];
1971 this -> fBlockM2[i] = 0.0;
1972 }
1973 this -> fHardwareBlockSum += scale * input->fHardwareBlockSum;
1974 this -> fHardwareBlockSumM2 = 0.0;
1975 this -> fNumberOfSamples += input->fNumberOfSamples;
1976 this -> fSequenceNumber = 0;
1977 this -> fErrorFlag |= (input->fErrorFlag);
1978 } else if (input == NULL && value != NULL) {
1979 TString loc="Standard exception from QwVQWK_Channel::ScaledAdd "
1980 +value->GetElementName()+" "
1981 +this->GetElementName()+" are not of the same type";
1982 throw(std::invalid_argument(loc.Data()));
1983 }
1984}

References fBlock, fBlockM2, fBlocksPerEvent, VQwDataElement::fErrorFlag, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fSequenceNumber, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ SetCalibrationToVolts()

void QwVQWK_Channel::SetCalibrationToVolts ( )
inline

Definition at line 278 of file QwVQWK_Channel.h.

References kVQWK_VoltsPerBit, and VQwHardwareChannel::SetCalibrationFactor().

+ Here is the call graph for this function:

◆ SetDefaultSampleSize()

void QwVQWK_Channel::SetDefaultSampleSize ( size_t num_samples_map)
inline

Definition at line 125 of file QwVQWK_Channel.h.

125 {
126 // This will be checked against the no.of samples read by the module
127 fNumberOfSamples_map = num_samples_map;
128 };

References fNumberOfSamples_map.

Referenced by LoadChannelParameters().

+ Here is the caller graph for this function:

◆ SetEventData()

void QwVQWK_Channel::SetEventData ( Double_t * block,
UInt_t sequencenumber = 0 )

Definition at line 374 of file QwVQWK_Channel.cc.

375{
376 fHardwareBlockSum = 0.0;
377 fHardwareBlockSumM2 = 0.0; // second moment is zero for single events
378 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
379 fBlock[i] = block[i];
380 fBlockM2[i] = 0.0; // second moment is zero for single events
381 fHardwareBlockSum += block[i];
382 }
384
385 fSequenceNumber = sequencenumber;
387
388// Double_t thispedestal = 0.0;
389// thispedestal = fPedestal * fNumberOfSamples;
390
392 return;
393}
void SetRawEventData() override

References fBlock, fBlockM2, fBlocksPerEvent, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fNumberOfSamples_map, fSequenceNumber, and SetRawEventData().

Referenced by SetHardwareSum().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ SetHardwareSum()

void QwVQWK_Channel::SetHardwareSum ( Double_t hwsum,
UInt_t sequencenumber = 0 )

TODO: SetHardwareSum should be removed, and SetEventData should be made protected.

Definition at line 358 of file QwVQWK_Channel.cc.

359{
360 Double_t* block = new Double_t[fBlocksPerEvent];
361 for (Int_t i = 0; i < fBlocksPerEvent; i++){
362 block[i] = hwsum / fBlocksPerEvent;
363 }
364 SetEventData(block);
365 delete[] block;
366 return;
367}
void SetEventData(Double_t *block, UInt_t sequencenumber=0)

References fBlocksPerEvent, and SetEventData().

+ Here is the call graph for this function:

◆ SetRawEventData()

void QwVQWK_Channel::SetRawEventData ( )
overridevirtual

Implements MQwMockable.

Definition at line 395 of file QwVQWK_Channel.cc.

395 {
398// Double_t hwsum_test = 0.0;
399// std::cout << "*******In QwVQWK_Channel::SetRawEventData for channel:\t" << this->GetElementName() << std::endl;
400 for (Int_t i = 0; i < fBlocksPerEvent; i++)
401 {
402 // The raw data is decoded ino calibrated values with the following (in ProcessEvent()):
403 // fBlock[i] = fCalibrationFactor * ( (1.0 * fBlock_raw[i] * fBlocksPerEvent / fNumberOfSamples) - fPedestal );
404 // We should invert that function here:
405/* if (fBlock[i]<-10.0 || fBlock[i]>+10.0)
406 QwError << "In QwVQWK_Channel::SetRawEventData for channel:\t" << this->GetElementName() << ", Block " << i << " is out of range (-10 V,+10V):" << fBlock[i] << QwLog::endl;
407*/
410 //hwsum_test +=fBlock[i] /(fBlocksPerEvent * 1.0);
411
412
413 // fBlock[i] = fCalibrationFactor * ((1.0 * fBlock_raw[i] * fBlocksPerEvent / fNumberOfSamples) - fPedestal);
414 // fHardwareBlockSum += fBlock[i];
415
416 /* std::cout << "\t fBlock[i] = " << std::setprecision(6) << fBlock[i] << "\n"
417 << "\t fCalibrationFactor = " << fCalibrationFactor << "\n"
418 << "\t fPedestal = " << fPedestal << "\n"
419 << "\t fNumberOfSamples = " << fNumberOfSamples << "\n"
420 << "\t fBlocksPerEvent = " << fBlocksPerEvent << "\n"
421 << "\t fBlock[i] / fCalibrationFactor + fPedestal = " << fBlock[i] / fCalibrationFactor + fPedestal << "\n"
422 << "\t That * fNumberOfSamples / (fBlocksPerEvent * 1) = " << (fBlock[i] / fCalibrationFactor + fPedestal) * fNumberOfSamples / (fBlocksPerEvent * 1.0) << "\n"
423 << "\t fBlock_raw[i] = " << fBlock_raw[i] << "\n"
424 << "\t fHardwareBlockSum_raw = " << fHardwareBlockSum_raw << "\n"
425 << std::endl;
426 */
427 }
428
429/* std::cout << "fBlock[0] = " << std::setprecision(16) << fBlock[0] << std::endl
430 << "fBlockraw[0] after calib: " << fCalibrationFactor * ((1.0 * fBlock_raw[0] * fBlocksPerEvent / fNumberOfSamples) - fPedestal) << std::endl;
431
432 std::cout << "fHardwareBlockSum = " << std::setprecision(8) << fHardwareBlockSum << std::endl;
433 std::cout << "hwsum_test = " << std::setprecision(8) << hwsum_test << std::endl;
434 std::cout << "fHardwareBlockSum_raw = " << std::setprecision(8) << fHardwareBlockSum_raw << std::endl;
435 std::cout << "fHardwareBlockSum_Raw after calibration = " << fCalibrationFactor * ((1.0 * fHardwareBlockSum_raw / fNumberOfSamples) - fPedestal) << std::endl;
436*/
437
439
440 return;
441}

References fBlock, fBlock_raw, fBlocksPerEvent, VQwHardwareChannel::fCalibrationFactor, fHardwareBlockSum_raw, fNumberOfSamples, fNumberOfSamples_map, VQwHardwareChannel::fPedestal, and fSoftwareBlockSum_raw.

Referenced by SetEventData().

+ Here is the caller graph for this function:

◆ SetVQWKSaturationLimt()

void QwVQWK_Channel::SetVQWKSaturationLimt ( Double_t sat_volts = 8.5)
inline

Definition at line 223 of file QwVQWK_Channel.h.

223 {//Set the absolute staturation limit in volts.
224 fSaturationABSLimit=sat_volts;
225 }

References fSaturationABSLimit.

Referenced by QwVQWK_Channel(), and QwVQWK_Channel().

+ Here is the caller graph for this function:

◆ SmearByResolution()

void QwVQWK_Channel::SmearByResolution ( double resolution)
overridevirtual

Implements MQwMockable.

Definition at line 338 of file QwVQWK_Channel.cc.

338 {
339
340 fHardwareBlockSum = 0.0;
341 fHardwareBlockSumM2 = 0.0; // second moment is zero for single events
342 for (Int_t i = 0; i < fBlocksPerEvent; i++) {
343 // std::cout << i << " " << fBlock[i] << "->";
344 //std::cout << "resolution = " << resolution << "\t for channel \t" << GetElementName() << std::endl;
345 fBlock[i] += resolution*sqrt(fBlocksPerEvent) * GetRandomValue();
346 // std::cout << fBlock[i] << ": ";
347 fBlockM2[i] = 0.0; // second moment is zero for single events
349 }
350 // std::cout << std::endl;
352
354 // SetRawEventData();
355 return;
356}

References fBlock, fBlockM2, fBlocksPerEvent, fHardwareBlockSum, fHardwareBlockSumM2, fNumberOfSamples, fNumberOfSamples_map, and MQwMockable::GetRandomValue().

+ Here is the call graph for this function:

◆ SubtractValueFrom()

void QwVQWK_Channel::SubtractValueFrom ( const VQwHardwareChannel * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1197 of file QwVQWK_Channel.cc.

1198{
1199 const QwVQWK_Channel* tmpptr;
1200 tmpptr = dynamic_cast<const QwVQWK_Channel*>(valueptr);
1201 if (tmpptr!=NULL){
1202 *this -= *tmpptr;
1203 } else {
1204 TString loc="Standard exception from QwVQWK_Channel::SubtractValueFrom = "
1205 +valueptr->GetElementName()+" is an incompatible type.";
1206 throw std::invalid_argument(loc.Data());
1207 }
1208}

References VQwDataElement::GetElementName(), QwVQWK_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ Sum()

void QwVQWK_Channel::Sum ( const QwVQWK_Channel & value1,
const QwVQWK_Channel & value2 )

Definition at line 1367 of file QwVQWK_Channel.cc.

1368{
1369 *this = value1;
1370 *this += value2;
1371}

References QwVQWK_Channel().

Referenced by QwQPD::ProcessEvent().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Friends And Related Symbol Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream & stream,
const QwVQWK_Channel & channel )
friend

Definition at line 1783 of file QwVQWK_Channel.cc.

1784{
1785 stream << channel.GetHardwareSum();
1786 return stream;
1787}

References GetHardwareSum(), and QwVQWK_Channel().

Field Documentation

◆ bBlock

Bool_t QwVQWK_Channel::bBlock
private

Definition at line 397 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bBlock_raw

Bool_t QwVQWK_Channel::bBlock_raw
private

Definition at line 398 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bDEBUG

const Bool_t QwVQWK_Channel::bDEBUG =kFALSE
staticprivate

debugging display purposes

For VQWK data element trimming uses

Definition at line 392 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), FillTreeVector(), and MatchNumberOfSamples().

◆ bDevice_Error_Code

Bool_t QwVQWK_Channel::bDevice_Error_Code
private

Definition at line 400 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bHw_sum

Bool_t QwVQWK_Channel::bHw_sum
private

Definition at line 395 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bHw_sum_raw

Bool_t QwVQWK_Channel::bHw_sum_raw
private

Definition at line 396 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bNum_samples

Bool_t QwVQWK_Channel::bNum_samples
private

Definition at line 399 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ bSequence_number

Bool_t QwVQWK_Channel::bSequence_number
private

Definition at line 401 of file QwVQWK_Channel.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

◆ fADC_Same_NumEvt

Int_t QwVQWK_Channel::fADC_Same_NumEvt
private

Keep track of how many events with same ADC value returned.

Definition at line 382 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fBlock

◆ fBlock_raw

Int_t QwVQWK_Channel::fBlock_raw[4]
private

Array of the sub-block data as read from the module.

Definition at line 339 of file QwVQWK_Channel.h.

Referenced by ClearEventData(), EncodeEventData(), GetRawValue(), operator=(), PrintInfo(), ProcessEvBuffer(), ProcessEvent(), and SetRawEventData().

◆ fBlockError

Double_t QwVQWK_Channel::fBlockError[4]
private

Uncertainty on the sub-block.

Definition at line 354 of file QwVQWK_Channel.h.

Referenced by CalculateRunningAverage(), ClearEventData(), and GetValueError().

◆ fBlockM2

◆ fBlocksPerEvent

◆ fErrorCount_HWSat

Int_t QwVQWK_Channel::fErrorCount_HWSat
private

check to see ADC channel is saturated

Definition at line 368 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_SameHW

Int_t QwVQWK_Channel::fErrorCount_SameHW
private

check to see ADC returning same HW value

Definition at line 372 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_sample

Int_t QwVQWK_Channel::fErrorCount_sample
private

for sample size check

Definition at line 369 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_Sequence

Int_t QwVQWK_Channel::fErrorCount_Sequence
private

sequence number check

Definition at line 371 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_SW_HW

Int_t QwVQWK_Channel::fErrorCount_SW_HW
private

HW_sum==SW_sum check.

Definition at line 370 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_ZeroHW

Int_t QwVQWK_Channel::fErrorCount_ZeroHW
private

check to see ADC returning zero

Definition at line 373 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fHardwareBlockSum

◆ fHardwareBlockSum_raw

Int_t QwVQWK_Channel::fHardwareBlockSum_raw
private

Module-based sum of the four sub-blocks as read from the module.

Definition at line 340 of file QwVQWK_Channel.h.

Referenced by ClearEventData(), EncodeEventData(), GetRawValue(), operator=(), PrintInfo(), ProcessEvBuffer(), ProcessEvent(), and SetRawEventData().

◆ fHardwareBlockSumError

Double_t QwVQWK_Channel::fHardwareBlockSumError
private

Uncertainty on the hardware sum.

Definition at line 357 of file QwVQWK_Channel.h.

Referenced by AssignScaledValue(), CalculateRunningAverage(), ClearEventData(), GetValueError(), operator=(), and PrintInfo().

◆ fHardwareBlockSumM2

◆ fNumberOfSamples

◆ fNumberOfSamples_map

UInt_t QwVQWK_Channel::fNumberOfSamples_map
private

Number of samples in the expected to read through the module. This value is set in the QwBeamline map file.

Definition at line 364 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), CopyFrom(), ForceMapfileSampleSize(), InitializeChannel(), QwVQWK_Channel(), QwVQWK_Channel(), RandomizeEventData(), SetDefaultSampleSize(), SetEventData(), SetRawEventData(), and SmearByResolution().

◆ fNumEvtsWithEventCutsRejected

Int_t QwVQWK_Channel::fNumEvtsWithEventCutsRejected
private

Counts the Event cut rejected events.

Definition at line 375 of file QwVQWK_Channel.h.

Referenced by IncrementErrorCounters(), InitializeChannel(), and PrintErrorCounters().

◆ fPrev_HardwareBlockSum

Double_t QwVQWK_Channel::fPrev_HardwareBlockSum
private

Previous Module-based sum of the four sub-blocks.

Definition at line 385 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fPreviousSequenceNumber

UInt_t QwVQWK_Channel::fPreviousSequenceNumber
private

Previous event sequence number for this channel.

Definition at line 362 of file QwVQWK_Channel.h.

Referenced by InitializeChannel().

◆ fSaturationABSLimit

Double_t QwVQWK_Channel::fSaturationABSLimit
private

absolute value of the VQWK saturation volt

Definition at line 389 of file QwVQWK_Channel.h.

Referenced by CopyFrom(), GetVQWKSaturationLimt(), QwVQWK_Channel(), QwVQWK_Channel(), and SetVQWKSaturationLimt().

◆ fSequenceNo_Counter

Int_t QwVQWK_Channel::fSequenceNo_Counter
private

Internal counter to keep track of the sequence number.

Definition at line 384 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fSequenceNo_Prev

Int_t QwVQWK_Channel::fSequenceNo_Prev
private

Keep the sequence number of the last event.

Definition at line 383 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fSequenceNumber

◆ fSoftwareBlockSum_raw

Int_t QwVQWK_Channel::fSoftwareBlockSum_raw
private

Sum of the data in the four sub-blocks raw.

Definition at line 341 of file QwVQWK_Channel.h.

Referenced by ClearEventData(), GetRawSoftwareSum(), operator=(), PrintInfo(), ProcessEvBuffer(), and SetRawEventData().

◆ kDEBUG

const Bool_t QwVQWK_Channel::kDEBUG = kFALSE
staticprivate

Definition at line 317 of file QwVQWK_Channel.h.

Referenced by ConstructBranch(), and ConstructBranchAndVector().

◆ kMaxChannels

const Int_t QwVQWK_Channel::kMaxChannels = 8
staticprivate

Definition at line 319 of file QwVQWK_Channel.h.

Referenced by GetBufferOffset().

◆ kTimePerSample

const Double_t QwVQWK_Channel::kTimePerSample = (2.0/30.0) * Qw::us
static

Definition at line 62 of file QwVQWK_Channel.h.

Referenced by RandomizeEventData().

◆ kVQWK_VoltsPerBit

const Double_t QwVQWK_Channel::kVQWK_VoltsPerBit = (20./(1<<18))
staticprivate

Conversion factor to translate the average bit count in an ADC channel into average voltage. The base factor is roughly 76 uV per count, and zero counts corresponds to zero voltage. Store as the exact value for 20 V range, 18 bit ADC.

Definition at line 323 of file QwVQWK_Channel.h.

Referenced by ApplyHWChecks(), GetAverageVolts(), and SetCalibrationToVolts().

◆ kWordsPerChannel

const Int_t QwVQWK_Channel::kWordsPerChannel = 6
staticprivate

Definition at line 318 of file QwVQWK_Channel.h.

Referenced by GetBufferOffset(), and ProcessEvBuffer().


The documentation for this class was generated from the following files: