JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwADC18_Channel Class Reference

Concrete hardware channel for HAPPEX 18-bit ADC modules. More...

#include <QwADC18_Channel.h>

+ Inheritance diagram for QwADC18_Channel:
+ Collaboration diagram for QwADC18_Channel:

Public Member Functions

 QwADC18_Channel ()
 
 QwADC18_Channel (TString name, TString datatosave="raw")
 
 QwADC18_Channel (const QwADC18_Channel &value)
 
 QwADC18_Channel (const QwADC18_Channel &value, VQwDataElement::EDataToSave datatosave)
 
 ~QwADC18_Channel () override
 
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 SetEventData (Double_t value)
 
void SetRawEventData () override
 
void EncodeEventData (std::vector< UInt_t > &buffer) override
 Encode the event data into a CODA buffer.
 
Bool_t IsHeaderWord (UInt_t word) const
 Decode the event data from a CODA buffer.
 
Int_t ProcessDataWord (UInt_t word)
 
Int_t ProcessEvBuffer (UInt_t *buffer, UInt_t num_words_left, UInt_t index=0) override
 Process the CODA event buffer for this element.
 
void ProcessEvent () override
 Process the event data according to pedestal and calibration factor.
 
QwADC18_Channeloperator= (const QwADC18_Channel &value)
 
void AssignScaledValue (const QwADC18_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 QwADC18_Channel &value)
 
QwADC18_Channeloperator+= (const QwADC18_Channel &value)
 
QwADC18_Channeloperator-= (const QwADC18_Channel &value)
 
QwADC18_Channeloperator*= (const QwADC18_Channel &value)
 
VQwHardwareChanneloperator+= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator-= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator*= (const VQwHardwareChannel &input) override
 
VQwHardwareChanneloperator/= (const VQwHardwareChannel &input) override
 
const QwADC18_Channel operator+ (const QwADC18_Channel &value) const
 
const QwADC18_Channel operator- (const QwADC18_Channel &value) const
 
const QwADC18_Channel operator* (const QwADC18_Channel &value) const
 
void Sum (const QwADC18_Channel &value1, const QwADC18_Channel &value2)
 
void Difference (const QwADC18_Channel &value1, const QwADC18_Channel &value2)
 
void Ratio (const QwADC18_Channel &numer, const QwADC18_Channel &denom)
 
void Product (const QwADC18_Channel &value1, const QwADC18_Channel &value2)
 
void DivideBy (const QwADC18_Channel &denom)
 
void AddChannelOffset (Double_t Offset)
 
void Scale (Double_t Offset) override
 
void AccumulateRunningSum (const QwADC18_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 QwADC18_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 SetADC18SaturationLimt (Double_t sat_volts=8.5)
 
Double_t GetADC18SaturationLimt ()
 
Int_t ApplyHWChecks () override
 
void IncrementErrorCounters () 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
 
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
 
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 QwADC18_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)
 
virtual 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)
 
virtual 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 * Qw::us
 

Protected Member Functions

