QwAnalysis
QwSIS3320_Samples Class Reference

SIS3320 sampling ADC samples. More...

#include <QwSIS3320_Samples.h>

Inherits TObject.

Public Member Functions

 QwSIS3320_Samples (UInt_t nsamples=256)
 
virtual ~QwSIS3320_Samples ()
 
size_t GetMinIndex () const
 
size_t GetMaxIndex () const
 
QwSIS3320_Type GetMinSample () const
 
QwSIS3320_Type GetMaxSample () const
 
QwSIS3320_Type GetSum () const
 
QwSIS3320_Type GetSample (size_t i) const
 
QwSIS3320_Type GetPedestal () const
 
QwSIS3320_Type GetSumInTimeWindow (const UInt_t start, const UInt_t stop) const
 
UInt_t GetNumberOfDataWords () const
 
void SetNumberOfDataWords (const UInt_t &numwords)
 
UInt_t GetNumberOfSamples () const
 
void SetNumberOfSamples (const UInt_t nsamples)
 
UInt_t GetSamplePointer () const
 
void SetSamplePointer (const UInt_t samplepointer)
 
UInt_t GetSamplesPerWord () const
 
void SetSamplesPerWord (const UInt_t nsamples)
 
TGraph * GetGraph () const
 
void UpdateGraph ()
 
void ClearEventData ()
 
Int_t ProcessEvBuffer (UInt_t *buffer, UInt_t num_words_left, UInt_t subelement=0)
 
QwSIS3320_Samplesoperator/= (const Double_t &value)
 
QwSIS3320_Samplesoperator*= (const Double_t &value)
 
QwSIS3320_Samplesoperator+= (const Double_t &value)
 
QwSIS3320_Samplesoperator-= (const Double_t &value)
 
const QwSIS3320_Samples operator/ (const Double_t &value) const
 
const QwSIS3320_Samples operator* (const Double_t &value) const
 
const QwSIS3320_Samples operator+ (const Double_t &value) const
 
const QwSIS3320_Samples operator- (const Double_t &value) const
 
QwSIS3320_Samplesoperator= (const QwSIS3320_Samples &value)
 
QwSIS3320_Samplesoperator+= (const QwSIS3320_Samples &value)
 
QwSIS3320_Samplesoperator-= (const QwSIS3320_Samples &value)
 
const QwSIS3320_Samples operator+ (const QwSIS3320_Samples &value) const
 
const QwSIS3320_Samples operator- (const QwSIS3320_Samples &value) const
 

Private Member Functions

std::pair< size_t, QwSIS3320_TypeGetMin () const
 
std::pair< size_t, QwSIS3320_TypeGetMax () const
 
 ClassDef (QwSIS3320_Samples, 1)
 

Private Attributes

UInt_t fSamplesPerWord
 Number of 12-bit sample values per data word. More...
 
UInt_t fNumberOfDataWords
 Number of data words in this data element. More...
 
UInt_t fSamplePointer
 Sample position in buffer. More...
 
std::vector< QwSIS3320_TypefSamples
 Samples values. More...
 
TGraph * fGraph
 Graph of samples. More...
 
size_t fTreeArrayIndex
 Index of this data element in tree. More...
 
size_t fTreeArrayNumEntries
 Number of entries from this data element. More...
 

Static Private Attributes

static std::vector
< QwSIS3320_Type
fIndex
 Samples index. More...
 

Friends

std::ostream & operator<< (std::ostream &stream, const QwSIS3320_Samples &s)
 

Detailed Description

SIS3320 sampling ADC samples.

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

The QwSIS3320_Samples should allow convenient access to the sampling 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 36 of file QwSIS3320_Samples.h.

Constructor & Destructor Documentation

QwSIS3320_Samples::QwSIS3320_Samples ( UInt_t  nsamples = 256)
inline

Definition at line 40 of file QwSIS3320_Samples.h.

References fGraph, SetNumberOfSamples(), SetSamplePointer(), and SetSamplesPerWord().

40  {
41  fGraph = 0;
44  SetNumberOfSamples(nsamples);
45  };
TGraph * fGraph
Graph of samples.
void SetNumberOfSamples(const UInt_t nsamples)
void SetSamplesPerWord(const UInt_t nsamples)
void SetSamplePointer(const UInt_t samplepointer)

+ Here is the call graph for this function:

