QwAnalysis
QwSIS3320_Samples.cc
Go to the documentation of this file.
1 /**
2  * \file QwSIS3320_Samples.cc
3  *
4  * \brief Implementation of the SIS3320 sampling ADC samples
5  *
6  * \author W. Deconinck
7  * \date 2009-09-04 18:06:23
8  * \ingroup QwCompton
9  *
10  * The QwSIS3320_Samples should allow convenient access to the sampling data
11  * collected with the SIS3320 for the Compton photon detector. This class
12  * implements its own sum, difference, and ratio methods inherited from the
13  * general VQwDataElement.
14  *
15  */
16 
17 #include "QwSIS3320_Samples.h"
19 
20 // System headers
21 #include <numeric>
22 #include <algorithm>
23 
24 // Qweak headers
25 #include "QwLog.h"
26 
27 std::vector<QwSIS3320_Type> QwSIS3320_Samples::fIndex;
28 
29 Int_t QwSIS3320_Samples::ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t subelement)
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 }
59 
60 
61 std::pair<size_t,QwSIS3320_Type> QwSIS3320_Samples::GetMax() const
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 }
66 
67 std::pair<size_t,QwSIS3320_Type> QwSIS3320_Samples::GetMin() const
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 }
72 
74 {
75  return std::accumulate(fSamples.begin(), fSamples.end(), 0.0);
76 }
77 
78 
79 QwSIS3320_Type QwSIS3320_Samples::GetSumInTimeWindow(const UInt_t start, const UInt_t stop) const
80 {
81  if (start >= fSamples.size() || stop >= fSamples.size()) return 0;
82  return std::accumulate(&fSamples.at(start), &fSamples.at(stop), 0.0);
83 }
84 
85 
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 }
93 
94 
95 /**
96  * Addition of offset to sampled data
97  * @param value Right-hand side
98  * @return Left-hand side
99  */
100 const QwSIS3320_Samples QwSIS3320_Samples::operator+ (const Double_t &value) const
101 {
102  QwSIS3320_Samples result = *this;
103  result += value;
104  return result;
105 }
106 
107 /**
108  * Subtraction of offset from sampled data
109  * @param value Right-hand side
110  * @return Left-hand side
111  */
112 const QwSIS3320_Samples QwSIS3320_Samples::operator- (const Double_t &value) const
113 {
114  QwSIS3320_Samples result = *this;
115  result -= value;
116  return result;
117 }
118 
119 /**
120  * Multiplication of factor to sampled data
121  * @param value Right-hand side
122  * @return Left-hand side
123  */
124 const QwSIS3320_Samples QwSIS3320_Samples::operator* (const Double_t &value) const
125 {
126  QwSIS3320_Samples result = *this;
127  result *= value;
128  return result;
129 }
130 
131 /**
132  * Division of factor from sampled data (not division-by-zero safe)
133  * @param value Right-hand side
134  * @return Left-hand side
135  */
136 const QwSIS3320_Samples QwSIS3320_Samples::operator/ (const Double_t &value) const
137 {
138  QwSIS3320_Samples result = *this;
139  result /= value;
140  return result;
141 }
142 
143 /**
144  * Multiplication assignment of factor to sampled data
145  * @param value Right-hand side
146  * @return Left-hand side
147  */
149 {
150  for (size_t i = 0; i < fSamples.size(); i++)
151  this->fSamples.at(i) *= value;
152  return *this;
153 }
154 
155 /**
156  * Division assignment of factor from sampled data (not division-by-zero safe)
157  * @param value Right-hand side
158  * @return Left-hand side
159  */
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 }
169 
170 /**
171  * Addition assignment of offset to sampled data
172  * @param value Right-hand side
173  * @return Left-hand side
174  */
176 {
177  for (size_t i = 0; i < fSamples.size(); i++)
178  this->fSamples.at(i) += value;
179  return *this;
180 }
181 
182 /**
183  * Subtraction assignment of offset from sampled data
184  * @param value Right-hand side
185  * @return Left-hand side
186  */
188 {
189  for (size_t i = 0; i < fSamples.size(); i++)
190  this->fSamples.at(i) -= value;
191  return *this;
192 }
193 
194 /**
195  * Addition of sampled data
196  * @param value Right-hand side
197  * @return Left-hand side
198  */
200 {
201  QwSIS3320_Samples result = *this;
202  result += value;
203  return result;
204 }
205 
206 /**
207  * Subtraction of sampled data
208  * @param value Right-hand side
209  * @return Left-hand side
210  */
212 {
213  QwSIS3320_Samples result = *this;
214  result -= value;
215  return result;
216 }
217 
218 /**
219  * Assignment of sampled data
220  * @param value Right-hand side
221  * @return Left-hand side
222  */
224 {
228  fSamples = value.fSamples;
229  return *this;
230 }
231 
232 /**
233  * Addition assignment of sampled data
234  * @param value Right-hand side
235  * @return Left-hand side
236  */
238 {
239  for (size_t i = 0; i < fSamples.size(); i++)
240  this->fSamples.at(i) += value.fSamples.at(i);
241  return *this;
242 }
243 
244 /**
245  * Subtraction assignment of sampled data
246  * @param value Right-hand side
247  * @return Left-hand side
248  */
250 {
251  for (size_t i = 0; i < fSamples.size(); i++)
252  this->fSamples.at(i) -= value.fSamples.at(i);
253  return *this;
254 }
SIS3320 sampling ADC samples.
UInt_t fSamplePointer
Sample position in buffer.
const QwSIS3320_Samples operator-(const Double_t &value) const
QwSIS3320_Samples & operator-=(const Double_t &value)
std::vector< QwSIS3320_Type > fSamples
Samples values.
TGraph * fGraph
Graph of samples.
UInt_t fNumberOfDataWords
Number of data words in this data element.
static std::vector< QwSIS3320_Type > fIndex
Samples index.
A logfile class, based on an identical class in the Hermes analyzer.
std::pair< size_t, QwSIS3320_Type > GetMin() const
QwSIS3320_Samples & operator=(const QwSIS3320_Samples &value)
const QwSIS3320_Samples operator*(const Double_t &value) const
QwSIS3320_Type GetSumInTimeWindow(const UInt_t start, const UInt_t stop) const
static const double min
Definition: QwUnits.h:76
Double_t QwSIS3320_Type
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t subelement=0)
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
QwSIS3320_Samples & operator/=(const Double_t &value)
ClassImp(QwF1TDC)
QwSIS3320_Type GetSum() const
QwSIS3320_Samples & operator+=(const Double_t &value)
std::pair< size_t, QwSIS3320_Type > GetMax() const
const QwSIS3320_Samples operator/(const Double_t &value) const
QwSIS3320_Samples & operator*=(const Double_t &value)
const QwSIS3320_Samples operator+(const Double_t &value) const
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40
UInt_t fSamplesPerWord
Number of 12-bit sample values per data word.