QwAnalysis
QwSIS3320_Accumulator Class Reference

SIS3320 sampling ADC accumulator. More...

#include <QwSIS3320_Accumulator.h>

+ Inheritance diagram for QwSIS3320_Accumulator:
+ Collaboration diagram for QwSIS3320_Accumulator:

Public Member Functions

 QwSIS3320_Accumulator (TString name="")
 
virtual ~QwSIS3320_Accumulator ()
 
Int_t GetMinValue () const
 
Int_t GetMaxValue () const
 
Int_t GetMinTime () const
 
Int_t GetMaxTime () const
 
Double_t GetNumberOfSamples () const
 
Double_t GetAccumulatorSum () const
 
Double_t GetAccumulatorAvg () const
 
void ClearEventData ()
 Clear the event data in this element. More...
 
Int_t ProcessEvBuffer (UInt_t *buffer, UInt_t num_words_left, UInt_t subelement=0)
 Process the CODA event buffer for this element. More...
 
void ProcessEvent ()
 
const QwSIS3320_Accumulator operator/ (const Double_t &value) const
 
const QwSIS3320_Accumulator operator* (const Double_t &value) const
 
const QwSIS3320_Accumulator operator+ (const Double_t &value) const
 
const QwSIS3320_Accumulator operator- (const Double_t &value) const
 
QwSIS3320_Accumulatoroperator/= (const Double_t &value)
 
QwSIS3320_Accumulatoroperator*= (const Double_t &value)
 
QwSIS3320_Accumulatoroperator+= (const Double_t &value)
 
QwSIS3320_Accumulatoroperator-= (const Double_t &value)
 
const QwSIS3320_Accumulator operator+ (const QwSIS3320_Accumulator &value) const
 
const QwSIS3320_Accumulator operator- (const QwSIS3320_Accumulator &value) const
 
QwSIS3320_Accumulatoroperator= (const QwSIS3320_Accumulator &value)
 
QwSIS3320_Accumulatoroperator+= (const QwSIS3320_Accumulator &value)
 
QwSIS3320_Accumulatoroperator-= (const QwSIS3320_Accumulator &value)
 
void Sum (const QwSIS3320_Accumulator &value1, const QwSIS3320_Accumulator &value2)
 
void Difference (const QwSIS3320_Accumulator &value1, const QwSIS3320_Accumulator &value2)
 
void Ratio (const QwSIS3320_Accumulator &numer, const QwSIS3320_Accumulator &denom)
 
void ConstructHistograms (TDirectory *folder, TString &prefix)
 Construct the histograms for this data element. More...
 
void FillHistograms ()
 Fill the histograms for this data element. More...
 
void ConstructBranchAndVector (TTree *tree, TString &prefix, std::vector< Double_t > &values)
 
void FillTreeVector (std::vector< Double_t > &values) const
 
- Public Member Functions inherited from VQwDataElement
 VQwDataElement ()
 Default constructor. More...
 
 VQwDataElement (const VQwDataElement &value)
 Copy constructor. More...
 
virtual ~VQwDataElement ()
 Virtual destructor. More...
 
Bool_t IsNameEmpty () const
 Is the name of this element empty? More...
 
void SetElementName (const TString &name)
 Set the name of this element. More...
 
virtual const TString & GetElementName () const
 Get the name of this element. More...
 
virtual void LoadChannelParameters (QwParameterFile &paramfile)
 
size_t GetNumberOfDataWords ()
 Get the number of data words in this data element. More...
 
UInt_t GetGoodEventCount () const
 
virtual void AssignValueFrom (const VQwDataElement *valueptr)
 
virtual VQwDataElementoperator+= (const VQwDataElement &value)
 Addition-assignment operator. More...
 
virtual VQwDataElementoperator-= (const VQwDataElement &value)
 Subtraction-assignment operator. More...
 
virtual void Sum (const VQwDataElement &value1, const VQwDataElement &value2)
 Sum operator. More...
 
