QwAnalysis
QwSIS3320_Accumulator.cc
Go to the documentation of this file.
1 /**
2  * \file QwSIS3320_Accumulator.cc
3  *
4  * \brief Implementation of the SIS3320 sampling ADC accumulator
5  *
6  * \author W. Deconinck
7  * \date 2009-09-04 18:06:23
8  * \ingroup QwCompton
9  *
10  * The QwSIS3320_Accumulator should allow convenient access to the accumulator
11  * data 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_Accumulator.h"
18 
19 // Qweak headers
20 #include "QwLog.h"
21 #include "QwHistogramHelper.h"
22 
23 // Fields in accumulator data buffer
24 const unsigned int QwSIS3320_Accumulator::INDEX_NUM = 0;
25 const unsigned int QwSIS3320_Accumulator::INDEX_SUM = 1;
26 
27 Int_t QwSIS3320_Accumulator::ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t subelement)
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 }
50 
51 
52 /**
53  * Addition of offset to accumulated data
54  * @param value Right-hand side
55  * @return Left-hand side
56  */
57 const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator+ (const Double_t &value) const
58 {
59  QwSIS3320_Accumulator result = *this;
60  result += value;
61  return result;
62 }
63 
64 /**
65  * Subtraction of offset from accumulated data
66  * @param value Right-hand side
67  * @return Left-hand side
68  */
69 const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator- (const Double_t &value) const
70 {
71  QwSIS3320_Accumulator result = *this;
72  result -= value;
73  return result;
74 }
75 
76 /**
77  * Multiplication of factor to accumulated data
78  * @param value Right-hand side
79  * @return Left-hand side
80  */
81 const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator* (const Double_t &value) const
82 {
83  QwSIS3320_Accumulator result = *this;
84  result *= value;
85  return result;
86 }
87 
88 /**
89  * Division of factor from accumulated data (not division-by-zero safe)
90  * @param value Right-hand side
91  * @return Left-hand side
92  */
93 const QwSIS3320_Accumulator QwSIS3320_Accumulator::operator/ (const Double_t &value) const
94 {
95  QwSIS3320_Accumulator result = *this;
96  if (value != 0)
97  result /= value;
98  return result;
99 }
100 
101 /**
102  * Multiplication assignment of factor to accumulated data
103  * @param value Right-hand side
104  * @return Left-hand side
105  */
107 {
108  if (!IsNameEmpty()) {
109  this->fAccumulatorSum *= value;
110  }
111  return *this;
112 }
113 
114 /**
115  * Division assignment of factor from accumulated data (not division-by-zero safe)
116  * @param value Right-hand side
117  * @return Left-hand side
118  */
120 {
121  if (!IsNameEmpty()) {
122  if (value != 0.0)
123  this->fAccumulatorSum /= value;
124  }
125  return *this;
126 }
127 
128 /**
129  * Addition assignment of offset to accumulated data
130  * @param value Right-hand side
131  * @return Left-hand side
132  */
134 {
135  if (!IsNameEmpty()) {
136  this->fAccumulatorSum += value;
137  }
138  return *this;
139 }
140 
141 /**
142  * Subtraction assignment of offset from accumulated data
143  * @param value Right-hand side
144  * @return Left-hand side
145  */
147 {
148  if (!IsNameEmpty()) {
149  this->fAccumulatorSum -= value;
150  }
151  return *this;
152 }
153 
154 /**
155  * Addition of accumulated data
156  * @param value Right-hand side
157  * @return Left-hand side
158  */
160 {
161  QwSIS3320_Accumulator result = *this;
162  result += value;
163  return result;
164 }
165 
166 /**
167  * Subtraction of accumulated data
168  * @param value Right-hand side
169  * @return Left-hand side
170  */
172 {
173  QwSIS3320_Accumulator result = *this;
174  result -= value;
175  return result;
176 }
177 
178 /**
179  * Assignment of accumulated data
180  * @param value Right-hand side
181  * @return Left-hand side
182  */
184 {
185  if (!IsNameEmpty()) {
186  this->fAccumulatorSum = value.fAccumulatorSum;
187  this->fNumberOfSamples = value.fNumberOfSamples;
188  }
189  return *this;
190 }
191 
192 /**
193  * Addition assignment of accumulated data
194  * @param value Right-hand side
195  * @return Left-hand side
196  */
198 {
199  if (!IsNameEmpty()) {
200  this->fAccumulatorSum += value.fAccumulatorSum;
201  this->fNumberOfSamples += value.fNumberOfSamples;
202  }
203  return *this;
204 }
205 
206 /**
207  * Subtraction assignment of accumulated data
208  * @param value Right-hand side
209  * @return Left-hand side
210  */
212 {
213  if (!IsNameEmpty()) {
214  this->fAccumulatorSum -= value.fAccumulatorSum;
215  this->fNumberOfSamples -= value.fNumberOfSamples;
216  }
217  return *this;
218 }
219 
220 
222 {
223  *this = value1;
224  *this += value2;
225 }
226 
227 
229 {
230  *this = value1;
231  *this -= value2;
232 }
233 
234 
236 {
237  if (!IsNameEmpty()) {
238  if (denom.fAccumulatorSum != 0.0)
241  }
242 }
243 
244 
245 void QwSIS3320_Accumulator::ConstructBranchAndVector(TTree *tree, TString &prefix, std::vector<Double_t> &values)
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 }
262 
263 void QwSIS3320_Accumulator::FillTreeVector(std::vector<Double_t> &values) const
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 }
282 
283 
284 void QwSIS3320_Accumulator::ConstructHistograms(TDirectory *folder, TString &prefix)
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 }
305 
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 }
static const unsigned int INDEX_SUM
Double_t GetAccumulatorAvg() const
Double_t fAccumulatorSum
Accumulator sum.
static const unsigned int INDEX_NUM
Number of entries from this data element.
void ConstructBranchAndVector(TTree *tree, TString &prefix, std::vector< Double_t > &values)
void Difference(const QwSIS3320_Accumulator &value1, const QwSIS3320_Accumulator &value2)
Bool_t IsNameEmpty() const
Is the name of this element empty?
const QwSIS3320_Accumulator operator/(const Double_t &value) const
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
Definition: MQwHistograms.h:46
SIS3320 sampling ADC accumulator.
QwSIS3320_Accumulator & operator=(const QwSIS3320_Accumulator &value)
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t subelement=0)
Process the CODA event buffer for this element.
QwSIS3320_Accumulator & operator/=(const Double_t &value)
void FillTreeVector(std::vector< Double_t > &values) const
QwSIS3320_Accumulator & operator-=(const Double_t &value)
A logfile class, based on an identical class in the Hermes analyzer.
size_t fTreeArrayIndex
Time-based accumulator limits.
const QwSIS3320_Accumulator operator-(const Double_t &value) const
Double_t GetAccumulatorSum() const
QwHistogramHelper gQwHists
Globally defined instance of the QwHistogramHelper class.
Double_t fNumberOfSamples
Number of accumulated samples.
UInt_t fNumberOfDataWords
Number of raw data words in this data element.
QwSIS3320_Accumulator & operator*=(const Double_t &value)
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
const QwSIS3320_Accumulator operator*(const Double_t &value) const
virtual const TString & GetElementName() const
Get the name of this element.
void Ratio(const QwSIS3320_Accumulator &numer, const QwSIS3320_Accumulator &denom)
void Sum(const QwSIS3320_Accumulator &value1, const QwSIS3320_Accumulator &value2)
QwSIS3320_Accumulator & operator+=(const Double_t &value)
void ConstructHistograms(TDirectory *folder, TString &prefix)
Construct the histograms for this data element.
#define QwWarning
Predefined log drain for warnings.
Definition: QwLog.h:45
Double_t GetNumberOfSamples() const
size_t fTreeArrayNumEntries
Index of this data element in tree.
TH1F * Construct1DHist(const TString &inputfile, const TString &name_title)
void FillHistograms()
Fill the histograms for this data element.
const QwSIS3320_Accumulator operator+(const Double_t &value) const
#define QwError
Predefined log drain for errors.
Definition: QwLog.h:40