QwADC18_Channeloperator/= (const QwADC18_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 Attributes

UInt_t fDiff_Raw
 
UInt_t fBase_Raw
 
UInt_t fPeak_Raw
 
UInt_t fValue_Raw
 
Double_t fValue
 
Double_t fValueM2
 
Double_t fValueError
 
QwADC18_ChannelfRunningSum
 Pointer to the running sum for this channel.
 
QwADC18_ChannelfDAC
 Pointer to the DAC channel for this channel.
 

Static Private Attributes

static const Bool_t kDEBUG = kFALSE
 
static const Int_t kHeaderWordsPerBank = 1
 
static const Int_t kFooterWordsPerBank = 1
 
static const Int_t kHeaderWordsPerModule = 12
 
static const Int_t kFooterWordsPerModule = 1
 
static const Int_t kDataWordsPerChannel = 3
 
static const Int_t kMaxChannels = 4
 
static const UInt_t mask31x = 0x80000000
 
static const UInt_t mask3029x = 0x60000000
 
static const UInt_t mask2625x = 0x06000000
 
static const UInt_t mask2422x = 0x01c00000
 
static const UInt_t mask21x = 0x00200000
 
static const UInt_t mask200x = 0x001fffff
 
static const UInt_t mask2118x = 0x003c0000
 
static const UInt_t mask170x = 0x0003ffff
 
static const UInt_t mask150x = 0x0000ffff
 
ADC Calibration
static const Double_t kADC18_VoltsPerBit = (20./(1<<18))
 

Friends

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

Channel configuration data members

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 ADC18 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
 
static const Bool_t bDEBUG =kFALSE
 debugging display purposes
 

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 boost::mt19937 fRandomnessGenerator
 Internal randomness generator.
 
static boost::normal_distribution< double > fNormalDistribution
 Internal normal probability distribution.
 
static boost::variate_generator< boost::mt19937, boost::normal_distribution< double > > fNormalRandomVariable
 Internal normal random variable.
 

Detailed Description

Concrete hardware channel for HAPPEX 18-bit ADC modules.

Decodes and processes data from HAPPEX 18-bit ADC channels, providing access to raw and calibrated values, statistical moments, single-event cuts, and running statistics. Implements the dual-operator pattern for both type-specific and polymorphic operations.

Definition at line 34 of file QwADC18_Channel.h.

Constructor & Destructor Documentation

◆ QwADC18_Channel() [1/4]

QwADC18_Channel::QwADC18_Channel ( )
inline

Definition at line 61 of file QwADC18_Channel.h.

61 : MQwMockable() {
62 InitializeChannel("","");
63 SetADC18SaturationLimt(8.5);//FIXME set the default saturation limit
64 };
void SetADC18SaturationLimt(Double_t sat_volts=8.5)
void InitializeChannel(TString name, TString datatosave) override
Initialize the fields in this object.

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

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

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

◆ QwADC18_Channel() [2/4]

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

Definition at line 65 of file QwADC18_Channel.h.

65 : MQwMockable() {
66 InitializeChannel(name, datatosave);
67 SetADC18SaturationLimt(8.5);//FIXME set the default saturation limit
68 };

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

+ Here is the call graph for this function:

◆ QwADC18_Channel() [3/4]

QwADC18_Channel::QwADC18_Channel ( const QwADC18_Channel & value)
inline

Definition at line 69 of file QwADC18_Channel.h.

69 :
70 VQwHardwareChannel(value), MQwMockable(value),
73 {
74 *this = value;
75 };
Double_t fSaturationABSLimit
absolute value of the ADC18 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 fNumberOfSamples_map, fSaturationABSLimit, MQwMockable::MQwMockable(), QwADC18_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ QwADC18_Channel() [4/4]

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

Definition at line 76 of file QwADC18_Channel.h.

76 :
77 VQwHardwareChannel(value,datatosave), MQwMockable(value),
80 {
81 *this = value;
82 };

References fNumberOfSamples_map, fSaturationABSLimit, MQwMockable::MQwMockable(), QwADC18_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ ~QwADC18_Channel()

QwADC18_Channel::~QwADC18_Channel ( )
inlineoverride

Definition at line 83 of file QwADC18_Channel.h.

83{ };

Member Function Documentation

◆ AccumulateRunningSum() [1/2]

void QwADC18_Channel::AccumulateRunningSum ( const QwADC18_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 1108 of file QwADC18_Channel.cc.

1109{
1110
1111 if(count==0){
1112 count = value.fGoodEventCount;
1113 }
1114
1115 // Moment calculations
1116 Int_t n1 = fGoodEventCount;
1117 Int_t n2 = count;
1118
1119 // If there are no good events, check whether device HW is good
1120 if (n2 == 0 && value.fErrorFlag == 0) {
1121 n2 = 1;
1122 }
1123
1124 if (ErrorMask == kPreserveError){
1125 //n = 1;
1126 if (n2 == 0) {
1127 n2 = 1;
1128 }
1129 if (count == -1) {
1130 n2 = -1;
1131 }
1132 }
1133
1134 Int_t n = n1 + n2;
1135
1136 // Set up variables
1137 Double_t M11 = fValue;
1138 Double_t M12 = value.fValue;
1139 Double_t M22 = value.fValueM2;
1140 if (n2 == 0) {
1141 // no good events for addition
1142 return;
1143 } else if (n2 == 1) {
1144 // simple version for addition of single event
1146 fValue += (M12 - M11) / n;
1147 fValueM2 += (M12 - M11) * (M12 - fValue); // note: using updated mean
1148 } else if (n2 > 1) {
1149 // general version for addition of multi-event sets
1150 fGoodEventCount += n2;
1151 fValue += n2 * (M12 - M11) / n;
1152 fValueM2 += M22 + n1 * n2 * (M12 - M11) * (M12 - M11) / n;
1153 }
1154
1155 // Nanny
1156 if (fValue != fValue)
1157 QwWarning << "Angry Nanny: NaN detected in " << GetElementName() << QwLog::endl;
1158}
#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
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(), VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fValue, fValueM2, VQwDataElement::GetElementName(), kPreserveError, QwADC18_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 QwADC18_Channel::AccumulateRunningSum ( const VQwHardwareChannel * value,
Int_t count = 0,
Int_t ErrorMask = 0xFFFFFFF )
inlineoverridevirtual

Reimplemented from VQwHardwareChannel.

Definition at line 168 of file QwADC18_Channel.h.

168 {
169 const QwADC18_Channel *tmp_ptr = dynamic_cast<const QwADC18_Channel*>(value);
170 if (tmp_ptr != NULL) {
171 AccumulateRunningSum(*tmp_ptr, count, ErrorMask);
172 } else {
173 throw std::invalid_argument("Standard exception from QwADC18_Channel::AccumulateRunningSum: incompatible hardware channel type");
174 }
175 };
void AccumulateRunningSum(const QwADC18_Channel &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)

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

+ Here is the call graph for this function:

◆ AddChannelOffset()

void QwADC18_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 1039 of file QwADC18_Channel.cc.

1040{
1041 if (!IsNameEmpty())
1042 {
1043 fValue += offset;
1044 }
1045}
Bool_t IsNameEmpty() const
Is the name of this element empty?

References fValue, and VQwDataElement::IsNameEmpty().

+ Here is the call graph for this function:

◆ AddValueFrom()

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

Implements VQwHardwareChannel.

Definition at line 777 of file QwADC18_Channel.cc.

778{
779 const QwADC18_Channel* tmpptr;
780 tmpptr = dynamic_cast<const QwADC18_Channel*>(valueptr);
781 if (tmpptr!=NULL){
782 *this += *tmpptr;
783 } else {
784 TString loc="Standard exception from QwADC18_Channel::AddValueFrom = "
785 +valueptr->GetElementName()+" is an incompatible type.";
786 throw std::invalid_argument(loc.Data());
787 }
788}

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

+ Here is the call graph for this function:

◆ ApplyHWChecks()

Int_t QwADC18_Channel::ApplyHWChecks ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 93 of file QwADC18_Channel.cc.

94{
95 Bool_t bStatus;
96 if (bEVENTCUTMODE>0){//Global switch to ON/OFF event cuts set at the event cut file
97
98 if (bDEBUG)
99 QwWarning<<" QwQWVK_Channel "<<GetElementName()<<" "<<GetNumberOfSamples()<<QwLog::endl;
100
101 // Sample size check
102 bStatus = MatchNumberOfSamples(fNumberOfSamples_map);//compare the default sample size with no.of samples read by the module
103 if (!bStatus) {
105 }
106
107 //check sequence number
109 if (fSequenceNo_Counter==0 || GetSequenceNumber()==0){//starting the data run
111 }
112
113 if (!MatchSequenceNumber(fSequenceNo_Prev)){//we have a sequence number error
115 if (bDEBUG) QwWarning<<" QwQWVK_Channel "<<GetElementName()<<" Sequence number previous value = "<<fSequenceNo_Prev<<" Current value= "<< GetSequenceNumber()<<QwLog::endl;
116 }
117 }
118 else {
119 fGoodEventCount = 1;
120 fErrorFlag = 0;
121 }
122
123 return fErrorFlag;
124}
static const UInt_t kErrorFlag_sample
Definition QwTypes.h:160
static const UInt_t kErrorFlag_Sequence
Definition QwTypes.h:162
Int_t fSequenceNo_Counter
Internal counter to keep track of the sequence number.
Int_t fSequenceNo_Prev
Keep the sequence number of the last event.
size_t GetSequenceNumber() const
Bool_t MatchNumberOfSamples(size_t numsamp)
static const Bool_t bDEBUG
debugging display purposes
Bool_t MatchSequenceNumber(size_t seqnum)
size_t GetNumberOfSamples() const

References bDEBUG, VQwHardwareChannel::bEVENTCUTMODE, QwLog::endl(), VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fNumberOfSamples_map, fSequenceNo_Counter, fSequenceNo_Prev, VQwDataElement::GetElementName(), GetNumberOfSamples(), GetSequenceNumber(), kErrorFlag_sample, kErrorFlag_Sequence, MatchNumberOfSamples(), MatchSequenceNumber(), and QwWarning.

+ Here is the call graph for this function:

◆ ApplySingleEventCuts() [1/2]

Bool_t QwADC18_Channel::ApplySingleEventCuts ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1271 of file QwADC18_Channel.cc.

1272{
1273 Bool_t status;
1274
1275 if (bEVENTCUTMODE>=2){//Global switch to ON/OFF event cuts set at the event cut file
1276
1277 if (fULimit < fLLimit){
1278 status=kTRUE;
1279 } else if (GetValue()<=fULimit && GetValue()>=fLLimit){
1280 if ((fErrorFlag)==0)
1281 status=kTRUE;
1282 else
1283 status=kFALSE;//If the device HW is failed
1284 }
1285 else{
1286 if (GetValue()> fULimit)
1288 else
1290 status=kFALSE;
1291 }
1292
1293 if (bEVENTCUTMODE==3){
1294 status=kTRUE; //Update the event cut fail flag but pass the event.
1295 }
1296
1297
1298 }
1299 else{
1300 status=kTRUE;
1301 //fErrorFlag=0;//we need to keep the device error codes
1302 }
1303
1304 return status;
1305}
static const UInt_t kErrorFlag_EventCut_L
Definition QwTypes.h:165
static const UInt_t kErrorFlag_EventCut_U
Definition QwTypes.h:166
Double_t GetValue() const

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

+ Here is the call graph for this function:

◆ ApplySingleEventCuts() [2/2]

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

Definition at line 1255 of file QwADC18_Channel.cc.

1256{
1257 Bool_t status = kFALSE;
1258
1259 if (UL < LL){
1260 status=kTRUE;
1261 } else if (GetValue()<=UL && GetValue()>=LL){
1262 if ((fErrorFlag & kPreserveError)!=0)
1263 status=kTRUE;
1264 else
1265 status=kFALSE;//If the device HW is failed
1266 }
1267 std::cout<<(this->fErrorFlag & kPreserveError)<<std::endl;
1268 return status;
1269}

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

+ Here is the call graph for this function:

◆ ArcTan()

void QwADC18_Channel::ArcTan ( const QwADC18_Channel & value)

Definition at line 1013 of file QwADC18_Channel.cc.

1014{
1015 if (!IsNameEmpty()) {
1016 this->fValue = atan(value.fValue);
1017 }
1018 // FIXME stuff missing here
1019}

References fValue, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ AssignScaledValue()

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

Definition at line 751 of file QwADC18_Channel.cc.

753{
754 if(this == &value) return;
755
756 if (!IsNameEmpty()) {
757 this->fValue = value.fValue * scale;
758 this->fValueError = value.fValueError;
759 this->fValueM2 = value.fValueM2 * scale * scale;
760 this->fErrorFlag = value.fErrorFlag;//error code is updated.
761 this->fGoodEventCount = value.fGoodEventCount;
762 }
763}

References VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fValue, fValueError, fValueM2, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ AssignValueFrom()

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

Implements VQwHardwareChannel.

Definition at line 765 of file QwADC18_Channel.cc.

766{
767 const QwADC18_Channel* tmpptr;
768 tmpptr = dynamic_cast<const QwADC18_Channel*>(valueptr);
769 if (tmpptr!=NULL){
770 *this = *tmpptr;
771 } else {
772 TString loc="Standard exception from QwADC18_Channel::AssignValueFrom = "
773 +valueptr->GetElementName()+" is an incompatible type.";
774 throw std::invalid_argument(loc.Data());
775 }
776}

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

+ Here is the call graph for this function:

◆ Blind() [1/2]

void QwADC18_Channel::Blind ( const QwBlinder * blinder)

Blind this channel as an asymmetry.

Blind this channel as an asymmetry

Parameters
blinderBlinder

Definition at line 1199 of file QwADC18_Channel.cc.

1200{
1201 if (!IsNameEmpty()) {
1202 if (blinder->IsBlinderOkay() && ((fErrorFlag)==0) ){
1203 blinder->BlindValue(fValue);
1204 } else {
1207 }
1208 }
1209}
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(), VQwDataElement::fErrorFlag, fValue, QwBlinder::IsBlinderOkay(), VQwDataElement::IsNameEmpty(), QwBlinder::kValue_BlinderFail, and QwBlinder::ModifyThisErrorCode().

+ Here is the call graph for this function:

◆ Blind() [2/2]

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

Blind this channel as a difference.

Blind this channel as a difference with specified yield

Parameters
blinderBlinder
yieldCorresponding yield

Definition at line 1216 of file QwADC18_Channel.cc.

1217{
1218 if (!IsNameEmpty()) {
1219 if (blinder->IsBlinderOkay() && ((fErrorFlag) ==0) ){
1220 blinder->BlindValue(fValue, yield.fValue);
1221 } else {
1222 blinder->ModifyThisErrorCode(fErrorFlag);//update the HW error code
1224 }
1225 }
1226}

References QwBlinder::BlindValue(), VQwDataElement::fErrorFlag, fValue, QwBlinder::IsBlinderOkay(), VQwDataElement::IsNameEmpty(), QwBlinder::kValue_BlinderFail, QwBlinder::ModifyThisErrorCode(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ CalculateRunningAverage()

void QwADC18_Channel::CalculateRunningAverage ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1161 of file QwADC18_Channel.cc.

1162{
1163 // See notes in QwVQWK_Channel; we are using:
1164 // error = sqrt(M2)/n,
1165 // or alternately we could use the unbiased estimator for both
1166 // the sigma and error on the mean:
1167 // error = sqrt(M2)/(n-1)
1168 if(fGoodEventCount > 0) {
1170 } else {
1171 fValueError = 0.0;
1172 }
1173}

References VQwDataElement::fGoodEventCount, fValueError, and fValueM2.

◆ ClearEventData()

void QwADC18_Channel::ClearEventData ( )
overridevirtual

Clear the event data in this element.

Reimplemented from VQwDataElement.

Definition at line 234 of file QwADC18_Channel.cc.

235{
236 fDiff_Raw = 0;
237 fBase_Raw = 0;
238 fPeak_Raw = 0;
239 fValue = 0.0;
240 fValueM2 = 0.0;
241 fValueError = 0.0;
242 fSequenceNumber = 0;
244 fGoodEventCount = 0;
245 fErrorFlag = 0;
246}
UInt_t fNumberOfSamples
Number of samples read through the module.
UInt_t fSequenceNumber
Event sequence number for this channel.

References fBase_Raw, fDiff_Raw, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fNumberOfSamples, fPeak_Raw, fSequenceNumber, fValue, fValueError, and fValueM2.

Referenced by InitializeChannel().

+ 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 * QwADC18_Channel::Clone ( VQwDataElement::EDataToSave datatosave) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 87 of file QwADC18_Channel.h.

87 {
88 return new QwADC18_Channel(*this,datatosave);
89 };

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

+ Here is the call graph for this function:

◆ ConstructBranch()

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

Implements VQwHardwareChannel.

Definition at line 603 of file QwADC18_Channel.cc.

604{
605 if (IsNameEmpty()){
606 // This channel is not used, so skip setting up the tree.
607 } else {
608 TString basename = prefix + GetElementName();
609 tree->Branch(basename, &fValue, basename+"/D");
610 }
611}

References fValue, VQwDataElement::GetElementName(), and VQwDataElement::IsNameEmpty().

+ Here is the call graph for this function:

◆ ConstructBranchAndVector()

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

Implements VQwHardwareChannel.

Definition at line 572 of file QwADC18_Channel.cc.

573{
574 if (IsNameEmpty()){
575 // This channel is not used, so skip setting up the tree.
576 } else {
577 // Decide what to store based on prefix
578 SetDataToSaveByPrefix(prefix);
579
580 TString basename = prefix(0, (prefix.First("|") >= 0)? prefix.First("|"): prefix.Length()) + GetElementName();
581 fTreeArrayIndex = values.size();
582
583 values.push_back("value", 'D');
584 if (fDataToSave == kMoments) {
585 values.push_back("value_m2", 'D');
586 values.push_back("value_err", 'D');
587 }
588
589 values.push_back("Device_Error_Code", 'i');
590 if (fDataToSave == kRaw){
591 values.push_back("raw", 'I');
592 values.push_back("diff", 'I');
593 values.push_back("peak", 'I');
594 values.push_back("base", 'I');
595 }
596
598 if (gQwHists.MatchDeviceParamsFromList(basename.Data()))
599 tree->Branch(basename, &(values[fTreeArrayIndex]), values.LeafList(fTreeArrayIndex).c_str());
600 }
601}
QwHistogramHelper gQwHists
Globally defined instance of the QwHistogramHelper class.
Bool_t MatchDeviceParamsFromList(const std::string &devicename)
size_type size() const noexcept
Definition QwRootFile.h:81
std::string LeafList(size_type start_index=0) const
Definition QwRootFile.h:228
void push_back(const std::string &name, const char type='D')
Definition QwRootFile.h:195
void SetDataToSaveByPrefix(const TString &prefix)
Set the flag indicating if raw or derived values are in this data element based on prefix.

References VQwHardwareChannel::fDataToSave, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, VQwDataElement::GetElementName(), gQwHists, VQwDataElement::IsNameEmpty(), VQwDataElement::kMoments, VQwDataElement::kRaw, QwRootTreeBranchVector::LeafList(), QwRootTreeBranchVector::push_back(), VQwHardwareChannel::SetDataToSaveByPrefix(), and QwRootTreeBranchVector::size().

+ Here is the call graph for this function:

◆ ConstructHistograms()

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

Construct the histograms for this data element.

Implements VQwDataElement.

Definition at line 511 of file QwADC18_Channel.cc.

512{
513 // If we have defined a subdirectory in the ROOT file, then change into it.
514 if (folder != NULL) folder->cd();
515
516 if (IsNameEmpty()){
517 // This channel is not used, so skip filling the histograms.
518 } else {
519 // Now create the histograms.
520 if (prefix == TString("asym_")
521 || prefix == TString("diff_")
522 || prefix == TString("yield_"))
524
525 TString basename, fullname;
526 basename = prefix + GetElementName();
527
528 if(fDataToSave==kRaw)
529 {
530 fHistograms.resize(2, NULL);
531 size_t index=0;
532 fHistograms[index++] = gQwHists.Construct1DHist(basename);
533 fHistograms[index++] = gQwHists.Construct1DHist(basename+Form("_raw"));
534 }
535 else if(fDataToSave==kDerived)
536 {
537 fHistograms.resize(1, NULL);
538 Int_t index=0;
539 fHistograms[index++] = gQwHists.Construct1DHist(basename);
540 }
541 else
542 {
543 // this is not recognized
544 }
545 }
546}
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
TH1F * Construct1DHist(const TString &inputfile, const TString &name_title)

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

+ Here is the call graph for this function:

◆ DeaccumulateRunningSum() [1/2]

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

Definition at line 177 of file QwADC18_Channel.h.

177 {
178 AccumulateRunningSum(value, -1, ErrorMask);
179 };

References AccumulateRunningSum(), and QwADC18_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 255 of file VQwHardwareChannel.h.

255 {
256 AccumulateRunningSum(value, -1, ErrorMask);
257 };

◆ Difference()

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

Definition at line 948 of file QwADC18_Channel.cc.

949{
950 *this = value1;
951 *this -= value2;
952}

References QwADC18_Channel().

+ Here is the call graph for this function:

◆ DivideBy() [1/2]

void QwADC18_Channel::DivideBy ( const QwADC18_Channel & denom)

Definition at line 1056 of file QwADC18_Channel.cc.

1057{
1058 *this /= denom;
1059}

References QwADC18_Channel().

+ Here is the call graph for this function:

◆ DivideBy() [2/2]

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

Implements VQwHardwareChannel.

Definition at line 813 of file QwADC18_Channel.cc.

814{
815 const QwADC18_Channel* tmpptr;
816 tmpptr = dynamic_cast<const QwADC18_Channel*>(valueptr);
817 if (tmpptr!=NULL){
818 *this /= *tmpptr;
819 } else {
820 TString loc="Standard exception from QwADC18_Channel::DivideBy = "
821 +valueptr->GetElementName()+" is an incompatible type.";
822 throw std::invalid_argument(loc.Data());
823 }
824}

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

+ Here is the call graph for this function:

◆ EncodeEventData()

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

Encode the event data into a CODA buffer.

Implements MQwMockable.

Definition at line 295 of file QwADC18_Channel.cc.

296{
297 UInt_t localbuf[kDataWordsPerChannel] = {0};
298
299 if (IsNameEmpty()) {
300 // This channel is not used, but is present in the data stream.
301 // Fill in with zero.
302 localbuf[0] = 0;
303 } else {
304 localbuf[0] |= mask31x & (0 << 31);
305 localbuf[0] |= mask3029x & (0 << 29);
306 localbuf[0] |= mask2625x & (0 << 25);
307 UInt_t exp_dtype = 1;
308 switch (exp_dtype) {
309 case 0:
310 break;
311 case 1:
312 case 2:
313 localbuf[0] |= mask2422x & (2 << 22);
314 localbuf[0] |= mask170x & fPeak_Raw;
315 break;
316 case 4:
317 break;
318 default:
319 QwError << "QwADC18_Channel::EncodeEventData: Unknown data type" << QwLog::endl;
320 }
321 }
322
323 for (Int_t i = 0; i < kDataWordsPerChannel; i++) {
324 buffer.push_back(localbuf[i]);
325 }
326}
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
static const UInt_t mask3029x
static const UInt_t mask2422x
static const UInt_t mask2625x
static const UInt_t mask170x
static const UInt_t mask31x
static const Int_t kDataWordsPerChannel

References QwLog::endl(), fPeak_Raw, VQwDataElement::IsNameEmpty(), kDataWordsPerChannel, mask170x, mask2422x, mask2625x, mask3029x, mask31x, and QwError.

+ Here is the call graph for this function:

◆ FillHistograms()

void QwADC18_Channel::FillHistograms ( )
overridevirtual

Fill the histograms for this data element.

Implements VQwDataElement.

Definition at line 548 of file QwADC18_Channel.cc.

549{
550 Int_t index=0;
551
552 if (IsNameEmpty())
553 {
554 // This channel is not used, so skip creating the histograms.
555 } else
556 {
557 if(fDataToSave==kRaw)
558 {
559 if (fHistograms[index] != NULL && (fErrorFlag)==0)
560 fHistograms[index++]->Fill(GetValue());
561 if (fHistograms[index] != NULL && (fErrorFlag)==0)
562 fHistograms[index++]->Fill(GetRawValue());
563 }
564 else if(fDataToSave==kDerived)
565 {
566 if (fHistograms[index] != NULL && (fErrorFlag)==0)
567 fHistograms[index++]->Fill(GetValue());
568 }
569 }
570}
Int_t GetRawValue() const

References VQwHardwareChannel::fDataToSave, VQwDataElement::fErrorFlag, MQwHistograms::fHistograms, GetRawValue(), GetValue(), VQwDataElement::IsNameEmpty(), VQwDataElement::kDerived, and VQwDataElement::kRaw.

+ Here is the call graph for this function:

◆ FillTreeVector()

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

Implements VQwHardwareChannel.

Definition at line 613 of file QwADC18_Channel.cc.

614{
615 if (IsNameEmpty()) {
616 // This channel is not used, so skip setting up the tree.
617 } else if (fTreeArrayNumEntries == 0) {
618 static bool warned = false;
619 if (!warned) {
620 QwError << "QwADC18_Channel::FillTreeVector: fTreeArrayNumEntries=="
621 << fTreeArrayNumEntries << " (no branch constructed?)" << QwLog::endl;
622 QwError << "Offending element is " << GetElementName() << QwLog::endl;
623 warned = true;
624 }
625 } else if (values.size() < fTreeArrayIndex+fTreeArrayNumEntries) {
626 QwError << "QwADC18_Channel::FillTreeVector: values.size()=="
627 << values.size() << " name: " << fElementName
628 << "; fTreeArrayIndex+fTreeArrayNumEntries=="
629 << fTreeArrayIndex << '+' << fTreeArrayNumEntries << '='
631 << QwLog::endl;
632 } else {
633 size_t index = fTreeArrayIndex;
634 values.SetValue(index++, this->fValue);
635 if (fDataToSave == kMoments) {
636 values.SetValue(index++, fValueM2);
637 values.SetValue(index++, fValueError);
638 }
639 values.SetValue(index++, this->fErrorFlag);
640 if(fDataToSave==kRaw){
641 values.SetValue(index++, this->fValue_Raw);
642 values.SetValue(index++, this->fDiff_Raw);
643 values.SetValue(index++, this->fPeak_Raw);
644 values.SetValue(index++, this->fBase_Raw);
645 }
646 }
647}
void SetValue(size_type index, Double_t val)
Definition QwRootFile.h:108
TString fElementName
Name of this data element.

References QwLog::endl(), fBase_Raw, VQwHardwareChannel::fDataToSave, fDiff_Raw, VQwDataElement::fElementName, VQwDataElement::fErrorFlag, fPeak_Raw, VQwHardwareChannel::fTreeArrayIndex, VQwHardwareChannel::fTreeArrayNumEntries, fValue, fValue_Raw, fValueError, fValueM2, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), VQwDataElement::kMoments, VQwDataElement::kRaw, QwError, QwRootTreeBranchVector::SetValue(), and QwRootTreeBranchVector::size().

+ Here is the call graph for this function:

◆ ForceMapfileSampleSize()

void QwADC18_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 112 of file QwADC18_Channel.h.

References fNumberOfSamples, and fNumberOfSamples_map.

◆ GetADC18SaturationLimt()

Double_t QwADC18_Channel::GetADC18SaturationLimt ( )
inline

Definition at line 210 of file QwADC18_Channel.h.

210 {
211 return fSaturationABSLimit;
212 }

References fSaturationABSLimit.

◆ GetAverageVolts()

Double_t QwADC18_Channel::GetAverageVolts ( ) const

Definition at line 490 of file QwADC18_Channel.cc.

491{
493}
static const Double_t kADC18_VoltsPerBit

References fNumberOfSamples, fValue, and kADC18_VoltsPerBit.

◆ GetBufferOffset()

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

Class: QwADC18_Channel Base class containing decoding functions for the HAPPEX 18-bit ADC The functions in this class will decode a single channel worth of ADC18_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 -1 (header)
Returns
The number of words offset to the beginning of this channel's data from the beginning of the ADC18 buffer.

Definition at line 62 of file QwADC18_Channel.cc.

63{
64 Int_t offset = -1;
65 if (moduleindex<0 ){
66 QwError << "QwADC18_Channel::GetBufferOffset: Invalid module index,"
67 << moduleindex
68 << ". Must be zero or greater."
69 << QwLog::endl;
70 } else if (channelindex<-1 || channelindex>kMaxChannels){
71 QwError << "QwADC18_Channel::GetBufferOffset: Invalid channel index,"
72 << channelindex
73 << ". Must be in range [-1," << kMaxChannels << "]."
74 << QwLog::endl;
75 } else {
76 if (channelindex == -1) {
77 // Header
78 offset = kHeaderWordsPerBank
80 } else {
81 // Data word
82 offset = kHeaderWordsPerBank
85 + channelindex * kDataWordsPerChannel;
86 }
87 }
88 return offset;
89}
static const Int_t kFooterWordsPerModule
static const Int_t kMaxChannels
static const Int_t kHeaderWordsPerBank
static const Int_t kHeaderWordsPerModule

References QwLog::endl(), kDataWordsPerChannel, kFooterWordsPerModule, kHeaderWordsPerBank, kHeaderWordsPerModule, kMaxChannels, and QwError.

Referenced by QwBeamDetectorID::QwBeamDetectorID().

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

◆ GetNumberOfSamples()

size_t QwADC18_Channel::GetNumberOfSamples ( ) const
inline

Definition at line 242 of file QwADC18_Channel.h.

242{return (fNumberOfSamples);};

References fNumberOfSamples.

Referenced by ApplyHWChecks().

+ 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 FillHistograms().

+ Here is the caller graph for this function:

◆ GetRawValue() [2/2]

Int_t QwADC18_Channel::GetRawValue ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 223 of file QwADC18_Channel.h.

223{ return fValue_Raw; };

References fValue_Raw.

◆ GetSequenceNumber()

size_t QwADC18_Channel::GetSequenceNumber ( ) const
inline

Definition at line 241 of file QwADC18_Channel.h.

241{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 ApplySingleEventCuts(), ApplySingleEventCuts(), FillHistograms(), and PrintValue().

+ Here is the caller graph for this function:

◆ GetValue() [2/2]

Double_t QwADC18_Channel::GetValue ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 224 of file QwADC18_Channel.h.

224{ return fValue; };

References fValue.

Referenced by operator<<.

+ 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);};
Double_t GetValueError() const

Referenced by PrintValue().

+ Here is the caller graph for this function:

◆ GetValueError() [2/2]

Double_t QwADC18_Channel::GetValueError ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 226 of file QwADC18_Channel.h.

226{ return fValueError; };

References fValueError.

◆ GetValueM2() [1/2]

Double_t VQwHardwareChannel::GetValueM2 ( ) const
inline

Definition at line 127 of file VQwHardwareChannel.h.

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

◆ GetValueM2() [2/2]

Double_t QwADC18_Channel::GetValueM2 ( size_t element) const
inlineoverridevirtual

Implements VQwHardwareChannel.

Definition at line 225 of file QwADC18_Channel.h.

225{ return fValueM2; };

References fValueM2.

◆ GetValueWidth() [1/2]

Double_t VQwHardwareChannel::GetValueWidth ( ) const
inline

Definition at line 129 of file VQwHardwareChannel.h.

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

Referenced by PrintValue().

+ 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 };
void RangeCheck(size_t element) const
Checks that the requested element is in range, to be used in accesses to subelements similar to std::...