virtual void Difference (const VQwDataElement &value1, const VQwDataElement &value2)
 Difference operator. More...
 
virtual void Ratio (const VQwDataElement &numer, const VQwDataElement &denom)
 Ratio operator. More...
 
virtual void PrintValue () const
 Print single line of value and error of this data element. More...
 
virtual void PrintInfo () const
 Print multiple lines of information about this data element. More...
 
virtual void SetSingleEventCuts (UInt_t errorflag, Double_t min, Double_t max, Double_t stability)
 set the upper and lower limits (fULimit and fLLimit), stability % and the error flag on this channel More...
 
virtual void PrintErrorCounters () const
 report number of events failed due to HW and event cut failure More...
 
virtual UInt_t GetEventcutErrorFlag ()
 return the error flag on this channel/device More...
 
virtual UInt_t UpdateErrorFlag ()
 Update the error flag based on the error flags of internally contained objects Return paramter is the "Eventcut Error Flag". More...
 
virtual void SetNeedsExternalClock (Bool_t needed)
 
virtual Bool_t NeedsExternalClock ()
 
virtual std::string GetExternalClockName ()
 
virtual void SetExternalClockPtr (const VQwHardwareChannel *clock)
 
virtual void SetExternalClockName (const std::string name)
 
virtual Double_t GetNormClockValue ()
 
TString GetSubsystemName () const
 Return the name of the inheriting subsystem name. More...
 
void SetSubsystemName (TString sysname)
 Set the name of the inheriting subsystem name. More...
 
TString GetModuleType () const
 Return the type of the beam instrument. More...
 
void SetModuleType (TString ModuleType)
 set the type of the beam instrument More...
 
- Public Member Functions inherited from MQwHistograms
void ShareHistograms (const MQwHistograms *source)
 Share histogram pointers between objects. More...
 

Protected Attributes

Double_t fNumberOfSamples
 Number of accumulated samples. More...
 
Double_t fAccumulatorSum
 Accumulator sum. More...
 
Double_t fAccumulatorAvg
 Accumulator average. More...
 
- Protected Attributes inherited from VQwDataElement
TString fElementName
 Name of this data element. More...
 
UInt_t fNumberOfDataWords
 Number of raw data words in this data element. More...
 
Int_t fGoodEventCount
 Number of good events accumulated in this element. More...
 
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) More...
 
UInt_t fErrorConfigFlag
 contains the global/local/stability flags More...
 
- Protected Attributes inherited from MQwHistograms
std::vector< TH1_ptrfHistograms
 Histograms associated with this data element. More...
 

Private Attributes

int fMinValue
 
int fMaxValue
 
int fMinTime
 Value-based accumulator limits. More...
 
int fMaxTime
 
size_t fTreeArrayIndex
 Time-based accumulator limits. More...
 
size_t fTreeArrayNumEntries
 Index of this data element in tree. More...
 

Static Private Attributes

static const unsigned int INDEX_NUM = 0
 Number of entries from this data element. More...
 
static const unsigned int INDEX_SUM = 1
 

Friends

std::ostream & operator<< (std::ostream &stream, const QwSIS3320_Accumulator &a)
 

Additional Inherited Members

- Public Types inherited from VQwDataElement
enum  EDataToSave { kRaw = 0, kDerived }
 
- Protected Member Functions inherited from VQwDataElement
void SetNumberOfDataWords (const UInt_t &numwords)
 Set the number of data words in this data element. More...
 
virtual VQwDataElementoperator= (const VQwDataElement &value)
 Arithmetic assignment operator: Should only copy event-based data. More...
 
virtual void UpdateErrorFlag (const UInt_t &error)
 
- Protected Member Functions inherited from MQwHistograms
 MQwHistograms ()
 Default constructor. More...
 
 MQwHistograms (const MQwHistograms &source)
 Copy constructor. More...
 
