Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaQWEAKHelicityReader.h
Go to the documentation of this file.
1#ifndef Podd_THaQWEAKHelicityReader_h_
2#define Podd_THaQWEAKHelicityReader_h_
3
5//
6// THaQWEAKHelicityReader
7//
8// Routines for decoding QWEAK helicity hardware
9//
11
12#include "Rtypes.h"
13
14class THaEvData;
15class TDatime;
16class TH1F;
17
19
20public:
22 virtual ~THaQWEAKHelicityReader() = default;
23
24 // when an event trigger the acquisition the reported helicity, T_settle and Pattern sync
25 //signals are recorded by the an Input Register.
26 UInt_t GetPatternTir() const { return fPatternTir; };
27 UInt_t GetHelicityTir() const { return fHelicityTir; };
28 UInt_t GetTSettleTir() const { return fTSettleTir; };
30 // in parallel, these signals are recorded in the ring buffer of the helicity gated scaler
31 // the maximum depth of this ring is for now 500.
32 UInt_t GetRingDepth() const { return fIRing; }; // how many events are in the ring
33 // I should have the depth of the ring defined in the data base
34 // let call it kHelRingDepth for now
35 void Print();
36
37 class ROCinfo {
38 public:
40 Bool_t valid() const;
41 UInt_t roc; // ROC to read out
42 UInt_t header; // Headers to search for (0 = ignore)
43 UInt_t index; // Index into buffer
44 };
45
46protected:
47
48 // Used by ReadDatabase
49 enum EROC { kHel = 0, kTime, kRing, kROC3 };
50 Int_t SetROCinfo( EROC which, UInt_t roc, UInt_t header, UInt_t index );
51
52 virtual void Clear( Option_t* opt="" );
53 virtual Int_t ReadData( const THaEvData& evdata );
54 Int_t ReadDatabase( const char* dbfilename, const char* prefix,
55 const TDatime& date, int debug_flag = 0 );
56 virtual void FillHisto();
57 void Begin();
58 void End();
59
66 // the ring map is defined in
67 // http://www.jlab.org/~adaq/halog/html/1011_archive/101101100943.html
68 // this map is valid after November 1rst
69 enum { kHelRingDepth = 500 };
77
78 // ROC information
79 // If a header is zero the index is taken to be from the start of
80 // the ROC (0 = first word of ROC), otherwise it's from the header
81 // (0 = first word after header).
82 ROCinfo fROCinfo[kROC3+1]; // for now only work with one roc (roc 11),
83 //should be extended to read helicity info from adc
84
85 Int_t fQWEAKDebug; // Debug level
86 Bool_t fHaveROCs; // Required ROCs are defined
87 Bool_t fNegGate; // Invert polarity of gate, TO DO implement this functionality
88 static const Int_t NHISTR = 12;
89 TH1F* fHistoR[12]; // Histograms
90
91private:
92
93 static UInt_t FindWord( const THaEvData& evdata, const ROCinfo& info );
94
95 ClassDef(THaQWEAKHelicityReader,0) // Helper class for reading QWEAK helicity data
96
97};
98
99#endif
100
int Int_t
unsigned int UInt_t
bool Bool_t
const UInt_t kMaxUInt
const char Option_t
#define ClassDef(name, id)
UInt_t fTimeStampRing[kHelRingDepth]
virtual ~THaQWEAKHelicityReader()=default
static UInt_t FindWord(const THaEvData &evdata, const ROCinfo &info)
UInt_t fPatternRing[kHelRingDepth]
UInt_t fT3Ring[kHelRingDepth]
UInt_t fHelicityRing[kHelRingDepth]
UInt_t fT10Ring[kHelRingDepth]
virtual void Clear(Option_t *opt="")
UInt_t fU3Ring[kHelRingDepth]
virtual Int_t ReadData(const THaEvData &evdata)
UInt_t fT5Ring[kHelRingDepth]
Int_t SetROCinfo(EROC which, UInt_t roc, UInt_t header, UInt_t index)
Int_t ReadDatabase(const char *dbfilename, const char *prefix, const TDatime &date, int debug_flag=0)