◆ IncrementErrorCounters()

void QwADC18_Channel::IncrementErrorCounters ( )
overridevirtual

Implements VQwHardwareChannel.

Definition at line 128 of file QwADC18_Channel.cc.

128 {
130 fErrorCount_sample++; //increment the hw error counter
132 fErrorCount_SW_HW++; //increment the hw error counter
134 fErrorCount_Sequence++; //increment the hw error counter
136 fErrorCount_SameHW++; //increment the hw error counter
138 fErrorCount_ZeroHW++; //increment the hw error counter
141 fNumEvtsWithEventCutsRejected++; //increment the event cut error counter
142 }
143}
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_SameHW
Definition QwTypes.h:163
Int_t fErrorCount_SW_HW
HW_sum==SW_sum check.
Int_t fErrorCount_SameHW
check to see ADC returning same HW value
Int_t fErrorCount_sample
for sample size check
Int_t fErrorCount_ZeroHW
check to see ADC returning zero
Int_t fNumEvtsWithEventCutsRejected
Counts the Event cut rejected events.
Int_t fErrorCount_Sequence
sequence number check

References 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, and kErrorFlag_ZeroHW.

◆ InitializeChannel() [1/2]

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

Initialize the fields in this object.

Implements VQwHardwareChannel.

Definition at line 147 of file QwADC18_Channel.cc.