virtual ~MQwHistograms ()
 Virtual destructor. More...
 
virtual MQwHistogramsoperator= (const MQwHistograms &value)
 
void Fill_Pointer (TH1_ptr hist_ptr, Double_t value)
 
void AddHistogram (TH1 *h)
 Register a histogram. More...
 

Detailed Description

SIS3320 sampling ADC accumulator.

Author
W. Deconinck
Date
2009-09-04 18:06:23

The QwSIS3320_Accumulator should allow convenient access to the accumulator data collected with the SIS3320 for the Compton photon detector. This class implements its own sum, difference, and ratio methods inherited from the general VQwDataElement.

Definition at line 34 of file QwSIS3320_Accumulator.h.

Constructor & Destructor Documentation

QwSIS3320_Accumulator::QwSIS3320_Accumulator ( TString  name = "")
inline

Definition at line 38 of file QwSIS3320_Accumulator.h.

References fMaxTime, fMaxValue, fMinTime, fMinValue, VQwDataElement::SetElementName(), and VQwDataElement::SetNumberOfDataWords().

38  {
39  SetElementName(name);
41  fMinValue = 0; fMaxValue = 0;
42  fMinTime = 0; fMaxTime = 0;
43  };
int fMinTime
Value-based accumulator limits.
void SetElementName(const TString &name)
Set the name of this element.
void SetNumberOfDataWords(const UInt_t &numwords)
Set the number of data words in this data element.

+ Here is the call graph for this function:

virtual QwSIS3320_Accumulator::~QwSIS3320_Accumulator ( )
inlinevirtual

Definition at line 44 of file QwSIS3320_Accumulator.h.

44 { };

Member Function Documentation

void QwSIS3320_Accumulator::ClearEventData ( )
inlinevirtual

Clear the event data in this element.

Reimplemented from VQwDataElement.

Definition at line 59 of file QwSIS3320_Accumulator.h.

References fAccumulatorSum, and fNumberOfSamples.

59 { fAccumulatorSum = 0; fNumberOfSamples = 0; };
Double_t fAccumulatorSum
Accumulator sum.
Double_t fNumberOfSamples
Number of accumulated samples.
void QwSIS3320_Accumulator::ConstructBranchAndVector ( TTree *  tree,
TString &  prefix,
std::vector< Double_t > &  values 
)

Definition at line 245 of file QwSIS3320_Accumulator.cc.

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

246 {
247  if (IsNameEmpty()) {
248  // This accumulator is not used, so skip setting up the tree.
249  } else {
250  TString basename = prefix + GetElementName();
251  fTreeArrayIndex = values.size();
252 
253  values.push_back(0.0);
254  TString list = "hw_sum/D";
255  values.push_back(0.0);
256  list += ":num_samples/D";
257 
258  fTreeArrayNumEntries = values.size() - fTreeArrayIndex;
259  tree->Branch(basename, &(values[fTreeArrayIndex]), list);
260  }
261 }
Bool_t IsNameEmpty() const
Is the name of this element empty?
size_t fTreeArrayIndex
Time-based accumulator limits.
virtual const TString & GetElementName() const
Get the name of this element.
size_t fTreeArrayNumEntries
Index of this data element in tree.

+ Here is the call graph for this function:

void QwSIS3320_Accumulator::ConstructHistograms ( TDirectory *  folder,
TString &  prefix 
)
virtual

Construct the histograms for this data element.

Implements VQwDataElement.

Definition at line 284 of file QwSIS3320_Accumulator.cc.

References QwHistogramHelper::Construct1DHist(), MQwHistograms::fHistograms, VQwDataElement::GetElementName(), gQwHists, and VQwDataElement::IsNameEmpty().

