JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwScaler.h
Go to the documentation of this file.
1/*!
2 * \file QwScaler.h
3 * \brief Scaler subsystem for counting and rate measurements
4 */
5
6#pragma once
7
8// System headers
9#include <vector>
10
11// ROOT headers
12#ifdef HAS_RNTUPLE_SUPPORT
13#include "ROOT/RNTupleModel.hxx"
14#include "ROOT/RField.hxx"
15#endif // HAS_RNTUPLE_SUPPORT
16
17// Qweak headers
18#include "VQwSubsystemParity.h"
19#include "QwScaler_Channel.h"
20
21
22/**
23 * \class QwScaler
24 * \ingroup QwAnalysis_BL
25 * \brief Subsystem managing scaler modules and derived rates
26 *
27 * Wraps hardware scaler channels, provides per-MPS processing, histogram
28 * and tree output, and utilities for normalization and cuts.
29 */
30class QwScaler: public VQwSubsystemParity, public MQwSubsystemCloneable<QwScaler>
31{
32 private:
33 /// Private default constructor (not implemented, will throw linker error on use)
35
36 public:
37
38 /// Constructor with name
39 QwScaler(const TString& name);
40 /// Copy constructor
41 QwScaler(const QwScaler& source)
42 : VQwSubsystem(source),VQwSubsystemParity(source)
43 {
44 fScaler.resize(source.fScaler.size());
45 for (size_t i = 0; i < fScaler.size(); i++) {
46 VQwScaler_Channel* scaler_tmp = 0;
47 if ((scaler_tmp = dynamic_cast<QwSIS3801D24_Channel*>(source.fScaler.at(i)))) {
48 QwSIS3801D24_Channel* scaler = dynamic_cast<QwSIS3801D24_Channel*>(source.fScaler.at(i));
49 fScaler.at(i) = new QwSIS3801D24_Channel(*scaler);
50 } else if ((scaler_tmp = dynamic_cast<QwSIS3801D32_Channel*>(source.fScaler.at(i)))) {
51 QwSIS3801D32_Channel* scaler = dynamic_cast<QwSIS3801D32_Channel*>(source.fScaler.at(i));
52 fScaler.at(i) = new QwSIS3801D32_Channel(*scaler);
53 }
54 }
55 }
56 /// Destructor
57 ~QwScaler() override;
58
59 // Handle command line options
60 static void DefineOptions(QwOptions &options);
61 void ProcessOptions(QwOptions &options) override;
62
63 Int_t LoadChannelMap(TString mapfile) override;
64 Int_t LoadInputParameters(TString pedestalfile) override;
65
66 void ClearEventData() override;
67
68 Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words) override;
69 Int_t ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override;
70 void ProcessEvent() override;
71
73 void ConstructHistograms(TDirectory *folder, TString &prefix) override;
74 void FillHistograms() override;
75
77 void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override;
78 void ConstructBranch(TTree *tree, TString& prefix) override { };
79 void ConstructBranch(TTree *tree, TString& prefix, QwParameterFile& trim_file) override { };
80 void FillTreeVector(QwRootTreeBranchVector &values) const override;
81
82 // RNTuple methods
83#ifdef HAS_RNTUPLE_SUPPORT
84 void ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) override;
85 void FillNTupleVector(std::vector<Double_t>& values) const override;
86#endif // HAS_RNTUPLE_SUPPORT
87
88 Bool_t Compare(VQwSubsystem *source);
89
90 VQwSubsystem& operator=(VQwSubsystem *value) override;
91 VQwSubsystem& operator+=(VQwSubsystem *value) override;
92 VQwSubsystem& operator-=(VQwSubsystem *value) override;
93 void Ratio(VQwSubsystem *value1, VQwSubsystem *value2) override;
94 void Scale(Double_t factor) override;
95
96 void AccumulateRunningSum(VQwSubsystem* value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override;
97 //remove one entry from the running sums for devices
98 void DeaccumulateRunningSum(VQwSubsystem* value, Int_t ErrorMask=0xFFFFFFF) override;
99 void CalculateRunningAverage() override;
100
101 Int_t LoadEventCuts(TString filename) override;
103 Bool_t ApplySingleEventCuts() override;
104
105 Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override{
106 //QwError << "************* test inside Scaler *****************" << QwLog::endl;
107 return kFALSE;
108 };
109
110 void IncrementErrorCounters() override;
111
112 void PrintErrorCounters() const override;
113 UInt_t GetEventcutErrorFlag() override;
114 //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
115 void UpdateErrorFlag(const VQwSubsystem *ev_error) override{
116 };
117
118 void PrintValue() const override;
119 void PrintInfo() const override;
120
122
123 Double_t GetDataForChannelInModule(Int_t modnum, Int_t channum) {
124 Int_t index = fModuleChannel_Map[std::pair<Int_t,Int_t>(modnum,channum)];
125 return fScaler.at(index)->GetValue();
126 }
127
128 Int_t GetChannelIndex(TString channelName, UInt_t module_number);
129
130 private:
131
132 // Number of good events
134
135 // Mapping from subbank to scaler channels
136 typedef std::map< Int_t, std::vector< std::vector<Int_t> > > Subbank_to_Scaler_Map_t;
138
139 // Mapping from module and channel number to scaler channel
140 typedef std::map< std::pair<Int_t,Int_t>, Int_t > Module_Channel_to_Scaler_Map_t;
142
143 // Mapping from name to scaler channel
144 typedef std::map< TString, Int_t> Name_to_Scaler_Map_t;
146
147 // Vector of scaler channels
148 std::vector< VQwScaler_Channel* > fScaler; // Raw channels
149 std::vector< UInt_t > fBufferOffset; // Offset in scaler buffer
150 std::vector< std::pair< VQwScaler_Channel*, double > > fNorm;
151};
Base and derived classes for scaler channel data handling.
class QwScaler_Channel< 0x00ffffff, 0 > QwSIS3801D24_Channel
class QwScaler_Channel< 0xffffffff, 0 > QwSIS3801D32_Channel
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
Abstract base class for scaler channels.
Base class for subsystems implementing container-delegation pattern.
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
static void DefineOptions()
Define options function (note: no virtual static functions in C++)
VQwSubsystem(const TString &name)
Constructor with name.
virtual void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)=0
Construct the branch and tree vector.
void ClearEventData() override
Definition QwScaler.cc:281
void Scale(Double_t factor) override
Definition QwScaler.cc:468
Int_t LoadChannelMap(TString mapfile) override
Definition QwScaler.cc:67
Int_t LoadInputParameters(TString pedestalfile) override
Definition QwScaler.cc:246
Bool_t Compare(VQwSubsystem *source)
Definition QwScaler.cc:558
Int_t LoadEventCuts(TString filename) override
Optional event cut file.
Definition QwScaler.cc:511
void AccumulateRunningSum(VQwSubsystem *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
Definition QwScaler.cc:478
void ProcessOptions(QwOptions &options) override
Definition QwScaler.cc:32
QwScaler(const QwScaler &source)
Copy constructor.
Definition QwScaler.h:41
Bool_t ApplySingleEventCuts() override
Apply the single event cuts.
Definition QwScaler.cc:516
void CalculateRunningAverage() override
Definition QwScaler.cc:504
void DeaccumulateRunningSum(VQwSubsystem *value, Int_t ErrorMask=0xFFFFFFF) override
Definition QwScaler.cc:491
VQwSubsystem & operator-=(VQwSubsystem *value) override
Definition QwScaler.cc:439
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
Name_to_Scaler_Map_t fName_Map
Definition QwScaler.h:145
void Ratio(VQwSubsystem *value1, VQwSubsystem *value2) override
Definition QwScaler.cc:454
Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override
Report the number of events failed due to HW and event cut failures.
Definition QwScaler.h:105
std::map< TString, Int_t > Name_to_Scaler_Map_t
Definition QwScaler.h:144
Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
Definition QwScaler.cc:299
std::vector< std::pair< VQwScaler_Channel *, double > > fNorm
Definition QwScaler.h:150
Double_t GetDataForChannelInModule(Int_t modnum, Int_t channum)
Definition QwScaler.h:123
void FillHistograms() override
Definition QwScaler.cc:368
Module_Channel_to_Scaler_Map_t fModuleChannel_Map
Definition QwScaler.h:141
Int_t ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
Definition QwScaler.cc:309
UInt_t GetEventcutErrorFlag() override
Return the error flag to the top level routines related to stability checks and ErrorFlag updates.
Definition QwScaler.cc:542
std::vector< UInt_t > fBufferOffset
Definition QwScaler.h:149
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
Definition QwScaler.cc:376
std::map< std::pair< Int_t, Int_t >, Int_t > Module_Channel_to_Scaler_Map_t
Definition QwScaler.h:140
void IncrementErrorCounters() override
Increment the error counters.
Definition QwScaler.cc:534
std::vector< VQwScaler_Channel * > fScaler
Definition QwScaler.h:148
QwScaler()
Private default constructor (not implemented, will throw linker error on use)
Bool_t SingleEventCuts()
VQwSubsystem & operator=(VQwSubsystem *value) override
Definition QwScaler.cc:410
void FillTreeVector(QwRootTreeBranchVector &values) const override
Definition QwScaler.cc:384
void PrintErrorCounters() const override
Definition QwScaler.cc:538
~QwScaler() override
Destructor.
Definition QwScaler.cc:50
Subbank_to_Scaler_Map_t fSubbank_Map
Definition QwScaler.h:137
VQwSubsystem & operator+=(VQwSubsystem *value) override
Definition QwScaler.cc:425
void PrintInfo() const override
Definition QwScaler.cc:580
void UpdateErrorFlag(const VQwSubsystem *ev_error) override
update the error flag in the subsystem level from the top level routines related to stability checks....
Definition QwScaler.h:115
std::map< Int_t, std::vector< std::vector< Int_t > > > Subbank_to_Scaler_Map_t
Definition QwScaler.h:136
void ConstructBranch(TTree *tree, TString &prefix) override
Construct the branch and tree vector.
Definition QwScaler.h:78
Int_t GetChannelIndex(TString channelName, UInt_t module_number)
Definition QwScaler.cc:547
Double_t * GetRawChannelArray()
void ProcessEvent() override
Definition QwScaler.cc:343
void ConstructBranch(TTree *tree, TString &prefix, QwParameterFile &trim_file) override
Construct the branch and tree vector based on the trim file.
Definition QwScaler.h:79
void PrintValue() const override
Definition QwScaler.cc:595
Int_t fGoodEventCount
Definition QwScaler.h:133
VQwSubsystemParity()
Private default constructor (not implemented, will throw linker error on use)