148{
149 SetElementName(name);
150 SetDataToSave(datatosave);
153
154 kFoundPedestal = 0;
155 kFoundGain = 0;
156
157 fPedestal = 0.0;
158 fCalibrationFactor = 1.0;
159
160 fTreeArrayIndex = 0;
162
164
168
169 // Use internal random variable by default
171
172 // Mock drifts
173 fMockDriftAmplitude.clear();
174 fMockDriftFrequency.clear();
175 fMockDriftPhase.clear();
176
177 // Mock asymmetries
178 fMockAsymmetry = 0.0;
179 fMockGaussianMean = 0.0;
180 fMockGaussianSigma = 0.0;
181
182 // Event cuts
183 fULimit=-1;
184 fLLimit=1;
186
187 fErrorFlag=0; //Initialize the error flag
188 fErrorConfigFlag=0; //Initialize the error config. flag
189
190 //init error counters//
197
198 fRunningSum = 0;
199
203
204 fGoodEventCount = 0;
205
206 bEVENTCUTMODE = 0;
207
208 //std::cout<< "name = "<<name<<" error count same _HW = "<<fErrorCount_SameHW <<std::endl;
209 return;
210}
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.
Int_t fErrorCount_HWSat
check to see ADC channel is saturated
Int_t fADC_Same_NumEvt
Keep track of how many events with same ADC value returned.
void ClearEventData() override
Clear the event data in this element.
QwADC18_Channel * fRunningSum
Pointer to the running sum for this channel.
UInt_t fPreviousSequenceNumber
Previous event sequence number for this channel.
UInt_t fErrorConfigFlag
contains the global/local/stability flags
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, 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, fPreviousSequenceNumber, fRunningSum, 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(), QwADC18_Channel(), and QwADC18_Channel().

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

◆ InitializeChannel() [2/2]

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

Initialize the fields in this object.

Implements VQwHardwareChannel.

Definition at line 214 of file QwADC18_Channel.cc.

214 {
215 InitializeChannel(name,datatosave);
216 SetSubsystemName(subsystem);
217 SetModuleType(instrumenttype);
218 //PrintInfo();
219}
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:

◆ IsHeaderWord()

Bool_t QwADC18_Channel::IsHeaderWord ( UInt_t word) const

Decode the event data from a CODA buffer.

Definition at line 328 of file QwADC18_Channel.cc.

329{
330 return ((rawd & mask31x) != 0);
331}

References mask31x.

Referenced by ProcessEvBuffer().

+ Here is the caller graph for this function:

◆ LoadChannelParameters()

void QwADC18_Channel::LoadChannelParameters ( QwParameterFile & paramfile)
overridevirtual

Reimplemented from VQwDataElement.

Definition at line 221 of file QwADC18_Channel.cc.

221 {
222 UInt_t value = 0;
223 if (paramfile.ReturnValue("sample_size",value)){
225 } else {
226 QwWarning << "ADC18 Channel "
227 << GetElementName()
228 << " cannot set the default sample size."
229 << QwLog::endl;
230 }
231};
void SetDefaultSampleSize(size_t num_samples_map)
Bool_t ReturnValue(const std::string keyname, T &retvalue)

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.

24 {
25 Bool_t ldebug=kFALSE;
26 Double_t asym=0.0, mean=0.0, sigma=0.0;
27 Double_t amplitude=0.0, phase=0.0, frequency=0.0;
28
29 //Check to see if this line contains "drift"
30 if (paramfile.GetLine().find("drift")!=std::string::npos){
31 // "drift" appears somewhere in the line. Assume it is the next token and move on...
32 paramfile.GetNextToken(); //Throw away this token. Now the next three should be the drift parameters.
33 //read 3 parameters
34 amplitude = paramfile.GetTypedNextToken<Double_t>();
35 phase = paramfile.GetTypedNextToken<Double_t>();
36 frequency = paramfile.GetTypedNextToken<Double_t>(); // The frequency we read in should be in Hz.
37 this->AddRandomEventDriftParameters(amplitude, phase, frequency*Qw::Hz);
38 // std::cout << "In MQwMockable::LoadMockDataParameters: amp = " << amplitude << "\t phase = " << phase << "\t freq = " << frequency << std::endl;
39 }
40 else {
41 asym = paramfile.GetTypedNextToken<Double_t>();
42 mean = paramfile.GetTypedNextToken<Double_t>();
43 sigma = paramfile.GetTypedNextToken<Double_t>();
44 if (ldebug==1) {
45 std::cout << "#################### \n";
46 std::cout << "asym, mean, sigma \n" << std::endl;
47 std::cout << asym << " / "
48 << mean << " / "
49 << sigma << " / "
50 << std::endl;
51 }
52 this->SetRandomEventParameters(mean, sigma);
53 this->SetRandomEventAsymmetry(asym);
54 }
55}
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.