285 {
286  // If we have defined a subdirectory in the ROOT file, then change into it.
287  if (folder != NULL) folder->cd();
288 
289  if (IsNameEmpty()) {
290 
291  // This channel is not used, so skip filling the histograms.
292 
293  } else {
294 
295  TString basename = prefix + GetElementName();
296 
297  fHistograms.resize(3, NULL);
298  size_t index = 0;
299  fHistograms[index++] = gQwHists.Construct1DHist(basename + Form("_nevents"));
300  fHistograms[index++] = gQwHists.Construct1DHist(basename + Form("_sum"));
301  fHistograms[index++] = gQwHists.Construct1DHist(basename + Form("_avg"));
302 
303  }
304 }
Bool_t IsNameEmpty() const
Is the name of this element empty?
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
Definition: MQwHistograms.h:46
QwHistogramHelper gQwHists
Globally defined instance of the QwHistogramHelper class.
virtual const TString & GetElementName() const
Get the name of this element.
TH1F * Construct1DHist(const TString &inputfile, const TString &name_title)

+ Here is the call graph for this function:

void QwSIS3320_Accumulator::Difference ( const QwSIS3320_Accumulator value1,
const QwSIS3320_Accumulator value2 
)

Definition at line 228 of file QwSIS3320_Accumulator.cc.

229 {
230  *this = value1;
231  *this -= value2;
232 }
void QwSIS3320_Accumulator::FillHistograms ( )
virtual

Fill the histograms for this data element.

Implements VQwDataElement.

Definition at line 306 of file QwSIS3320_Accumulator.cc.

References MQwHistograms::fHistograms, GetAccumulatorAvg(), GetAccumulatorSum(), GetNumberOfSamples(), and VQwDataElement::IsNameEmpty().

307 {
308  Int_t index = 0;
309  if (IsNameEmpty()) {
310 
311  // This channel is not used, so skip creating the histograms.
312 
313  } else {
314 
315  if (fHistograms[index] != NULL)
316  fHistograms[index++]->Fill(GetNumberOfSamples());
317  if (fHistograms[index] != NULL)
318  fHistograms[index+1]->Fill(GetAccumulatorSum());
319  if (fHistograms[index] != NULL)
320  fHistograms[index+1]->Fill(GetAccumulatorAvg());
321 
322  }
323 }
Double_t GetAccumulatorAvg() const
Bool_t IsNameEmpty() const
Is the name of this element empty?
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
Definition: MQwHistograms.h:46
Double_t GetAccumulatorSum() const
Double_t GetNumberOfSamples() const

+ Here is the call graph for this function:

void QwSIS3320_Accumulator::FillTreeVector ( std::vector< Double_t > &  values) const

Definition at line 263 of file QwSIS3320_Accumulator.cc.

References QwLog::endl(), fTreeArrayIndex, fTreeArrayNumEntries, GetAccumulatorSum(), GetNumberOfSamples(), VQwDataElement::IsNameEmpty(), and QwError.

264 {
265  if (IsNameEmpty()) {
266  // This accumulator is not used, so skip filling the tree vector.
267  } else if (fTreeArrayNumEntries <= 0) {
268  QwError << "QwSIS3320_Accumulator::FillTreeVector: fTreeArrayNumEntries == "
270  } else if (values.size() < fTreeArrayIndex + fTreeArrayNumEntries) {
271  QwError << "QwSIS3320_Accumulator::FillTreeVector: values.size() == "
272  << values.size()
273  << "; fTreeArrayIndex + fTreeArrayNumEntries == "
275  << QwLog::endl;
276  } else {
277  size_t index = fTreeArrayIndex;
278  values[index++] = GetAccumulatorSum();
279  values[index++] = GetNumberOfSamples();
280  }
281 }
Bool_t IsNameEmpty() const
Is the name of this element empty?
size_t fTreeArrayIndex
Time-based accumulator limits.
Double_t GetAccumulatorSum() const
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
Double_t GetNumberOfSamples() const
size_t fTreeArrayNumEntries
Index of this data element in tree.
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40

+ Here is the call graph for this function:

