JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwBeamMod.h
Go to the documentation of this file.
1/*!
2 * \file QwBeamMod.h
3 * \brief Beam modulation subsystem for parity analysis
4 * \author Joshua Hoskins
5 * \date 2010-05-25
6 */
7
8#pragma once
9
10// System headers
11#include <vector>
12
13// ROOT headers
14#include "TTree.h"
15#ifdef HAS_RNTUPLE_SUPPORT
16#include <ROOT/RNTupleModel.hxx>
17#endif // HAS_RNTUPLE_SUPPORT
18
19// Qweak headers
20#include "VQwSubsystemParity.h"
21#include "QwVQWK_Channel.h"
22#include "QwWord.h"
23
24// Forward declarations
25class QwBeamMod;
26
27//enum EBeamInstrumentType{kBPMStripline = 0,
28// kBCM,
29// kCombinedBCM,
30// kCombinedBPM};
31
32// this emun vector needs to be coherent with the DetectorTypes declaration in the QwBeamLine constructor
33
34
35
36/*****************************************************************
37* Class:
38******************************************************************/
39/**
40 * \class QwModChannelID
41 * \ingroup QwAnalysis_BL
42 * \brief Mapping information for beam modulation channels
43 */
45{
46 public:
47
48 QwModChannelID(Int_t subbankid, Int_t wordssofar,TString name,
49 TString modtype ,QwBeamMod * obj);
50
51 QwModChannelID(Int_t subbankid, QwParameterFile &paramfile);
52
53
54
55/* QwModChannelID):fSubbankIndex(-1),fWordInSubbank(-1),fTypeID(-1),fIndex(-1), */
56/* fSubelement(999999),fmoduletype(""),fmodulename("") */
57/* {}; */
58
59 Int_t fSubbankIndex; //Generated from ROCID(readout CPU) & BankID(corresponds to internal headers to ID different types of data)
61 //first word reported for this channel in the subbank
62 //(eg VQWK channel report 6 words for each event, scalers only report one word per event)
63
64 // The first word of the subbank gets fWordInSubbank=0
65
66 TString fmoduletype; // eg: VQWK, SCALER
67 TString fmodulename;
68 Int_t modnum;
69 Int_t channum;
70
71 // TString fdetectortype;
72
74 Int_t fTypeID; // type of detector eg: lumi or stripline, etc..
75 Int_t fIndex; // index of this detector in the vector containing all the detector of same type
76
77 void Print();
78
79};
80
81
82
83/*****************************************************************
84* Class:
85******************************************************************/
86/**
87 * \class QwBeamMod
88 * \ingroup QwAnalysis_BL
89 * \brief Subsystem for beam modulation studies and FFB handling
90 *
91 * Decodes modulation ramp and pattern words, maintains channels affected
92 * by modulation, and computes relevant summaries for regression.
93 */
94class QwBeamMod: public VQwSubsystemParity, public MQwSubsystemCloneable<QwBeamMod> {
95
96 private:
97 ///
99
100 public:
101 /// Constructor with name
102 QwBeamMod(const TString& name)
103 : VQwSubsystem(name),VQwSubsystemParity(name)
104 {
106 fFFB_Flag=kTRUE;
109 fBmwObj_Index = -1;
110
111 // Initialize the bmwobj cuts with UL < LL to disable the cut
112 fBMWObj_LL = +1;
113 fBMWObj_UL = -1;
114 };
115 /// Copy constructor
116 QwBeamMod(const QwBeamMod& source)
117 : VQwSubsystem(source),VQwSubsystemParity(source),
118 fWord(source.fWord)
119 {
120 // std::cout<< "Here in the copy constructor" << std::endl;
121 this->fModChannel.reserve(source.fModChannel.size());
122 for(size_t i=0;i< source.fModChannel.size();i++) {
123 this->fModChannel.push_back(source.fModChannel[i]->Clone());
124 *(this->fModChannel[i]) = *(source.fModChannel[i]);
125 //source.fModChannel[i]->PrintValue();
126 }
127 }
128 /// Virtual destructor
129 ~QwBeamMod() override {};
130
131 /* derived from VQwSubsystem */
132
133 //Handle command line options
134 void ProcessOptions(QwOptions &options) override;
135 void AccumulateRunningSum(VQwSubsystem*, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override;
136 //remove one entry from the running sums for devices
137 void DeaccumulateRunningSum(VQwSubsystem* value, Int_t ErrorMask=0xFFFFFFF) override{
138 };
139
140 Int_t LoadChannelMap(TString mapfile) override;
141 void LoadEventCuts_Init() override;
142 void LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag) override;
143 void LoadEventCuts_Fin(Int_t &eventcut_flag) override;
144 Int_t LoadGeometry(TString mapfile);
145 Int_t LoadInputParameters(TString pedestalfile) override;
146
147
148 Bool_t ApplySingleEventCuts() override;//derived from VQwSubsystemParity
149 void IncrementErrorCounters() override;
150 void PrintErrorCounters() const override;// report number of events failed due to HW and event cut failures
151 UInt_t GetEventcutErrorFlag() override;//return the error flag
152
153 Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override;
154
155 //update the error flag in the subsystem level from the top level routines related to stability checks. This will uniquely update the errorflag at each channel based on the error flag in the corresponding channel in the ev_error subsystem
156 void UpdateErrorFlag(const VQwSubsystem *ev_error) override;
157
158 Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words) override;
159 Int_t ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words) override;
160// void PrintDetectorID();
161
162 void ClearEventData() override;
163
164 void ProcessEvent() override;
165 void ProcessEvent_2() override;
166
167 VQwSubsystem& operator= (VQwSubsystem *value) override;
168 VQwSubsystem& operator+= (VQwSubsystem *value) override;
169 VQwSubsystem& operator-= (VQwSubsystem *value) override;
170
171
172 void Ratio(VQwSubsystem *numer, VQwSubsystem *denom) override;
173
174 void Scale(Double_t factor) override;
175
176 void CalculateRunningAverage() override;
177 void PrintModChannelID();
178
180 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
181 void FillHistograms() override;
182
184 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override;
185 void ConstructBranch(TTree *tree, TString& prefix) override { };
186 void ConstructBranch(TTree *tree, TString& prefix, QwParameterFile& trim_file) override { };
187 void FillTreeVector(QwRootTreeBranchVector &values) const override;
188#ifdef HAS_RNTUPLE_SUPPORT
189 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) override;
190 void FillNTupleVector(std::vector<Double_t>& values) const override;
191#endif // HAS_RNTUPLE_SUPPORT
192
193#ifdef __USE_DATABASE__
194 void FillDB_MPS(QwParityDB *db, TString datatype);
195 void FillDB(QwParityDB *db, TString datatype);
196 void FillErrDB(QwParityDB *db, TString datatype);
197#endif // __USE_DATABASE__
198 void WritePromptSummary(QwPromptSummary *ps, TString type) override;
199
200 Bool_t Compare(VQwSubsystem *source);
201
202 void Print();
203
204 protected:
205 Int_t GetDetectorIndex(TString name);
207
208 std::vector <VQwHardwareChannel*> fModChannel;
209 std::vector <QwModChannelID> fModChannelID;
210 std::vector <QwWord> fWord;
211 std::vector < std::pair<Int_t, Int_t> > fWordsPerSubbank;
212
213
214
215 private:
216
221 Bool_t fFFB_Flag;
222 static const Bool_t bDEBUG=kFALSE;
223
230
231};
Word-level data manipulation and bit operations.
Decoding and management for VQWK ADC channels (6x32-bit datawords)
ULong64_t BankID_t
Definition QwTypes.h:21
UInt_t ROCID_t
Definition QwTypes.h:20
Virtual base class for parity analysis subsystems.
Command-line and configuration file options processor.
Definition QwOptions.h:141
Configuration file parser with flexible tokenization and search capabilities.
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:53
Base class for subsystems implementing container-delegation pattern.
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
VQwSubsystem(const TString &name)
Constructor with name.
virtual void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)=0
Construct the branch and tree vector.
Int_t fWordInSubbank
Definition QwBeamMod.h:60
QwModChannelID(Int_t subbankid, Int_t wordssofar, TString name, TString modtype, QwBeamMod *obj)
TString fmodulename
Definition QwBeamMod.h:67
Int_t kUnknownDeviceType
Definition QwBeamMod.h:73
Int_t fSubbankIndex
Definition QwBeamMod.h:59
TString fmoduletype
Definition QwBeamMod.h:66
Subsystem for beam modulation studies and FFB handling.
Definition QwBeamMod.h:94
void ClearEventData() override
Definition QwBeamMod.cc:579
void LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag) override
Definition QwBeamMod.cc:285
static const Bool_t bDEBUG
Definition QwBeamMod.h:222
void LoadEventCuts_Init() override
Definition QwBeamMod.cc:281
VQwSubsystem & operator-=(VQwSubsystem *value) override
Definition QwBeamMod.cc:663
Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override
Report the number of events failed due to HW and event cut failures.
Definition QwBeamMod.cc:325
void PrintErrorCounters() const override
Definition QwBeamMod.cc:521
Int_t fRampChannelIndex
Definition QwBeamMod.h:224
UInt_t fBmwObj_Index
Definition QwBeamMod.h:226
void PrintModChannelID()
Definition QwBeamMod.cc:882
VQwSubsystem & operator=(VQwSubsystem *value) override
Assignment Note: Must be called at the beginning of all subsystems routine call to operator=(VQwSubsy...
Definition QwBeamMod.cc:625
VQwSubsystem & operator+=(VQwSubsystem *value) override
Definition QwBeamMod.cc:646
Int_t fBMWObj_UL
Definition QwBeamMod.h:228
Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
Definition QwBeamMod.cc:572
Int_t fBMWObj_LL
Definition QwBeamMod.h:227
UInt_t fFFB_holdoff
Definition QwBeamMod.h:218
void ProcessOptions(QwOptions &options) override
Process the command line options.
Definition QwBeamMod.cc:42
void AccumulateRunningSum(VQwSubsystem *, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
Update the running sums for devices.
Definition QwBeamMod.cc:706
Int_t LoadChannelMap(TString mapfile) override
Mandatory map file definition.
Definition QwBeamMod.cc:46
void FillHistograms() override
Fill the histograms for this subsystem.
Definition QwBeamMod.cc:746
void ProcessEvent_2() override
Process the event data again, including data from other subsystems. Not all derived classes will requ...
Definition QwBeamMod.cc:566
Int_t LoadInputParameters(TString pedestalfile) override
Mandatory parameter file definition.
Definition QwBeamMod.cc:339
Int_t LoadGeometry(TString mapfile)
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
void ConstructBranch(TTree *tree, TString &prefix) override
Construct the branch and tree vector.
Definition QwBeamMod.h:185
QwBeamMod(const QwBeamMod &source)
Copy constructor.
Definition QwBeamMod.h:116
UInt_t fFFB_ErrorFlag
Definition QwBeamMod.h:220
void CalculateRunningAverage() override
Calculate the average for all good events.
Definition QwBeamMod.cc:704
void WritePromptSummary(QwPromptSummary *ps, TString type) override
Definition QwBeamMod.cc:976
void Ratio(VQwSubsystem *numer, VQwSubsystem *denom) override
Definition QwBeamMod.cc:680
void FillTreeVector(QwRootTreeBranchVector &values) const override
Fill the tree vector.
Definition QwBeamMod.cc:816
Int_t fPatternWordIndex
Definition QwBeamMod.h:225
std::vector< VQwHardwareChannel * > fModChannel
Definition QwBeamMod.h:208
UInt_t GetEventcutErrorFlag() override
Return the error flag to the top level routines related to stability checks and ErrorFlag updates.
Definition QwBeamMod.cc:542
~QwBeamMod() override
Virtual destructor.
Definition QwBeamMod.h:129
Int_t GetDetectorIndex(TString name)
Definition QwBeamMod.cc:591
void Scale(Double_t factor) override
Definition QwBeamMod.cc:698
Bool_t Compare(VQwSubsystem *source)
Definition QwBeamMod.cc:708
void ConstructBranch(TTree *tree, TString &prefix, QwParameterFile &trim_file) override
Construct the branch and tree vector based on the trim file.
Definition QwBeamMod.h:186
QwBeamMod(const TString &name)
Constructor with name.
Definition QwBeamMod.h:102
std::vector< QwWord > fWord
Definition QwBeamMod.h:210
Int_t ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
TODO: The non-event-type-aware ProcessEvBuffer routine should be replaced with the event-type-aware v...
Definition QwBeamMod.cc:384
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
Construct the branch and tree vector.
Definition QwBeamMod.cc:797
void ProcessEvent() override
Definition QwBeamMod.cc:556
UInt_t fFFB_holdoff_Counter
Definition QwBeamMod.h:219
void Print()
Definition QwBeamMod.cc:869
UInt_t fFFB_Index
Definition QwBeamMod.h:217
UInt_t fBmwObj_ErrorFlag
Definition QwBeamMod.h:229
void IncrementErrorCounters() override
Increment the error counters.
Definition QwBeamMod.cc:514
std::vector< std::pair< Int_t, Int_t > > fWordsPerSubbank
Definition QwBeamMod.h:211
std::vector< QwModChannelID > fModChannelID
Definition QwBeamMod.h:209
Bool_t ApplySingleEventCuts() override
Apply the single event cuts.
Definition QwBeamMod.cc:450
Int_t fTreeArrayIndex
Definition QwBeamMod.h:206
void DeaccumulateRunningSum(VQwSubsystem *value, Int_t ErrorMask=0xFFFFFFF) override
remove one entry from the running sums for devices
Definition QwBeamMod.h:137
Bool_t fFFB_Flag
Definition QwBeamMod.h:221
void LoadEventCuts_Fin(Int_t &eventcut_flag) override
Definition QwBeamMod.cc:318
VQwSubsystemParity()
Private default constructor (not implemented, will throw linker error on use)
virtual void FillDB_MPS(QwParityDB *, TString)
Fill the database with MPS-based variables Note that most subsystems don't need to do this.
virtual void FillDB(QwParityDB *, TString)
Fill the database.
virtual UInt_t UpdateErrorFlag()
Uses the error flags of contained data elements to update Returns the error flag to the top level rou...
virtual void FillErrDB(QwParityDB *, TString)