References MQwMockable::fMockDriftAmplitude, MQwMockable::fMockDriftFrequency, and MQwMockable::fMockDriftPhase.

◆ MatchNumberOfSamples()

Bool_t QwADC18_Channel::MatchNumberOfSamples ( size_t numsamp)

Definition at line 1237 of file QwADC18_Channel.cc.

1238{
1239 Bool_t status = kTRUE;
1240 if (!IsNameEmpty()){
1241 status = (fNumberOfSamples==numsamp);
1242 if (! status){
1243 if (bDEBUG)
1244 QwError << "QwADC18_Channel::MatchNumberOfSamples: Channel "
1245 << GetElementName()
1246 << " had fNumberOfSamples==" << fNumberOfSamples
1247 << " and was supposed to have " << numsamp
1248 << std::endl;
1249 }
1250 }
1251 return status;
1252}

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

Referenced by ApplyHWChecks().

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

◆ MatchSequenceNumber()

Bool_t QwADC18_Channel::MatchSequenceNumber ( size_t seqnum)

Definition at line 1228 of file QwADC18_Channel.cc.

1229{
1230 Bool_t status = kTRUE;
1231 if (!IsNameEmpty()){
1232 status = (fSequenceNumber == seqnum);
1233 }
1234 return status;
1235}

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 QwADC18_Channel::MultiplyBy ( const VQwHardwareChannel * valueptr)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 801 of file QwADC18_Channel.cc.

802{
803 const QwADC18_Channel* tmpptr;
804 tmpptr = dynamic_cast<const QwADC18_Channel*>(valueptr);
805 if (tmpptr!=NULL){
806 *this *= *tmpptr;
807 } else {
808 TString loc="Standard exception from QwADC18_Channel::MultiplyBy = "
809 +valueptr->GetElementName()+" is an incompatible type.";
810 throw std::invalid_argument(loc.Data());
811 }
812}

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

+ Here is the call graph for this function:

◆ operator*()

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

Definition at line 861 of file QwADC18_Channel.cc.

862{
863 QwADC18_Channel result = *this;
864 result *= value;
865 return result;
866}

References QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator*=() [1/2]

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

Definition at line 868 of file QwADC18_Channel.cc.

869{
870 if (!IsNameEmpty()){
871 this->fValue *= value.fValue;
872 this->fDiff_Raw = 0;
873 this->fPeak_Raw = 0;
874 this->fBase_Raw = 0;
875 this->fValueM2 = 0.0;
876 this->fErrorFlag |= (value.fErrorFlag);//error code is ORed.
877 }
878 return *this;
879}

References fBase_Raw, fDiff_Raw, VQwDataElement::fErrorFlag, fPeak_Raw, fValue, fValueM2, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator*=() [2/2]

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

Implements VQwHardwareChannel.

Definition at line 911 of file QwADC18_Channel.cc.

912{
913 const QwADC18_Channel* tmpptr;
914 tmpptr = dynamic_cast<const QwADC18_Channel*>(&source);
915 if (tmpptr!=NULL){
916 *this *= *tmpptr;
917 } else {
918 TString loc="Standard exception from QwADC18_Channel::operator*= "
919 +source.GetElementName()+" "
920 +this->GetElementName()+" are not of the same type";
921 throw(std::invalid_argument(loc.Data()));
922 }
923 return *this;
924}

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

+ Here is the call graph for this function:

◆ operator+()

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

Definition at line 827 of file QwADC18_Channel.cc.

828{
829 QwADC18_Channel result = *this;
830 result += value;
831 return result;
832}

References QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator+=() [1/2]

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

Definition at line 834 of file QwADC18_Channel.cc.

835{
836 if (!IsNameEmpty()) {
837 this->fValue += value.fValue;
838 this->fValueM2 = 0.0;
839 this->fErrorFlag |= value.fErrorFlag;//error code is ORed.
840 }
841 return *this;
842}

References VQwDataElement::fErrorFlag, fValue, fValueM2, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator+=() [2/2]

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

Implements VQwHardwareChannel.

Definition at line 881 of file QwADC18_Channel.cc.

882{
883 const QwADC18_Channel* tmpptr;
884 tmpptr = dynamic_cast<const QwADC18_Channel*>(&source);
885 if (tmpptr!=NULL){
886 *this += *tmpptr;
887 } else {
888 TString loc="Standard exception from QwADC18_Channel::operator+= "
889 +source.GetElementName()+" "
890 +this->GetElementName()+" are not of the same type";
891 throw(std::invalid_argument(loc.Data()));
892 }
893 return *this;
894}

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

+ Here is the call graph for this function:

◆ operator-()

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

Definition at line 844 of file QwADC18_Channel.cc.

845{
846 QwADC18_Channel result = *this;
847 result -= value;
848 return result;
849}

References QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator-=() [1/2]

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

Definition at line 851 of file QwADC18_Channel.cc.

852{
853 if (!IsNameEmpty()){
854 this->fValue -= value.fValue;
855 this->fValueM2 = 0.0;
856 this->fErrorFlag |= (value.fErrorFlag);//error code is ORed.
857 }
858 return *this;
859}

References VQwDataElement::fErrorFlag, fValue, fValueM2, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ operator-=() [2/2]

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

Implements VQwHardwareChannel.

Definition at line 896 of file QwADC18_Channel.cc.

897{
898 const QwADC18_Channel* tmpptr;
899 tmpptr = dynamic_cast<const QwADC18_Channel*>(&source);
900 if (tmpptr!=NULL){
901 *this -= *tmpptr;
902 } else {
903 TString loc="Standard exception from QwADC18_Channel::operator-= "
904 +source.GetElementName()+" "
905 +this->GetElementName()+" are not of the same type";
906 throw(std::invalid_argument(loc.Data()));
907 }
908 return *this;
909}

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

+ Here is the call graph for this function:

◆ operator/=() [1/2]

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

Definition at line 971 of file QwADC18_Channel.cc.

972{
973 // In this function, leave the "raw" variables untouched.
974 Double_t ratio;
975 Double_t variance;
976 if (!IsNameEmpty()) {
977 // The variances are calculated using the following formula:
978 // Var[ratio] = ratio^2 (Var[numer] / numer^2 + Var[denom] / denom^2)
979 //
980 // This requires that both the numerator and denominator are non-zero!
981 //
982 if (this->fValue != 0.0 && denom.fValue != 0.0){
983 ratio = (this->fValue) / (denom.fValue);
984 variance = ratio * ratio *
985 (this->fValueM2 / this->fValue / this->fValue
986 + denom.fValueM2 / denom.fValue / denom.fValue);
987 fValue = ratio;
988 fValueM2 = variance;
989 } else if (this->fValue == 0.0) {
990 fValue = 0.0;
991 fValueM2 = 0.0;
992 } else {
993 QwVerbose << "Attempting to divide by zero in "
995 fValue = 0.0;
996 fValueM2 = 0.0;
997 }
998
999 // Remaining variables
1000 // Don't change fGoodEventCount.
1001 // 'OR' the device error codes together.
1002 fErrorFlag |= denom.fErrorFlag;
1003 }
1004
1005 // Nanny
1006 if (fValue != fValue)
1007 QwWarning << "Angry Nanny: NaN detected in " << GetElementName() << QwLog::endl;
1008 return *this;
1009}
#define QwVerbose
Predefined log drain for verbose messages.
Definition QwLog.h:54

References QwLog::endl(), VQwDataElement::fErrorFlag, fValue, fValueM2, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), QwADC18_Channel(), QwVerbose, and QwWarning.

+ Here is the call graph for this function:

◆ operator/=() [2/2]

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

Implements VQwHardwareChannel.

Definition at line 926 of file QwADC18_Channel.cc.

927{
928 const QwADC18_Channel* tmpptr;
929 tmpptr = dynamic_cast<const QwADC18_Channel*>(&source);
930 if (tmpptr!=NULL){
931 *this /= *tmpptr;
932 } else {
933 TString loc="Standard exception from QwADC18_Channel::operator/= "
934 +source.GetElementName()+" "
935 +this->GetElementName()+" are not of the same type";
936 throw(std::invalid_argument(loc.Data()));
937 }
938 return *this;
939}

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

+ Here is the call graph for this function:

◆ operator=()

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

Definition at line 734 of file QwADC18_Channel.cc.

735{
736 if (this == &value) return *this;
737
738 if (!IsNameEmpty()) {
740 this->fDiff_Raw = value.fDiff_Raw;
741 this->fPeak_Raw = value.fPeak_Raw;
742 this->fBase_Raw = value.fBase_Raw;
743 this->fValue_Raw = value.fValue_Raw;
744 this->fValue = value.fValue;
745 this->fValueError = value.fValueError;
746 this->fValueM2 = value.fValueM2;
747 }
748 return *this;
749}
VQwHardwareChannel & operator=(const VQwHardwareChannel &value)
Arithmetic assignment operator: Should only copy event-based data.