virtual QwSIS3320_Samples::~QwSIS3320_Samples ( )
inlinevirtual

Definition at line 46 of file QwSIS3320_Samples.h.

References fGraph.

46  {
47  if (fGraph) delete fGraph;
48  };
TGraph * fGraph
Graph of samples.

Member Function Documentation

QwSIS3320_Samples::ClassDef ( QwSIS3320_Samples  ,
 
)
private
void QwSIS3320_Samples::ClearEventData ( )
inline

Definition at line 93 of file QwSIS3320_Samples.h.

References fSamples.

Referenced by QwSIS3320_Channel::ClearEventData(), and QwSIS3320_Channel::InitializeChannel().

93 { fSamples.clear(); };
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

TGraph* QwSIS3320_Samples::GetGraph ( ) const
inline

Definition at line 87 of file QwSIS3320_Samples.h.

References fGraph.

87 { return fGraph; };
TGraph * fGraph
Graph of samples.
std::pair< size_t, QwSIS3320_Type > QwSIS3320_Samples::GetMax ( ) const
private

Definition at line 61 of file QwSIS3320_Samples.cc.

References fSamples.

Referenced by GetMaxIndex(), and GetMaxSample().

62 {
63  std::vector<QwSIS3320_Type>::const_iterator max = std::max_element(fSamples.begin(), fSamples.end());
64  return std::pair<size_t,QwSIS3320_Type>(max - fSamples.begin(),*max);
65 }
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

size_t QwSIS3320_Samples::GetMaxIndex ( ) const
inline

Definition at line 51 of file QwSIS3320_Samples.h.

References GetMax().

51 { return GetMax().first; };
std::pair< size_t, QwSIS3320_Type > GetMax() const

+ Here is the call graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetMaxSample ( ) const
inline

Definition at line 53 of file QwSIS3320_Samples.h.

References GetMax().

53 { return GetMax().second; };
std::pair< size_t, QwSIS3320_Type > GetMax() const

+ Here is the call graph for this function:

std::pair< size_t, QwSIS3320_Type > QwSIS3320_Samples::GetMin ( ) const
private

Definition at line 67 of file QwSIS3320_Samples.cc.

References fSamples, and Qw::min.

Referenced by GetMinIndex(), and GetMinSample().