Double_t QwSIS3320_Accumulator::GetAccumulatorAvg ( ) const
inline

Definition at line 53 of file QwSIS3320_Accumulator.h.

References fAccumulatorSum, and fNumberOfSamples.

Referenced by FillHistograms().

53  {
54  if (fAccumulatorSum > 0)
56  else return 0.0;
57  };
Double_t fAccumulatorSum
Accumulator sum.
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the caller graph for this function:

Double_t QwSIS3320_Accumulator::GetAccumulatorSum ( ) const
inline

Definition at line 52 of file QwSIS3320_Accumulator.h.

References fAccumulatorSum.

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

52 { return fAccumulatorSum; };
Double_t fAccumulatorSum
Accumulator sum.

+ Here is the caller graph for this function:

Int_t QwSIS3320_Accumulator::GetMaxTime ( ) const
inline

Definition at line 49 of file QwSIS3320_Accumulator.h.

References fMaxTime.

49 { return fMaxTime; };
Int_t QwSIS3320_Accumulator::GetMaxValue ( ) const
inline

Definition at line 47 of file QwSIS3320_Accumulator.h.

References fMaxValue.

47 { return fMaxValue; };
Int_t QwSIS3320_Accumulator::GetMinTime ( ) const
inline

Definition at line 48 of file QwSIS3320_Accumulator.h.

References fMinTime.

48 { return fMinTime; };
int fMinTime
Value-based accumulator limits.
Int_t QwSIS3320_Accumulator::GetMinValue ( ) const
inline

Definition at line 46 of file QwSIS3320_Accumulator.h.

References fMinValue.

46 { return fMinValue; };
Double_t QwSIS3320_Accumulator::GetNumberOfSamples ( ) const
inline

Definition at line 51 of file QwSIS3320_Accumulator.h.

References fNumberOfSamples.

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

51 { return fNumberOfSamples; };
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the caller graph for this function:

const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator* ( const Double_t &  value) const

Multiplication of factor to accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 81 of file QwSIS3320_Accumulator.cc.

82 {
83  QwSIS3320_Accumulator result = *this;
84  result *= value;
85  return result;
86 }
SIS3320 sampling ADC accumulator.
QwSIS3320_Accumulator & QwSIS3320_Accumulator::operator*= ( const Double_t &  value)

Multiplication assignment of factor to accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 106 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, and VQwDataElement::IsNameEmpty().