References fBase_Raw, fDiff_Raw, fPeak_Raw, fValue, fValue_Raw, fValueError, fValueM2, VQwDataElement::IsNameEmpty(), VQwHardwareChannel::operator=(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ PrintErrorCounterHead()

void QwADC18_Channel::PrintErrorCounterHead ( )
static

Definition at line 1307 of file QwADC18_Channel.cc.

1308{
1309 TString message;
1310 message = Form("%30s","Device name");
1311 message += Form("%9s", "HW Sat");
1312 message += Form("%9s", "Sample");
1313 message += Form("%9s", "SW_HW");
1314 message += Form("%9s", "Sequence");
1315 message += Form("%9s", "SameHW");
1316 message += Form("%9s", "ZeroHW");
1317 message += Form("%9s", "EventCut");
1318 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1319 QwMessage << message << QwLog::endl;
1320 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1321}
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49

References QwLog::endl(), and QwMessage.

+ Here is the call graph for this function:

◆ PrintErrorCounters()

void QwADC18_Channel::PrintErrorCounters ( ) const
overridevirtual

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

Reimplemented from VQwDataElement.

Definition at line 1328 of file QwADC18_Channel.cc.

1329{
1330 TString message;
1332 message = Form("%30s", GetElementName().Data());
1333 message += Form("%9d", fErrorCount_HWSat);
1334 message += Form("%9d", fErrorCount_sample);
1335 message += Form("%9d", fErrorCount_SW_HW);
1336 message += Form("%9d", fErrorCount_Sequence);
1337 message += Form("%9d", fErrorCount_SameHW);
1338 message += Form("%9d", fErrorCount_ZeroHW);
1339 message += Form("%9d", fNumEvtsWithEventCutsRejected);
1340
1341 if((fDataToSave == kRaw) && (!kFoundPedestal||!kFoundGain)){
1342 message += " >>>>> No Pedestal or Gain in map file";
1343 }
1344
1345 QwMessage << message << QwLog::endl;
1346 }
1347}

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 QwADC18_Channel::PrintErrorCounterTail ( )
static

Definition at line 1323 of file QwADC18_Channel.cc.

1324{
1325 QwMessage << "---------------------------------------------------------------------------------------------" << QwLog::endl;
1326}

References QwLog::endl(), and QwMessage.

+ Here is the call graph for this function:

◆ PrintInfo()

void QwADC18_Channel::PrintInfo ( ) const
overridevirtual

Print multiple lines of information about this data element.

Reimplemented from VQwDataElement.

Definition at line 495 of file QwADC18_Channel.cc.

496{
497 QwMessage<<"***************************************"<<QwLog::endl;
498 QwMessage<<"Subsystem "<<GetSubsystemName()<<QwLog::endl;
499 QwMessage<<"Beam Instrument Type: "<<GetModuleType()<<QwLog::endl;
500 QwMessage<<"QwADC18 channel: "<<GetElementName()<<QwLog::endl;
501 QwMessage<<"fPedestal= "<< fPedestal<<QwLog::endl;
502 QwMessage<<"fCalibrationFactor= "<<fCalibrationFactor<<QwLog::endl;
503 QwMessage<<"fSequenceNumber= "<<fSequenceNumber<<QwLog::endl;
504 QwMessage<<"fNumberOfSamples= "<<fNumberOfSamples<<QwLog::endl;
505 QwMessage<<"fDiff_Raw= "<<fDiff_Raw<<QwLog::endl;
506 QwMessage<<"fPeak_Raw= "<<fPeak_Raw<<QwLog::endl;
507 QwMessage<<"fBase_Raw= "<<fBase_Raw<<QwLog::endl;
508 QwMessage<<"fValue = "<<std::setprecision(8) <<fValue << QwLog::endl;
509}
TString GetSubsystemName() const
Return the name of the inheriting subsystem name.
TString GetModuleType() const
Return the type of the beam instrument.

References QwLog::endl(), fBase_Raw, VQwHardwareChannel::fCalibrationFactor, fDiff_Raw, fNumberOfSamples, fPeak_Raw, VQwHardwareChannel::fPedestal, fSequenceNumber, fValue, VQwDataElement::GetElementName(), VQwDataElement::GetModuleType(), VQwDataElement::GetSubsystemName(), and QwMessage.

+ Here is the call graph for this function:

◆ PrintValue()

void QwADC18_Channel::PrintValue ( ) const
overridevirtual

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

Reimplemented from VQwDataElement.

Definition at line 1176 of file QwADC18_Channel.cc.

1177{
1178 QwMessage << std::setprecision(8)
1179 << std::setw(18) << std::left << GetSubsystemName() << ""
1180 << std::setw(18) << std::left << GetModuleType() << ""
1181 << std::setw(18) << std::left << GetElementName() << ""
1182 << std::setw(12) << std::left << GetValue() << " +/- "
1183 << std::setw(12) << std::left << GetValueError() << " sig "
1184 << std::setw(12) << std::left << GetValueWidth() << ""
1185 << std::setw(12) << std::left << fGoodEventCount << ""
1186 << QwLog::endl;
1187}

References QwLog::endl(), VQwDataElement::fGoodEventCount, VQwDataElement::GetElementName(), VQwDataElement::GetModuleType(), VQwDataElement::GetSubsystemName(), GetValue(), GetValueError(), GetValueWidth(), and QwMessage.

+ Here is the call graph for this function:

◆ ProcessDataWord()

Int_t QwADC18_Channel::ProcessDataWord ( UInt_t word)

Definition at line 333 of file QwADC18_Channel.cc.

334{
335 // "Actual" values from data word
336 UInt_t act_dtype = (rawd & mask2422x) >> 22;
337 [[maybe_unused]]
338 UInt_t act_chan = (act_dtype != 4) ?
339 ((rawd & mask3029x) >> 29) : 0;
340 UInt_t act_dvalue = (rawd & mask2625x) >> 25;
341 UInt_t act_snum = (act_dtype == 1 || act_dtype == 2) ?
342 ((rawd & mask2118x) >> 18) : 0;
343
344 // Interpret by data type
345 UInt_t value_raw = 0;
346 switch (act_dtype) {
347 case 0: // Diff word
348 static UInt_t prev_dvalue = act_dvalue;
349 if (act_dvalue != prev_dvalue) {
350 QwError << "QwADC18_Channel::ProcessEvBuffer: Number of samples changed " << act_dvalue << " " << prev_dvalue << QwLog::endl;
351 return 0;
352 }
353 value_raw = rawd & mask200x;
354 if (rawd & mask21x) value_raw = -((~value_raw & 0x1fffffff) + 1);
355 fNumberOfSamples = (1 << act_dvalue);
356 return value_raw;
357 break;
358 case 1: // Peak word
359 case 2: // Base word
360 if (act_snum != fNumberOfSamples) {
361 QwError << "QwADC18_Channel::ProcessEvBuffer: Number of samples changed " << act_snum << " " << fNumberOfSamples << QwLog::endl;
362 return 0;
363 }
364 if (act_dvalue != 0) {
365 QwError << "QwADC18_Channel::ProcessEvBuffer: Divider value non-zero 0x" << std::hex << rawd << std::dec << QwLog::endl;
366 return 0;
367 }
368 return rawd & mask170x;
369 break;
370 case 4: // DAC word
371 if (act_dvalue != 0) {
372 QwError << "QwADC18_Channel::ProcessEvBuffer: Divider value non-zero 0x" << std::hex << rawd << std::dec << QwLog::endl;
373 return 0;
374 }
375 return rawd & mask150x;
376 break;
377 default:
378 QwError << "QwADC18_Channel::ProcessEvBuffer: Unknown data type 0x" << std::hex << rawd << std::dec << QwLog::endl;
379 return 0;
380 }
381}
static const UInt_t mask200x
static const UInt_t mask2118x
static const UInt_t mask21x
static const UInt_t mask150x

References QwLog::endl(), fNumberOfSamples, mask150x, mask170x, mask200x, mask2118x, mask21x, mask2422x, mask2625x, mask3029x, and QwError.

Referenced by ProcessEvBuffer().

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

◆ ProcessEvBuffer()

Int_t QwADC18_Channel::ProcessEvBuffer ( UInt_t * buffer,
UInt_t num_words_left,
UInt_t subelement = 0 )
overridevirtual

Process the CODA event buffer for this element.

Implements VQwDataElement.

Definition at line 384 of file QwADC18_Channel.cc.

385{
386 Bool_t debug = false;
387
388 // Print buffer
389 if (debug) {
390 QwOut << GetElementName() << " : " << QwLog::endl << std::hex;
391 UInt_t n = 25;
392 for (size_t i = 0; i < num_words_left; i++) {
393 QwOut << "0x" << std::setfill('0') << std::setw(8) << buffer[i] << " ";
394 if (i % n == n - 1) QwOut << QwLog::endl;
395 }
396 QwOut << std::dec << std::setfill(' ') << std::setw(0) << QwLog::endl;
397 }
398
399 UInt_t words_read = 0;
400 if (IsNameEmpty()) {
401 // This channel is not used, but is present in the data stream.
402 // Skip over this data.
403 words_read = kDataWordsPerChannel;
404 } else if (num_words_left >= fNumberOfDataWords) {
405
406 // Is this a header word?
407 if (IsHeaderWord(buffer[0])) {
408
409 // Debug output
410 if (debug) {
411 QwOut << " : header " << std::hex;
412 UInt_t n = kHeaderWordsPerModule;
413 for (size_t i = 0; i < n && i < num_words_left; i++) {
414 QwOut << "0x" << std::setfill('0') << std::setw(8) << buffer[i] << " ";
415 }
416 QwOut << std::dec << std::setfill(' ') << std::setw(0) << QwLog::endl;
417 }
418
419 // Check if enough words left
420 if (num_words_left < kHeaderWordsPerModule) {
421 QwError << "QwADC18_Channel::ProcessEvBuffer: Not enough words left!" << QwLog::endl;
422 return num_words_left;
423 }
424
425 // FIXME Catch 0xfa180bad
426
427 // Header word: read DAC value
429
430 words_read = kHeaderWordsPerModule;
431
432 } else {
433
434 // Debug output
435 if (debug) {
436 QwOut << " : channel " << std::hex;
437 UInt_t n = kDataWordsPerChannel;
438 for (size_t i = 0; i < n && i < num_words_left; i++) {
439 QwOut << "0x" << std::setfill('0') << std::setw(8) << buffer[i] << " ";
440 }
441 QwOut << std::dec << std::setfill(' ') << std::setw(0) << QwLog::endl;
442 }
443
444 // Check if enough words left
445 if (num_words_left < kDataWordsPerChannel) {
446 QwError << "QwADC18_Channel::ProcessEvBuffer: Not enough words left!" << QwLog::endl;
447 return num_words_left;
448 }
449
450 // Data channel words: read diff, peak, base
451 fDiff_Raw = ProcessDataWord(buffer[0]);
452 fPeak_Raw = ProcessDataWord(buffer[1]);
453 fBase_Raw = ProcessDataWord(buffer[2]);
454
455 words_read = kDataWordsPerChannel;
456 }
457
458 } else {
459 QwError << "QwADC18_Channel::ProcessEvBuffer: Not enough words!" << QwLog::endl;
460 }
461
462 return words_read;
463}
#define QwOut
Predefined log drain for explicit output.
Definition QwLog.h:34
Bool_t IsHeaderWord(UInt_t word) const
Decode the event data from a CODA buffer.
Int_t ProcessDataWord(UInt_t word)
UInt_t fNumberOfDataWords
Number of raw data words in this data element.

