JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwPMT_Channel.h
Go to the documentation of this file.
1/*!
2 * \file QwPMT_Channel.h
3 * \brief PMT channel data element for tracking subsystem
4 * \author P. M. King
5 * \date 2007-05-08
6 */
7
8#pragma once
9
10#include <vector>
11#include "TTree.h"
12
13//jpan: Mersenne Twistor: A 623-diminsionally equidistributed
14//uniform pseudorandom number generator
15#include "TRandom3.h"
16
17#include "VQwDataElement.h"
18
19
20/**
21 * \class QwPMT_Channel
22 * \ingroup QwTracking
23 * \brief PMT channel data element for photomultiplier tube readout
24 *
25 * Handles data from photomultiplier tube channels including raw values,
26 * calibration, and basic data element operations for PMT-based detectors
27 * in the tracking system.
28 */
30 /******************************************************************
31 * Class: QwPMT_Channel
32 *
33 ******************************************************************/
34 public:
35 /// Default constructor
37 // Prepare the random number generator.
38 gRandom->SetSeed();
39 };
40 /// Copy constructor
42 : VQwDataElement(source)
43 { };
44 /// Constructor with name
45 QwPMT_Channel(TString name) {
47 };
48 /// Virtual destructor
49 ~QwPMT_Channel() override { };
50
51 void InitializeChannel(TString name){
52 SetElementName(name);
54 };
55
56 void ClearEventData() override;
57 void RandomizeEventData(int helicity, int SlotNum, int ChanNum);
58 void EncodeEventData(std::vector<UInt_t> &TrigBuffer);
59 Int_t ProcessEvBuffer(UInt_t* buffer, UInt_t num_words_left, UInt_t subelement=0) override{return 0;};
60
61 void SetValue(Double_t data) { fValue = data; };
62 Double_t GetValue() const { return fValue; };
63 void SetSubbankID(const Int_t bank_index) { fCrate = bank_index; };
64 void SetModule(const Int_t slot_num) { fModule = slot_num; };
65 Int_t GetSubbankID() const { return fCrate; };
66 Int_t GetModule() const { return fModule; };
67
68 void ProcessEvent();
69
73 void Sum(const QwPMT_Channel &value1, const QwPMT_Channel &value2);
74 void Difference(const QwPMT_Channel &value1, const QwPMT_Channel &value2);
75
76 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
77 void FillHistograms() override;
78
79 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values);
80 void FillTreeVector(QwRootTreeBranchVector &values) const;
81
82 void PrintValue() const override;
83 void PrintInfo() const override;
84 void PrintErrorCounters() const override {};
85
86
87 protected:
88
89
90 private:
91 static const Bool_t kDEBUG;
92
93 /* ADC Calibration */
94 static const Double_t kPMT_VoltsPerBit;
95
96 /* Channel information data members */
97
98
99 /* Ntuple array indices */
102
103 /* Event data members */
104 Double_t fValue;
105
106 Int_t fCrate; ///< ROC number
107 Int_t fModule; ///< slot number
108};
Definition of the pure virtual base class of all data elements.
const TString QwBPMStripline< T >::subelement[4]
PMT channel data element for photomultiplier tube readout.
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 SetSubbankID(const Int_t bank_index)
QwPMT_Channel(TString name)
Constructor with name.
void SetValue(Double_t data)
void ClearEventData() override
Clear the event-scoped ADC word value.
void PrintErrorCounters() const override
report number of events failed due to HW and event cut failure
void Sum(const QwPMT_Channel &value1, const QwPMT_Channel &value2)
Int_t fCrate
ROC number.
Int_t fModule
slot number
void ConstructHistograms(TDirectory *folder, TString &prefix) override
Create histograms for this channel within an optional folder.
Double_t GetValue() const
Int_t GetModule() 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.
void InitializeChannel(TString name)
QwPMT_Channel & operator+=(const QwPMT_Channel &value)
void PrintValue() const override
Print a compact value summary for this PMT channel.
Int_t ProcessEvBuffer(UInt_t *buffer, UInt_t num_words_left, UInt_t subelement=0) override
Process the CODA event buffer for this element.
void SetModule(const Int_t slot_num)
QwPMT_Channel()
Default constructor.
void Difference(const QwPMT_Channel &value1, const QwPMT_Channel &value2)
size_t fTreeArrayIndex
QwPMT_Channel(const QwPMT_Channel &source)
Copy constructor.
QwPMT_Channel & operator=(const QwPMT_Channel &value)
Copy-assign from another PMT channel (event-scoped data).
~QwPMT_Channel() override
Virtual destructor.
Int_t GetSubbankID() const
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:53
VQwDataElement()
Default constructor.
void SetElementName(const TString &name)
Set the name of this element.