JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwPMT_Channel.cc
Go to the documentation of this file.
1/*!
2 * \file QwPMT_Channel.cc
3 * \brief Implementation for PMT channel data element
4 * \author P. M. King
5 * \date 2009-03-07
6 */
7
8#include "QwPMT_Channel.h"
9
10// Qweak headers
11#include "QwLog.h"
12#include "QwHistogramHelper.h"
13#include "QwRootFile.h"
14
15const Bool_t QwPMT_Channel::kDEBUG = kFALSE;
16
17
18/*! Conversion factor to translate the average bit count in an ADC
19 * channel into average voltage.
20 * The base factor is roughly 76 uV per count, and zero counts corresponds
21 * to zero voltage.
22 * Store as the exact value for 20 V range, 18 bit ADC.
23 */
24const Double_t QwPMT_Channel::kPMT_VoltsPerBit = (20./(1<<18));
25
26
27/** \brief Clear the event-scoped ADC word value. */
31
32/**
33 * \brief Generate a mock ADC word for testing.
34 * \param helicity Helicity state indicator (unused).
35 * \param SlotNum V775 slot number to encode in the word.
36 * \param ChanNum V775 channel number to encode in the word.
37 */
38void QwPMT_Channel::RandomizeEventData(int helicity, int SlotNum, int ChanNum){
39
40 Double_t mean = 1500.0;
41 Double_t sigma = 300.0;
42 UInt_t fV775Dataword = abs( (Int_t)gRandom->Gaus(mean,sigma) );
43
44 UInt_t fV775SlotNumber = SlotNum;
45 UInt_t fV775ChannelNumber = ChanNum;
46 const UInt_t fV775DataValidBit = 0x00004000;
47
48 UInt_t word = fV775Dataword | (fV775SlotNumber<<27);
49 word = word | (fV775ChannelNumber<<16) | fV775DataValidBit;
50 fValue = word;
51}
52
53/** \brief Encode this channel's word into the trigger buffer. */
54void QwPMT_Channel::EncodeEventData(std::vector<UInt_t> &TrigBuffer)
55{
56// std::cout<<"QwPMT_Channel::EncodeEventData() not fully implemented yet."<<std::endl;
57
58 Long_t localbuf;
59
60 if (IsNameEmpty()) {
61 // This channel is not used, but is present in the data stream.
62 // Skip over this data.
63 } else {
64 localbuf = (Long_t) (this->fValue);
65 TrigBuffer.push_back(localbuf);
66 }
67
68}
69
70/** \brief Process the event (no-op for simple PMT channel). */
75
76
77/**
78 * \brief Create histograms for this channel within an optional folder.
79 * \param folder Optional ROOT folder to cd() into.
80 * \param prefix Histogram name prefix.
81 */
82void QwPMT_Channel::ConstructHistograms(TDirectory *folder, TString &prefix)
83{
84 // If we have defined a subdirectory in the ROOT file, then change into it.
85 if (folder != NULL) folder->cd();
86
87 if (GetElementName() == "") {
88 // This channel is not used, so skip filling the histograms.
89 } else {
90 // Now create the histograms.
91 TString basename, fullname;
92 basename = prefix + GetElementName();
93
94 fHistograms.resize(1, NULL);
95 size_t index = 0;
96 fHistograms[index] = gQwHists.Construct1DHist(basename);
97 index++;
98 }
99}
100
101/** \brief Fill histograms for this channel if present. */
103{
104 size_t index = 0;
105 if (GetElementName() == "") {
106 // This channel is not used, so skip creating the histograms.
107 } else {
108 if (fHistograms[index] != NULL)
109 fHistograms[index]->Fill(fValue);
110 index++;
111 }
112}
113
114/**
115 * \brief Construct a ROOT branch and append a value slot to the vector.
116 * \param tree Output tree.
117 * \param prefix Branch name prefix.
118 * \param values Output value vector to be appended.
119 */
120void QwPMT_Channel::ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)
121{
122 if (GetElementName() == "") {
123 // This channel is not used, so skip setting up the tree.
124 } else {
125 TString basename = prefix + GetElementName();
126 fTreeArrayIndex = values.size();
127
128 values.push_back(basename.Data(), 'D');
129
131 tree->Branch(basename, &(values[fTreeArrayIndex]), values.LeafList(fTreeArrayIndex).c_str());
132 }
133}
134
135/** \brief Write this channel's value into the tree vector slot. */
137{
138 if (GetElementName()==""){
139 // This channel is not used, so skip filling the tree vector.
140 } else if (fTreeArrayNumEntries<=0){
141 std::cerr << "QwPMT_Channel::FillTreeVector: fTreeArrayNumEntries=="
142 << fTreeArrayNumEntries << std::endl;
143 } else if (values.size() < fTreeArrayIndex+fTreeArrayNumEntries){
144 std::cerr << "QwPMT_Channel::FillTreeVector: values.size()=="
145 << values.size()
146 << "; fTreeArrayIndex+fTreeArrayNumEntries=="
148 << std::endl;
149 } else {
150 size_t index=fTreeArrayIndex;
151 values.SetValue(index++, this->fValue);
152 }
153}
154
155
156
157/** \brief Copy-assign from another PMT channel (event-scoped data). */
159 if (this != &value) {
160 if (GetElementName()!=""){
162 this->fValue = value.fValue;
163 }
164 }
165 return *this;
166}
167
169 if (GetElementName()!=""){
170 this->fValue += value.fValue;
171 }
172 return *this;
173}
174
176 if (GetElementName()!=""){
177 this->fValue -= value.fValue;
178 }
179 return *this;
180}
181
182void QwPMT_Channel::Sum(const QwPMT_Channel &value1, const QwPMT_Channel &value2) {
183 *this = value1;
184 *this += value2;
185}
186
187void QwPMT_Channel::Difference(const QwPMT_Channel &value1, const QwPMT_Channel &value2) {
188 *this = value1;
189 *this -= value2;
190}
191
192/** \brief Print a compact value summary for this PMT channel. */
194{
195 QwMessage << std::setprecision(4)
196 << std::setw(18) << std::left << GetElementName() << ", "
197 << std::setw(15) << std::left << GetValue()
198 << QwLog::endl;
199}
200
201/** \brief Print a placeholder info line for this PMT channel. */
203{
204 std::cout << "QwPMT_Channel::Print() not implemented yet." << std::endl;
205}
Helper functions and utilities for ROOT histogram management.
QwHistogramHelper gQwHists
Globally defined instance of the QwHistogramHelper class.
A logfile class, based on an identical class in the Hermes analyzer.
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
ROOT file and tree management wrapper classes.
PMT channel data element for tracking subsystem.
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)
Construct a ROOT branch and append a value slot to the vector.
void ProcessEvent()
Process the event (no-op for simple PMT channel).
void EncodeEventData(std::vector< UInt_t > &TrigBuffer)
Encode this channel's word into the trigger buffer.
void RandomizeEventData(int helicity, int SlotNum, int ChanNum)
Generate a mock ADC word for testing.
static const Double_t kPMT_VoltsPerBit
QwPMT_Channel & operator-=(const QwPMT_Channel &value)
static const Bool_t kDEBUG
void FillHistograms() override
Fill histograms for this channel if present.
void ClearEventData() override
Clear the event-scoped ADC word value.
void Sum(const QwPMT_Channel &value1, const QwPMT_Channel &value2)
void ConstructHistograms(TDirectory *folder, TString &prefix) override
Create histograms for this channel within an optional folder.
Double_t GetValue() const
void PrintInfo() const override
Print a placeholder info line for this PMT channel.
size_t fTreeArrayNumEntries
void FillTreeVector(QwRootTreeBranchVector &values) const
Write this channel's value into the tree vector slot.
QwPMT_Channel & operator+=(const QwPMT_Channel &value)
void PrintValue() const override
Print a compact value summary for this PMT channel.
QwPMT_Channel()
Default constructor.
void Difference(const QwPMT_Channel &value1, const QwPMT_Channel &value2)
size_t fTreeArrayIndex
QwPMT_Channel & operator=(const QwPMT_Channel &value)
Copy-assign from another PMT channel (event-scoped data).
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:53
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 SetValue(size_type index, Double_t val)
Definition QwRootFile.h:108
virtual const TString & GetElementName() const
Get the name of this element.
Bool_t IsNameEmpty() const
Is the name of this element empty?
VQwDataElement & operator=(const VQwDataElement &value)
Arithmetic assignment operator: Should only copy event-based data.