References QwLog::endl(), fBase_Raw, fDiff_Raw, VQwHardwareChannel::fNumberOfDataWords, fPeak_Raw, fValue_Raw, VQwDataElement::GetElementName(), IsHeaderWord(), VQwDataElement::IsNameEmpty(), kDataWordsPerChannel, kHeaderWordsPerModule, ProcessDataWord(), QwError, and QwOut.

+ Here is the call graph for this function:

◆ ProcessEvent()

void QwADC18_Channel::ProcessEvent ( )
overridevirtual

Process the event data according to pedestal and calibration factor.

Implements VQwHardwareChannel.

Definition at line 467 of file QwADC18_Channel.cc.

468{
469 if (fNumberOfSamples == 0 && fDiff_Raw == 0) {
470 // There isn't valid data for this channel. Just flag it and move on.
471 fValue = 0.0;
472 fValueM2 = 0.0;
474 } else if (fNumberOfSamples == 0) {
475 // This is probably a more serious problem.
476 QwWarning << "QwADC18_Channel::ProcessEvent: Channel "
477 << GetElementName()
478 << " has fNumberOfSamples == 0 but has valid data in the value filed."
479 << "Flag this as an error."
480 << QwLog::endl;
481 fValue = 0.0;
482 fValueM2 = 0.0;
484 } else {
486 fValueM2 = 0.0; // second moment is zero for single events
487 }
488}

References QwLog::endl(), VQwHardwareChannel::fCalibrationFactor, fDiff_Raw, VQwDataElement::fErrorFlag, fNumberOfSamples, VQwHardwareChannel::fPedestal, fValue, fValueM2, VQwDataElement::GetElementName(), kErrorFlag_sample, and QwWarning.

+ Here is the call graph for this function:

◆ Product()

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

Definition at line 1022 of file QwADC18_Channel.cc.

1023{
1024 if (!IsNameEmpty()){
1025 fValue = value1.fValue * value2.fValue;
1026 fDiff_Raw = 0;
1027 fPeak_Raw = 0;
1028 fBase_Raw = 0;
1029
1030 // Remaining variables
1032 fErrorFlag = (value1.fErrorFlag | value2.fErrorFlag);//error code is ORed.
1033 }
1034}

References fBase_Raw, fDiff_Raw, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fPeak_Raw, fValue, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ RandomizeEventData()

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

Internally generate random event data.

Implements MQwMockable.

Definition at line 248 of file QwADC18_Channel.cc.

249{
250 // Calculate drift (if time is not specified, it stays constant at zero)
251 Double_t drift = 0.0;
252 for (UInt_t i = 0; i < fMockDriftFrequency.size(); i++) {
253 drift += fMockDriftAmplitude[i] * sin(2.0 * Qw::pi * fMockDriftFrequency[i] * time + fMockDriftPhase[i]);
254 }
255
256 Double_t value = fMockGaussianMean * (1 + helicity * fMockAsymmetry)
258 + drift;
259
260 fValue = value;
261 fSequenceNumber = 0;
263}
static const double pi
Angles: base unit is radian.
Definition QwUnits.h:107
Double_t GetRandomValue()

References MQwMockable::fMockAsymmetry, MQwMockable::fMockDriftAmplitude, MQwMockable::fMockDriftFrequency, MQwMockable::fMockDriftPhase, MQwMockable::fMockGaussianMean, MQwMockable::fMockGaussianSigma, fNumberOfSamples, fNumberOfSamples_map, fSequenceNumber, fValue, MQwMockable::GetRandomValue(), and Qw::pi.

+ Here is the call graph for this function:

◆ Ratio()

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

Definition at line 954 of file QwADC18_Channel.cc.

955{
956 if (!IsNameEmpty()) {
957 *this = numer;
958 *this /= denom;
959
960 // Set the raw values to zero.
961 fDiff_Raw = 0;
962 fPeak_Raw = 0;
963 fBase_Raw = 0;
964
965 // Remaining variables
967 fErrorFlag = (numer.fErrorFlag|denom.fErrorFlag);//error code is ORed.
968 }
969}

References fBase_Raw, fDiff_Raw, VQwDataElement::fErrorFlag, VQwDataElement::fGoodEventCount, fPeak_Raw, VQwDataElement::IsNameEmpty(), and QwADC18_Channel().

+ Here is the call graph for this function:

◆ Scale()

void QwADC18_Channel::Scale ( Double_t Offset)
overridevirtual

Implements VQwHardwareChannel.

Definition at line 1047 of file QwADC18_Channel.cc.

1048{
1049 if (!IsNameEmpty())
1050 {
1051 fValue *= scale;
1052 fValueM2 *= scale * scale;
1053 }
1054}

References fValue, fValueM2, and VQwDataElement::IsNameEmpty().

+ Here is the call graph for this function:

◆ ScaledAdd()

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

Implements VQwHardwareChannel.

Definition at line 1349 of file QwADC18_Channel.cc.

1350{
1351 const QwADC18_Channel* input = dynamic_cast<const QwADC18_Channel*>(value);
1352
1353 // follows same steps as += but w/ scaling factor
1354 if(input!=NULL && !IsNameEmpty()){
1355 // QwWarning << "Adding " << input->GetElementName()
1356 // << " to " << GetElementName()
1357 // << " with scale factor " << scale
1358 // << QwLog::endl;
1359 // PrintValue();
1360 // input->PrintValue();
1361 this->fDiff_Raw = 0;
1362 this->fPeak_Raw = 0;
1363 this->fBase_Raw = 0;
1364 this->fValue += scale * input->fValue;
1365 this->fValueM2 = 0.0;
1366 this->fNumberOfSamples += input->fNumberOfSamples;
1367 this->fSequenceNumber = 0;
1368 this->fErrorFlag |= (input->fErrorFlag);
1369 } else if (input == NULL && value != NULL) {
1370 TString loc="Standard exception from QwADC18_Channel::ScaledAdd "
1371 +value->GetElementName()+" "
1372 +this->GetElementName()+" are not of the same type";
1373 throw(std::invalid_argument(loc.Data()));
1374 }
1375}

References fBase_Raw, fDiff_Raw, VQwDataElement::fErrorFlag, fNumberOfSamples, fPeak_Raw, fSequenceNumber, fValue, fValueM2, VQwDataElement::GetElementName(), VQwDataElement::IsNameEmpty(), QwADC18_Channel(), and VQwHardwareChannel::VQwHardwareChannel().

+ Here is the call graph for this function:

◆ SetADC18SaturationLimt()

void QwADC18_Channel::SetADC18SaturationLimt ( Double_t sat_volts = 8.5)
inline

Definition at line 205 of file QwADC18_Channel.h.

205 {
206 fSaturationABSLimit = sat_volts;
207 }

References fSaturationABSLimit.

Referenced by QwADC18_Channel(), and QwADC18_Channel().

+ Here is the caller graph for this function:

◆ SetCalibrationToVolts()

void QwADC18_Channel::SetCalibrationToVolts ( )
inline

Definition at line 244 of file QwADC18_Channel.h.

References kADC18_VoltsPerBit, and VQwHardwareChannel::SetCalibrationFactor().

+ Here is the call graph for this function:

◆ SetDefaultSampleSize()

void QwADC18_Channel::SetDefaultSampleSize ( size_t num_samples_map)
inline

Definition at line 100 of file QwADC18_Channel.h.

100 {
101 // This will be checked against the no.of samples read by the module
102 fNumberOfSamples_map = num_samples_map;
103 };

References fNumberOfSamples_map.

Referenced by LoadChannelParameters().

+ Here is the caller graph for this function:

◆ SetEventData()

void QwADC18_Channel::SetEventData ( Double_t value)

Definition at line 277 of file QwADC18_Channel.cc.

278{
279 fValue = value;
280 fDiff_Raw = (UInt_t) value;
281 fPeak_Raw = (UInt_t) value;
282 fBase_Raw = (UInt_t) 0;
283 fSequenceNumber = 0;
285}

References fBase_Raw, fDiff_Raw, fNumberOfSamples, fNumberOfSamples_map, fPeak_Raw, fSequenceNumber, and fValue.

◆ SetRawEventData()

void QwADC18_Channel::SetRawEventData ( )
overridevirtual

Implements MQwMockable.

Definition at line 287 of file QwADC18_Channel.cc.

References VQwHardwareChannel::fCalibrationFactor, fDiff_Raw, fNumberOfSamples, fNumberOfSamples_map, fPeak_Raw, VQwHardwareChannel::fPedestal, and fValue.

Referenced by SmearByResolution().

+ Here is the caller graph for this function:

◆ SmearByResolution()

void QwADC18_Channel::SmearByResolution ( double resolution)
overridevirtual

Implements MQwMockable.

Definition at line 265 of file QwADC18_Channel.cc.

266{
267 fValue = resolution*GetRandomValue();
268 fValueM2 = 0.0;
271}
void SetRawEventData() override

References fNumberOfSamples, fNumberOfSamples_map, fValue, fValueM2, MQwMockable::GetRandomValue(), and SetRawEventData().

+ Here is the call graph for this function:

◆ SubtractValueFrom()

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

Implements VQwHardwareChannel.

Definition at line 789 of file QwADC18_Channel.cc.