68 {
69  std::vector<QwSIS3320_Type>::const_iterator min = std::min_element(fSamples.begin(), fSamples.end());
70  return std::pair<size_t,QwSIS3320_Type>(min - fSamples.begin(),*min);
71 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
static const double min
Definition: QwUnits.h:76

+ Here is the caller graph for this function:

size_t QwSIS3320_Samples::GetMinIndex ( ) const
inline

Definition at line 50 of file QwSIS3320_Samples.h.

References GetMin().

50 { return GetMin().first; };
std::pair< size_t, QwSIS3320_Type > GetMin() const

+ Here is the call graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetMinSample ( ) const
inline

Definition at line 52 of file QwSIS3320_Samples.h.

References GetMin().

52 { return GetMin().second; };
std::pair< size_t, QwSIS3320_Type > GetMin() const

+ Here is the call graph for this function:

UInt_t QwSIS3320_Samples::GetNumberOfDataWords ( ) const
inline

Definition at line 60 of file QwSIS3320_Samples.h.

References fNumberOfDataWords.

60 { return fNumberOfDataWords; };
UInt_t fNumberOfDataWords
Number of data words in this data element.
UInt_t QwSIS3320_Samples::GetNumberOfSamples ( ) const
inline

Definition at line 65 of file QwSIS3320_Samples.h.

References fSamples.

Referenced by operator<<(), SetNumberOfSamples(), and SetSamplesPerWord().

65 { return fSamples.size(); };
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetPedestal ( ) const
inline

Definition at line 57 of file QwSIS3320_Samples.h.

References GetSample().

57 { return GetSample(0); };
QwSIS3320_Type GetSample(size_t i) const

+ Here is the call graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetSample ( size_t  i) const
inline

Definition at line 56 of file QwSIS3320_Samples.h.

References fSamples.

Referenced by QwSIS3320_Channel::FillHistograms(), GetPedestal(), and operator<<().

56 { return fSamples.at(i); };
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

UInt_t QwSIS3320_Samples::GetSamplePointer ( ) const
inline

Definition at line 75 of file QwSIS3320_Samples.h.

References fSamplePointer.

75 { return fSamplePointer; };
UInt_t fSamplePointer
Sample position in buffer.
UInt_t QwSIS3320_Samples::GetSamplesPerWord ( ) const
inline

Definition at line 80 of file QwSIS3320_Samples.h.

References fSamplesPerWord.

Referenced by SetNumberOfSamples(), and SetSamplesPerWord().

80 { return fSamplesPerWord; };
UInt_t fSamplesPerWord
Number of 12-bit sample values per data word.

+ Here is the caller graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetSum ( ) const

Definition at line 73 of file QwSIS3320_Samples.cc.

References fSamples.

Referenced by QwSIS3320_Channel::FillHistograms(), and QwSIS3320_Channel::PrintInfo().

74 {
75  return std::accumulate(fSamples.begin(), fSamples.end(), 0.0);
76 }
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

QwSIS3320_Type QwSIS3320_Samples::GetSumInTimeWindow ( const UInt_t  start,
const UInt_t  stop 
) const

Definition at line 79 of file QwSIS3320_Samples.cc.

References fSamples.

Referenced by QwSIS3320_Channel::PrintInfo().

80 {
81  if (start >= fSamples.size() || stop >= fSamples.size()) return 0;
82  return std::accumulate(&fSamples.at(start), &fSamples.at(stop), 0.0);
83 }
std::vector< QwSIS3320_Type > fSamples
Samples values.

+ Here is the caller graph for this function:

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

Multiplication of factor to sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 124 of file QwSIS3320_Samples.cc.

125 {
126  QwSIS3320_Samples result = *this;
127  result *= value;
128  return result;
129 }
SIS3320 sampling ADC samples.
QwSIS3320_Samples & QwSIS3320_Samples::operator*= ( const Double_t &  value)

Multiplication assignment of factor to sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 148 of file QwSIS3320_Samples.cc.

References fSamples.

149 {
150  for (size_t i = 0; i < fSamples.size(); i++)
151  this->fSamples.at(i) *= value;
152  return *this;
153 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
const QwSIS3320_Samples QwSIS3320_Samples::operator+ ( const Double_t &  value) const

Addition of offset to sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 100 of file QwSIS3320_Samples.cc.

101 {
102  QwSIS3320_Samples result = *this;
103  result += value;
104  return result;
105 }
SIS3320 sampling ADC samples.
const QwSIS3320_Samples QwSIS3320_Samples::operator+ ( const QwSIS3320_Samples value) const

Addition of sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 199 of file QwSIS3320_Samples.cc.

200 {
201  QwSIS3320_Samples result = *this;
202  result += value;
203  return result;
204 }
SIS3320 sampling ADC samples.
QwSIS3320_Samples & QwSIS3320_Samples::operator+= ( const Double_t &  value)

Addition assignment of offset to sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 175 of file QwSIS3320_Samples.cc.

References fSamples.

176 {
177  for (size_t i = 0; i < fSamples.size(); i++)
178  this->fSamples.at(i) += value;
179  return *this;
180 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
QwSIS3320_Samples & QwSIS3320_Samples::operator+= ( const QwSIS3320_Samples value)

Addition assignment of sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 237 of file QwSIS3320_Samples.cc.

References fSamples.

238 {
239  for (size_t i = 0; i < fSamples.size(); i++)
240  this->fSamples.at(i) += value.fSamples.at(i);
241  return *this;
242 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
const QwSIS3320_Samples QwSIS3320_Samples::operator- ( const Double_t &  value) const

Subtraction of offset from sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 112 of file QwSIS3320_Samples.cc.

113 {
114  QwSIS3320_Samples result = *this;
115  result -= value;
116  return result;
117 }
SIS3320 sampling ADC samples.
const QwSIS3320_Samples QwSIS3320_Samples::operator- ( const QwSIS3320_Samples value) const

Subtraction of sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 211 of file QwSIS3320_Samples.cc.

212 {
213  QwSIS3320_Samples result = *this;
214  result -= value;
215  return result;
216 }
SIS3320 sampling ADC samples.
QwSIS3320_Samples & QwSIS3320_Samples::operator-= ( const Double_t &  value)

Subtraction assignment of offset from sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 187 of file QwSIS3320_Samples.cc.

References fSamples.

188 {
189  for (size_t i = 0; i < fSamples.size(); i++)
190  this->fSamples.at(i) -= value;
191  return *this;
192 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
QwSIS3320_Samples & QwSIS3320_Samples::operator-= ( const QwSIS3320_Samples value)

Subtraction assignment of sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 249 of file QwSIS3320_Samples.cc.

References fSamples.

250 {
251  for (size_t i = 0; i < fSamples.size(); i++)
252  this->fSamples.at(i) -= value.fSamples.at(i);
253  return *this;
254 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
const QwSIS3320_Samples QwSIS3320_Samples::operator/ ( const Double_t &  value) const

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

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 136 of file QwSIS3320_Samples.cc.

137 {
138  QwSIS3320_Samples result = *this;
139  result /= value;
140  return result;
141 }
SIS3320 sampling ADC samples.
QwSIS3320_Samples & QwSIS3320_Samples::operator/= ( const Double_t &  value)

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

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 160 of file QwSIS3320_Samples.cc.

References fSamples.

161 {
162  for (size_t i = 0; i < fSamples.size(); i++)
163  if (value != 0.0)
164  this->fSamples.at(i) /= value;
165  else
166  this->fSamples.at(i) = 0.0;
167  return *this;
168 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
QwSIS3320_Samples & QwSIS3320_Samples::operator= ( const QwSIS3320_Samples value)

Assignment of sampled data

Parameters
valueRight-hand side
Returns
Left-hand side

Definition at line 223 of file QwSIS3320_Samples.cc.

References fNumberOfDataWords, fSamplePointer, fSamples, and fSamplesPerWord.

224 {
228  fSamples = value.fSamples;
229  return *this;
230 }
UInt_t fSamplePointer
Sample position in buffer.
std::vector< QwSIS3320_Type > fSamples
Samples values.
UInt_t fNumberOfDataWords
Number of data words in this data element.
UInt_t fSamplesPerWord
Number of 12-bit sample values per data word.
Int_t QwSIS3320_Samples::ProcessEvBuffer ( UInt_t *  buffer,
UInt_t  num_words_left,
UInt_t  subelement = 0 
)

Definition at line 29 of file QwSIS3320_Samples.cc.

References QwLog::endl(), fNumberOfDataWords, fSamples, fSamplesPerWord, and QwError.

30 {
31  UInt_t words_read = 0;
32  UInt_t index = 0;
33  if (num_words_left >= fNumberOfDataWords) {
34  // Read all words
35  for (size_t i = 0; i < fNumberOfDataWords; i++) {
36  // Shift bits as necessary
37  switch (fSamplesPerWord) {
38  case 1:
39  fSamples[index++] = buffer[i];
40  break;
41  case 2:
42  fSamples[index++] = buffer[i] & 0xFFFF; // lowest 16 bits
43  fSamples[index++] = (buffer[i] >> 16) & 0xFFFF; // highest 16 bits
44  break;
45  default:
46  QwError << "QwSIS3320_Samples: Illegal number of samples per word!" << QwLog::endl;
47  words_read = 0;
48  return words_read;
49  }
50  }
51  words_read = fNumberOfDataWords;
52 
53  } else {
54  QwError << "QwSIS3320_Samples: Not enough words while processing buffer!" << QwLog::endl;
55  }
56 
57  return words_read;
58 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
UInt_t fNumberOfDataWords
Number of data words in this data element.
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40
UInt_t fSamplesPerWord
Number of 12-bit sample values per data word.

+ Here is the call graph for this function:

void QwSIS3320_Samples::SetNumberOfDataWords ( const UInt_t &  numwords)
inline

Definition at line 61 of file QwSIS3320_Samples.h.

References fNumberOfDataWords.

Referenced by SetNumberOfSamples(), and SetSamplesPerWord().

61  {
62  fNumberOfDataWords = numwords;
63  };
UInt_t fNumberOfDataWords
Number of data words in this data element.

+ Here is the caller graph for this function:

void QwSIS3320_Samples::SetNumberOfSamples ( const UInt_t  nsamples)
inline

Definition at line 66 of file QwSIS3320_Samples.h.

References fIndex, fSamples, GetNumberOfSamples(), GetSamplesPerWord(), and SetNumberOfDataWords().

Referenced by QwSIS3320_Channel::ProcessEvent(), and QwSIS3320_Samples().

66  {
67  // Initialize index vector
68  fIndex.resize(nsamples);
69  for (size_t i = 0; i < fIndex.size(); i++) fIndex[i] = i;
70  // Initialize sample vector
71  fSamples.resize(nsamples);
73  };
UInt_t GetSamplesPerWord() const
UInt_t GetNumberOfSamples() const
std::vector< QwSIS3320_Type > fSamples
Samples values.
static std::vector< QwSIS3320_Type > fIndex
Samples index.
void SetNumberOfDataWords(const UInt_t &numwords)

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwSIS3320_Samples::SetSamplePointer ( const UInt_t  samplepointer)
inline

Definition at line 76 of file QwSIS3320_Samples.h.

References fSamplePointer.

Referenced by QwSIS3320_Samples().

76  {
77  fSamplePointer = samplepointer;
78  };
UInt_t fSamplePointer
Sample position in buffer.

+ Here is the caller graph for this function:

void QwSIS3320_Samples::SetSamplesPerWord ( const UInt_t  nsamples)
inline

Definition at line 81 of file QwSIS3320_Samples.h.

References fSamplesPerWord, GetNumberOfSamples(), GetSamplesPerWord(), and SetNumberOfDataWords().

Referenced by QwSIS3320_Samples().

81  {
82  fSamplesPerWord = nsamples;
84  };
UInt_t GetSamplesPerWord() const
UInt_t GetNumberOfSamples() const
void SetNumberOfDataWords(const UInt_t &numwords)
UInt_t fSamplesPerWord
Number of 12-bit sample values per data word.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwSIS3320_Samples::UpdateGraph ( )

Definition at line 86 of file QwSIS3320_Samples.cc.

References fGraph, fIndex, and fSamples.

87 {
88  const size_t n = fSamples.size();
89  const QwSIS3320_Type* vx = &(*fIndex.begin());
90  const QwSIS3320_Type* vy = &(*fSamples.begin());
91  fGraph = new TGraph(n, vx, vy);
92 }
std::vector< QwSIS3320_Type > fSamples
Samples values.
TGraph * fGraph
Graph of samples.
static std::vector< QwSIS3320_Type > fIndex
Samples index.
Double_t QwSIS3320_Type

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  stream,
const QwSIS3320_Samples s 
)
friend

Definition at line 144 of file QwSIS3320_Samples.h.

145 {
146  for (size_t i = 0; i < s.GetNumberOfSamples(); i++)
147  stream << s.GetSample(i) << " ";
148  return stream;
149 }
UInt_t GetNumberOfSamples() const
QwSIS3320_Type GetSample(size_t i) const

Field Documentation

TGraph* QwSIS3320_Samples::fGraph
private

Graph of samples.

Definition at line 134 of file QwSIS3320_Samples.h.

Referenced by GetGraph(), QwSIS3320_Samples(), UpdateGraph(), and ~QwSIS3320_Samples().

std::vector< QwSIS3320_Type > QwSIS3320_Samples::fIndex
staticprivate

Samples index.

Definition at line 130 of file QwSIS3320_Samples.h.

Referenced by SetNumberOfSamples(), and UpdateGraph().

UInt_t QwSIS3320_Samples::fNumberOfDataWords
private

Number of data words in this data element.

Definition at line 125 of file QwSIS3320_Samples.h.

Referenced by GetNumberOfDataWords(), operator=(), ProcessEvBuffer(), and SetNumberOfDataWords().

UInt_t QwSIS3320_Samples::fSamplePointer
private

Sample position in buffer.

Definition at line 127 of file QwSIS3320_Samples.h.

Referenced by GetSamplePointer(), operator=(), and SetSamplePointer().

UInt_t QwSIS3320_Samples::fSamplesPerWord
private

Number of 12-bit sample values per data word.

Definition at line 123 of file QwSIS3320_Samples.h.

Referenced by GetSamplesPerWord(), operator=(), ProcessEvBuffer(), and SetSamplesPerWord().

size_t QwSIS3320_Samples::fTreeArrayIndex
private

Index of this data element in tree.

Definition at line 137 of file QwSIS3320_Samples.h.

size_t QwSIS3320_Samples::fTreeArrayNumEntries
private

Number of entries from this data element.

Definition at line 138 of file QwSIS3320_Samples.h.


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