107 {
108  if (!IsNameEmpty()) {
109  this->fAccumulatorSum *= value;
110  }
111  return *this;
112 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?

+ Here is the call graph for this function:

const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator+ ( const Double_t &  value) const

Addition of offset to accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 57 of file QwSIS3320_Accumulator.cc.

58 {
59  QwSIS3320_Accumulator result = *this;
60  result += value;
61  return result;
62 }
SIS3320 sampling ADC accumulator.
const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator+ ( const QwSIS3320_Accumulator value) const

Addition of accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 159 of file QwSIS3320_Accumulator.cc.

160 {
161  QwSIS3320_Accumulator result = *this;
162  result += value;
163  return result;
164 }
SIS3320 sampling ADC accumulator.
QwSIS3320_Accumulator & QwSIS3320_Accumulator::operator+= ( const Double_t &  value)

Addition assignment of offset to accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 133 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, and VQwDataElement::IsNameEmpty().

134 {
135  if (!IsNameEmpty()) {
136  this->fAccumulatorSum += value;
137  }
138  return *this;
139 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?

+ Here is the call graph for this function:

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

Addition assignment of accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 197 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, fNumberOfSamples, and VQwDataElement::IsNameEmpty().

198 {
199  if (!IsNameEmpty()) {
200  this->fAccumulatorSum += value.fAccumulatorSum;
201  this->fNumberOfSamples += value.fNumberOfSamples;
202  }
203  return *this;
204 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the call graph for this function:

const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator- ( const Double_t &  value) const

Subtraction of offset from accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 69 of file QwSIS3320_Accumulator.cc.

70 {
71  QwSIS3320_Accumulator result = *this;
72  result -= value;
73  return result;
74 }
SIS3320 sampling ADC accumulator.
const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator- ( const QwSIS3320_Accumulator value) const

Subtraction of accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 171 of file QwSIS3320_Accumulator.cc.

172 {
173  QwSIS3320_Accumulator result = *this;
174  result -= value;
175  return result;
176 }
SIS3320 sampling ADC accumulator.
QwSIS3320_Accumulator & QwSIS3320_Accumulator::operator-= ( const Double_t &  value)

Subtraction assignment of offset from accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 146 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, and VQwDataElement::IsNameEmpty().

147 {
148  if (!IsNameEmpty()) {
149  this->fAccumulatorSum -= value;
150  }
151  return *this;
152 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?

+ Here is the call graph for this function:

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

Subtraction assignment of accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 211 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, fNumberOfSamples, and VQwDataElement::IsNameEmpty().

212 {
213  if (!IsNameEmpty()) {
214  this->fAccumulatorSum -= value.fAccumulatorSum;
215  this->fNumberOfSamples -= value.fNumberOfSamples;
216  }
217  return *this;
218 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the call graph for this function:

const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator/ ( const Double_t &  value) const

Division of factor from accumulated data (not division-by-zero safe)

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 93 of file QwSIS3320_Accumulator.cc.

94 {
95  QwSIS3320_Accumulator result = *this;
96  if (value != 0)
97  result /= value;
98  return result;
99 }
SIS3320 sampling ADC accumulator.
QwSIS3320_Accumulator & QwSIS3320_Accumulator::operator/= ( const Double_t &  value)

Division assignment of factor from accumulated data (not division-by-zero safe)

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 119 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, and VQwDataElement::IsNameEmpty().

120 {
121  if (!IsNameEmpty()) {
122  if (value != 0.0)
123  this->fAccumulatorSum /= value;
124  }
125  return *this;
126 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?

+ Here is the call graph for this function:

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

Assignment of accumulated data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 183 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, fNumberOfSamples, and VQwDataElement::IsNameEmpty().

184 {
185  if (!IsNameEmpty()) {
186  this->fAccumulatorSum = value.fAccumulatorSum;
187  this->fNumberOfSamples = value.fNumberOfSamples;
188  }
189  return *this;
190 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the call graph for this function:

Int_t QwSIS3320_Accumulator::ProcessEvBuffer ( UInt_t *  buffer,
UInt_t  num_words_left,
UInt_t  subelement = 0 
)
virtual

Process the CODA event buffer for this element.

Implements VQwDataElement.

Reimplemented in QwSIS3320_LogicalAccumulator.

Definition at line 27 of file QwSIS3320_Accumulator.cc.

References QwLog::endl(), fAccumulatorSum, VQwDataElement::fNumberOfDataWords, fNumberOfSamples, and QwWarning.

28 {
29  UInt_t words_read = 0;
30  Long64_t hw_sum;
31  Long64_t num_samples;
32  if (num_words_left >= fNumberOfDataWords) {
33 
34  // Read all words
35  num_samples = buffer[0];
36  hw_sum = buffer[2]; // first assign Int_t to Long_t, so that we
37  hw_sum <<= sizeof(Int_t); // can shift it into the higher registers
38  hw_sum += buffer[1];
39  words_read = fNumberOfDataWords;
40 
41  fNumberOfSamples = num_samples;
42  fAccumulatorSum = hw_sum;
43 
44  } else {
45  QwWarning << "QwSIS3320_Accumulator::ProcessEvBuffer: Not enough words!" << QwLog::endl;
46  }
47 
48  return words_read;
49 }
Double_t fAccumulatorSum
Accumulator sum.
Double_t fNumberOfSamples
Number of accumulated samples.
UInt_t fNumberOfDataWords
Number of raw data words in this data element.
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
#define QwWarning
Predefined log drain for warnings.
Definition: QwLog.h:45

+ Here is the call graph for this function:

void QwSIS3320_Accumulator::ProcessEvent ( )
inline

Definition at line 61 of file QwSIS3320_Accumulator.h.

61 { };
void QwSIS3320_Accumulator::Ratio ( const QwSIS3320_Accumulator numer,
const QwSIS3320_Accumulator denom 
)

Definition at line 235 of file QwSIS3320_Accumulator.cc.

References fAccumulatorSum, fNumberOfSamples, and VQwDataElement::IsNameEmpty().

236 {
237  if (!IsNameEmpty()) {
238  if (denom.fAccumulatorSum != 0.0)
241  }
242 }
Double_t fAccumulatorSum
Accumulator sum.
Bool_t IsNameEmpty() const
Is the name of this element empty?
Double_t fNumberOfSamples
Number of accumulated samples.

+ Here is the call graph for this function:

void QwSIS3320_Accumulator::Sum ( const QwSIS3320_Accumulator value1,
const QwSIS3320_Accumulator value2 
)

Definition at line 221 of file QwSIS3320_Accumulator.cc.

222 {
223  *this = value1;
224  *this += value2;
225 }

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  stream,
const QwSIS3320_Accumulator a 
)
friend

Definition at line 112 of file QwSIS3320_Accumulator.h.

112  {
113  stream << a.GetAccumulatorSum() << " (" << a.GetNumberOfSamples() << ")";
114  return stream;
115 }
Double_t GetAccumulatorSum() const
Double_t GetNumberOfSamples() const

Field Documentation

Double_t QwSIS3320_Accumulator::fAccumulatorAvg
protected

Accumulator average.

Definition at line 95 of file QwSIS3320_Accumulator.h.

Double_t QwSIS3320_Accumulator::fAccumulatorSum
protected
int QwSIS3320_Accumulator::fMaxTime
private

Definition at line 100 of file QwSIS3320_Accumulator.h.

Referenced by GetMaxTime(), and QwSIS3320_Accumulator().

int QwSIS3320_Accumulator::fMaxValue
private

Definition at line 99 of file QwSIS3320_Accumulator.h.

Referenced by GetMaxValue(), and QwSIS3320_Accumulator().

int QwSIS3320_Accumulator::fMinTime
private

Value-based accumulator limits.

Definition at line 100 of file QwSIS3320_Accumulator.h.

Referenced by GetMinTime(), and QwSIS3320_Accumulator().

int QwSIS3320_Accumulator::fMinValue
private

Definition at line 99 of file QwSIS3320_Accumulator.h.

Referenced by GetMinValue(), and QwSIS3320_Accumulator().

Double_t QwSIS3320_Accumulator::fNumberOfSamples
protected

Number of accumulated samples.

Definition at line 93 of file QwSIS3320_Accumulator.h.

Referenced by ClearEventData(), GetAccumulatorAvg(), GetNumberOfSamples(), operator+=(), operator-=(), operator=(), ProcessEvBuffer(), and Ratio().

size_t QwSIS3320_Accumulator::fTreeArrayIndex
private

Time-based accumulator limits.

Definition at line 103 of file QwSIS3320_Accumulator.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

size_t QwSIS3320_Accumulator::fTreeArrayNumEntries
private

Index of this data element in tree.

Definition at line 104 of file QwSIS3320_Accumulator.h.

Referenced by ConstructBranchAndVector(), and FillTreeVector().

const unsigned int QwSIS3320_Accumulator::INDEX_NUM = 0
staticprivate

Number of entries from this data element.

Definition at line 106 of file QwSIS3320_Accumulator.h.

const unsigned int QwSIS3320_Accumulator::INDEX_SUM = 1
staticprivate

Definition at line 107 of file QwSIS3320_Accumulator.h.


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