790{
791 const QwADC18_Channel* tmpptr;
792 tmpptr = dynamic_cast<const QwADC18_Channel*>(valueptr);
793 if (tmpptr!=NULL){
794 *this -= *tmpptr;
795 } else {
796 TString loc="Standard exception from QwADC18_Channel::SubtractValueFrom = "
797 +valueptr->GetElementName()+" is an incompatible type.";
798 throw std::invalid_argument(loc.Data());
799 }
800}

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

+ Here is the call graph for this function:

◆ Sum()

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

Definition at line 942 of file QwADC18_Channel.cc.

943{
944 *this = value1;
945 *this += value2;
946}

References QwADC18_Channel().

+ Here is the call graph for this function:

Friends And Related Symbol Documentation

◆ operator<<

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

Definition at line 1189 of file QwADC18_Channel.cc.

1190{
1191 stream << channel.GetValue();
1192 return stream;
1193}
Double_t GetValue(size_t element) const override

References GetValue(), and QwADC18_Channel().

Field Documentation

◆ bBlock

Bool_t QwADC18_Channel::bBlock
private

Definition at line 346 of file QwADC18_Channel.h.

◆ bBlock_raw

Bool_t QwADC18_Channel::bBlock_raw
private

Definition at line 347 of file QwADC18_Channel.h.

◆ bDEBUG

const Bool_t QwADC18_Channel::bDEBUG =kFALSE
staticprivate

debugging display purposes

For ADC18 data element trimming uses

Definition at line 341 of file QwADC18_Channel.h.

Referenced by ApplyHWChecks(), and MatchNumberOfSamples().

◆ bDevice_Error_Code

Bool_t QwADC18_Channel::bDevice_Error_Code
private

Definition at line 349 of file QwADC18_Channel.h.

◆ bHw_sum

Bool_t QwADC18_Channel::bHw_sum
private

Definition at line 344 of file QwADC18_Channel.h.

◆ bHw_sum_raw

Bool_t QwADC18_Channel::bHw_sum_raw
private

Definition at line 345 of file QwADC18_Channel.h.

◆ bNum_samples

Bool_t QwADC18_Channel::bNum_samples
private

Definition at line 348 of file QwADC18_Channel.h.

◆ bSequence_number

Bool_t QwADC18_Channel::bSequence_number
private

Definition at line 350 of file QwADC18_Channel.h.

◆ fADC_Same_NumEvt

Int_t QwADC18_Channel::fADC_Same_NumEvt
private

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

Definition at line 331 of file QwADC18_Channel.h.

Referenced by InitializeChannel().

◆ fBase_Raw

UInt_t QwADC18_Channel::fBase_Raw
private

◆ fDAC

QwADC18_Channel* QwADC18_Channel::fDAC
private

Pointer to the DAC channel for this channel.

Definition at line 301 of file QwADC18_Channel.h.

◆ fDiff_Raw

◆ fErrorCount_HWSat

Int_t QwADC18_Channel::fErrorCount_HWSat
private

check to see ADC channel is saturated

Definition at line 320 of file QwADC18_Channel.h.

Referenced by InitializeChannel(), and PrintErrorCounters().

◆ fErrorCount_SameHW

Int_t QwADC18_Channel::fErrorCount_SameHW
private

check to see ADC returning same HW value

Definition at line 324 of file QwADC18_Channel.h.

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

◆ fErrorCount_sample

Int_t QwADC18_Channel::fErrorCount_sample
private

for sample size check

Definition at line 321 of file QwADC18_Channel.h.

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

◆ fErrorCount_Sequence

Int_t QwADC18_Channel::fErrorCount_Sequence
private

sequence number check

Definition at line 323 of file QwADC18_Channel.h.

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

◆ fErrorCount_SW_HW

Int_t QwADC18_Channel::fErrorCount_SW_HW
private

HW_sum==SW_sum check.

Definition at line 322 of file QwADC18_Channel.h.

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

◆ fErrorCount_ZeroHW

Int_t QwADC18_Channel::fErrorCount_ZeroHW
private

check to see ADC returning zero

Definition at line 325 of file QwADC18_Channel.h.

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

◆ fNumberOfSamples

◆ fNumberOfSamples_map

UInt_t QwADC18_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 317 of file QwADC18_Channel.h.

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

◆ fNumEvtsWithEventCutsRejected

Int_t QwADC18_Channel::fNumEvtsWithEventCutsRejected
private

Counts the Event cut rejected events.

Definition at line 327 of file QwADC18_Channel.h.

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

◆ fPeak_Raw

◆ fPrev_HardwareBlockSum

Double_t QwADC18_Channel::fPrev_HardwareBlockSum
private

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

Definition at line 334 of file QwADC18_Channel.h.

◆ fPreviousSequenceNumber

UInt_t QwADC18_Channel::fPreviousSequenceNumber
private

Previous event sequence number for this channel.

Definition at line 315 of file QwADC18_Channel.h.

Referenced by InitializeChannel().

◆ fRunningSum

QwADC18_Channel* QwADC18_Channel::fRunningSum
private

Pointer to the running sum for this channel.

Definition at line 298 of file QwADC18_Channel.h.

Referenced by InitializeChannel().

◆ fSaturationABSLimit

Double_t QwADC18_Channel::fSaturationABSLimit
private

absolute value of the ADC18 saturation volt

Definition at line 338 of file QwADC18_Channel.h.

Referenced by GetADC18SaturationLimt(), QwADC18_Channel(), QwADC18_Channel(), and SetADC18SaturationLimt().

◆ fSequenceNo_Counter

Int_t QwADC18_Channel::fSequenceNo_Counter
private

Internal counter to keep track of the sequence number.

Definition at line 333 of file QwADC18_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fSequenceNo_Prev

Int_t QwADC18_Channel::fSequenceNo_Prev
private

Keep the sequence number of the last event.

Definition at line 332 of file QwADC18_Channel.h.

Referenced by ApplyHWChecks(), and InitializeChannel().

◆ fSequenceNumber

UInt_t QwADC18_Channel::fSequenceNumber
private

Event sequence number for this channel.

Definition at line 314 of file QwADC18_Channel.h.

Referenced by ClearEventData(), GetSequenceNumber(), MatchSequenceNumber(), PrintInfo(), RandomizeEventData(), ScaledAdd(), and SetEventData().

◆ fValue

◆ fValue_Raw

UInt_t QwADC18_Channel::fValue_Raw
private

Definition at line 270 of file QwADC18_Channel.h.

Referenced by FillTreeVector(), GetRawValue(), operator=(), and ProcessEvBuffer().

◆ fValueError

Double_t QwADC18_Channel::fValueError
private

◆ fValueM2

◆ kADC18_VoltsPerBit

const Double_t QwADC18_Channel::kADC18_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 305 of file QwADC18_Channel.h.

Referenced by GetAverageVolts(), and SetCalibrationToVolts().

◆ kDataWordsPerChannel

const Int_t QwADC18_Channel::kDataWordsPerChannel = 3
staticprivate

Definition at line 282 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), GetBufferOffset(), and ProcessEvBuffer().

◆ kDEBUG

const Bool_t QwADC18_Channel::kDEBUG = kFALSE
staticprivate

Definition at line 277 of file QwADC18_Channel.h.

◆ kFooterWordsPerBank

const Int_t QwADC18_Channel::kFooterWordsPerBank = 1
staticprivate

Definition at line 279 of file QwADC18_Channel.h.

◆ kFooterWordsPerModule

const Int_t QwADC18_Channel::kFooterWordsPerModule = 1
staticprivate

Definition at line 281 of file QwADC18_Channel.h.

Referenced by GetBufferOffset().

◆ kHeaderWordsPerBank

const Int_t QwADC18_Channel::kHeaderWordsPerBank = 1
staticprivate

Definition at line 278 of file QwADC18_Channel.h.

Referenced by GetBufferOffset().

◆ kHeaderWordsPerModule

const Int_t QwADC18_Channel::kHeaderWordsPerModule = 12
staticprivate

Definition at line 280 of file QwADC18_Channel.h.

Referenced by GetBufferOffset(), and ProcessEvBuffer().

◆ kMaxChannels

const Int_t QwADC18_Channel::kMaxChannels = 4
staticprivate

Definition at line 283 of file QwADC18_Channel.h.

Referenced by GetBufferOffset().

◆ kTimePerSample

const Double_t QwADC18_Channel::kTimePerSample = 2.0 * Qw::us
static

Definition at line 47 of file QwADC18_Channel.h.

◆ mask150x

const UInt_t QwADC18_Channel::mask150x = 0x0000ffff
staticprivate

Definition at line 294 of file QwADC18_Channel.h.

Referenced by ProcessDataWord().

◆ mask170x

const UInt_t QwADC18_Channel::mask170x = 0x0003ffff
staticprivate

Definition at line 293 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), and ProcessDataWord().

◆ mask200x

const UInt_t QwADC18_Channel::mask200x = 0x001fffff
staticprivate

Definition at line 291 of file QwADC18_Channel.h.

Referenced by ProcessDataWord().

◆ mask2118x

const UInt_t QwADC18_Channel::mask2118x = 0x003c0000
staticprivate

Definition at line 292 of file QwADC18_Channel.h.

Referenced by ProcessDataWord().

◆ mask21x

const UInt_t QwADC18_Channel::mask21x = 0x00200000
staticprivate

Definition at line 290 of file QwADC18_Channel.h.

Referenced by ProcessDataWord().

◆ mask2422x

const UInt_t QwADC18_Channel::mask2422x = 0x01c00000
staticprivate

Definition at line 289 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), and ProcessDataWord().

◆ mask2625x

const UInt_t QwADC18_Channel::mask2625x = 0x06000000
staticprivate

Definition at line 288 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), and ProcessDataWord().

◆ mask3029x

const UInt_t QwADC18_Channel::mask3029x = 0x60000000
staticprivate

Definition at line 287 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), and ProcessDataWord().

◆ mask31x

const UInt_t QwADC18_Channel::mask31x = 0x80000000
staticprivate

Definition at line 286 of file QwADC18_Channel.h.

Referenced by EncodeEventData(), and IsHeaderWord().


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