JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
LRBCorrector Class Reference

Linear-regression corrector applying per-burst slopes to data. More...

#include <LRBCorrector.h>

+ Inheritance diagram for LRBCorrector:
+ Collaboration diagram for LRBCorrector:

Public Member Functions

 LRBCorrector (const TString &name)
 Constructor with name.
 
void ParseConfigFile (QwParameterFile &file) override
 
Int_t LoadChannelMap (const std::string &mapfile) override
 
Int_t ConnectChannels (QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff) override
 
void ProcessData () override
 
void UpdateBurstCounter (Short_t burstcounter) override
 
- Public Member Functions inherited from VQwDataHandler
 VQwDataHandler (const TString &name)
 
 VQwDataHandler (const VQwDataHandler &source)
 
void SetPointer (QwHelicityPattern *ptr)
 
void SetPointer (QwSubsystemArrayParity *ptr)
 
virtual Int_t ConnectChannels (QwSubsystemArrayParity &, QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff)
 
virtual Int_t ConnectChannels (QwSubsystemArrayParity &)
 
Int_t ConnectChannels (QwHelicityPattern &helicitypattern)
 
virtual void FinishDataHandler ()
 
 ~VQwDataHandler () override
 
TString GetName ()
 
virtual void ClearEventData ()
 
void InitRunningSum ()
 
void AccumulateRunningSum ()
 
virtual void AccumulateRunningSum (VQwDataHandler &value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF)
 
void CalculateRunningAverage ()
 
void PrintValue () const
 
void WritePromptSummary (QwPromptSummary *ps, TString type)
 
virtual void ConstructTreeBranches (QwRootFile *treerootfile, const std::string &treeprefix="", const std::string &branchprefix="")
 
virtual void FillTreeBranches (QwRootFile *treerootfile)
 
virtual void ConstructNTupleFields (QwRootFile *treerootfile, const std::string &treeprefix="", const std::string &branchprefix="")
 RNTuple methods.
 
virtual void FillNTupleFields (QwRootFile *treerootfile)
 
virtual void ConstructHistograms (TDirectory *, TString &)
 Construct the histograms in a folder with a prefix.
 
virtual void FillHistograms ()
 Fill the histograms.
 
void FillTreeVector (QwRootTreeBranchVector &values) const
 
void ConstructBranchAndVector (TTree *tree, TString &prefix, QwRootTreeBranchVector &values)
 
void SetRunLabel (TString x)
 
Int_t LoadChannelMap ()
 
Bool_t PublishInternalValues () const override
 Publish all variables of the subsystem.
 
Bool_t PublishByRequest (TString device_name) override
 Try to publish an internal variable matching the submitted name.
 
- Public Member Functions inherited from MQwPublishable_child< QwDataHandlerArray, VQwDataHandler >
 MQwPublishable_child ()
 Default constructor Initializes the child object and sets up self-reference for publishing.
 
 MQwPublishable_child (const MQwPublishable_child &source)
 Copy constructor.
 
virtual ~MQwPublishable_child ()
 Virtual destructor.
 
void SetParent (QwDataHandlerArray *parent)
 Set the parent container for this child object.
 
QwDataHandlerArrayGetParent () const
 Get the parent container for this child object.
 
- Public Member Functions inherited from MQwCloneable< VQwDataHandler, LRBCorrector >
 ~MQwCloneable () override
 Virtual destructor.
 
VQwDataHandlerClone () const override
 Concrete clone method.
 
const VQwFactory< VQwDataHandler > * Factory () const override
 Factory getter.
 
- Public Member Functions inherited from VQwCloneable< VQwDataHandler >
virtual ~VQwCloneable ()
 Virtual destructor.
 
std::string GetClassName () const
 Get demangled name of this class.
 

Protected Member Functions

 LRBCorrector ()
 
- Protected Member Functions inherited from VQwDataHandler
 VQwDataHandler ()
 
void SetEventcutErrorFlagPointer (const UInt_t *errorflagptr)
 
UInt_t GetEventcutErrorFlag () const
 
std::pair< EQwHandleType, std::string > ParseHandledVariable (const std::string &variable)
 
void CalcOneOutput (const VQwHardwareChannel *dv, VQwHardwareChannel *output, std::vector< const VQwHardwareChannel * > &ivs, std::vector< Double_t > &sens)
 
- Protected Member Functions inherited from MQwPublishable_child< QwDataHandlerArray, VQwDataHandler >
Bool_t RequestExternalValue (const TString &name, VQwHardwareChannel *value) const
 Retrieve the variable name from other subsystem arrays Get the value corresponding to some variable name from a different data array.
 
const VQwHardwareChannelRequestExternalPointer (const TString &name) const
 Retrieve a pointer to an external variable by name Requests a direct pointer to a variable from sibling subsystems via the parent container.
 
Bool_t PublishInternalValue (const TString name, const TString desc, const VQwHardwareChannel *element) const
 Publish a variable from this child into the parent container.
 

Protected Attributes

std::string fAlphaFileBase
 
std::string fAlphaFileSuff
 
std::string fAlphaFilePath
 
std::vector< EQwHandleTypefIndependentType
 
std::vector< std::string > fIndependentName
 
std::vector< std::string > fIndependentFull
 
std::vector< const VQwHardwareChannel * > fIndependentVar
 
std::vector< Double_t > fIndependentValues
 
Short_t fLastCycle
 
std::map< Short_t, std::vector< std::vector< Double_t > > > fSensitivity
 
- Protected Attributes inherited from VQwDataHandler
Int_t fPriority
 
Short_t fBurstCounter
 When a datahandler array is processed, handlers with lower priority will be processed before handlers with higher priority.
 
TString fName
 
std::string fMapFile
 
std::string fTreeName
 
std::string fTreeComment
 
std::string fPrefix
 
TString run_label
 
const UInt_t * fErrorFlagPtr
 Error flag pointer.
 
QwSubsystemArrayParityfSubsystemArray
 Single event pointer.
 
QwHelicityPatternfHelicityPattern
 Helicity pattern pointer.
 
std::vector< std::string > fDependentFull
 
std::vector< EQwHandleTypefDependentType
 
std::vector< std::string > fDependentName
 
std::vector< const VQwHardwareChannel * > fDependentVar
 
std::vector< Double_t > fDependentValues
 
std::vector< VQwHardwareChannel * > fOutputVar
 
std::vector< Double_t > fOutputValues
 
std::vector< std::vector< TString > > fPublishList
 
std::string ParseSeparator
 
Bool_t fKeepRunningSum
 
Bool_t fRunningsumFillsTree
 
VQwDataHandlerfRunningsum
 

Additional Inherited Members

- Public Types inherited from VQwDataHandler
enum  EQwHandleType {
  kHandleTypeUnknown =0 , kHandleTypeMps , kHandleTypeAsym , kHandleTypeDiff ,
  kHandleTypeYield
}
 
typedef std::vector< VQwHardwareChannel * >::iterator Iterator_HdwChan
 
typedef std::vector< VQwHardwareChannel * >::const_iterator ConstIterator_HdwChan
 
- Static Public Member Functions inherited from MQwCloneable< VQwDataHandler, LRBCorrector >
static VQwDataHandlerCreate (const std::string &name)
 Object creation.
 
static LRBCorrectorCast (LRBCorrector *type)
 Object dynamic cast.
 

Detailed Description

Linear-regression corrector applying per-burst slopes to data.

Loads cycle-dependent sensitivities and applies linear regression corrections to monitored channels, selecting the appropriate set based on the current burst counter.

Definition at line 23 of file LRBCorrector.h.

Constructor & Destructor Documentation

◆ LRBCorrector() [1/2]

LRBCorrector::LRBCorrector ( const TString & name)

Constructor with name.

Definition at line 31 of file LRBCorrector.cc.

31 :VQwDataHandler(name),
32 fLastCycle(0)
33{
34 ParseSeparator = "_";
35 fKeepRunningSum = kTRUE;
36}
Short_t fLastCycle
std::string ParseSeparator

References VQwDataHandler::fKeepRunningSum, fLastCycle, VQwDataHandler::ParseSeparator, and VQwDataHandler::VQwDataHandler().

+ Here is the call graph for this function:

◆ LRBCorrector() [2/2]

LRBCorrector::LRBCorrector ( )
inlineprotected

Definition at line 53 of file LRBCorrector.h.

53{ }

Member Function Documentation

◆ ConnectChannels()

Int_t LRBCorrector::ConnectChannels ( QwSubsystemArrayParity & asym,
QwSubsystemArrayParity & diff )
overridevirtual

Connect to external pointers for dependent variables from asym/diff subsystem arrays, creating corresponding output variables as clones.

Fill vector of pointers to the relevant data elements

Reimplemented from VQwDataHandler.

Definition at line 182 of file LRBCorrector.cc.

185{
187
188 // Add independent variables
189 for (size_t iv = 0; iv < fIndependentName.size(); iv++) {
190 // Get the independent variables
191 const VQwHardwareChannel* iv_ptr = 0;
192 iv_ptr = this->RequestExternalPointer(fIndependentFull.at(iv));
193 if (iv_ptr==NULL){
194 switch (fIndependentType.at(iv)) {
195 case kHandleTypeAsym:
196 iv_ptr = asym.RequestExternalPointer(fIndependentName.at(iv));
197 break;
198 case kHandleTypeDiff:
199 iv_ptr = diff.RequestExternalPointer(fIndependentName.at(iv));
200 break;
201 default:
202 QwWarning << "Independent variable for corrector has unknown type."
203 << QwLog::endl;
204 break;
205 }
206 }
207 if (iv_ptr) {
208 //QwMessage << " iv: " << fIndependentName.at(iv) /*<< " (sens = " << fSensitivity.at(dv).at(iv) << ")"*/ << QwLog::endl;
209 fIndependentVar.push_back(iv_ptr);
210 } else {
211 QwWarning << "Independent variable " << fIndependentName.at(iv) << " could not be found."
212 << QwLog::endl;
213 }
214 }
215
216 QwMessage << "In LRBCorrector::ConnectChannels; Number of IVs: " << fIndependentVar.size()
217 << " Number of DVs: " << fDependentVar.size() << QwLog::endl;
218 return 0;
219}
#define QwWarning
Predefined log drain for warnings.
Definition QwLog.h:44
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
const VQwHardwareChannel * RequestExternalPointer(const TString &name) const
const VQwHardwareChannel * RequestExternalPointer(const TString &name) const
Retrieve a direct pointer to an external variable Searches for the named variable in external subsyst...
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
std::vector< std::string > fIndependentName
std::vector< EQwHandleType > fIndependentType
std::vector< const VQwHardwareChannel * > fIndependentVar
std::vector< std::string > fIndependentFull
std::vector< const VQwHardwareChannel * > fDependentVar
virtual Int_t ConnectChannels(QwSubsystemArrayParity &, QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff)

References VQwDataHandler::ConnectChannels(), QwLog::endl(), VQwDataHandler::fDependentVar, fIndependentFull, fIndependentName, fIndependentType, fIndependentVar, VQwDataHandler::kHandleTypeAsym, VQwDataHandler::kHandleTypeDiff, QwMessage, QwWarning, MQwPublishable< U, T >::RequestExternalPointer(), and MQwPublishable_child< QwDataHandlerArray, VQwDataHandler >::RequestExternalPointer().

+ Here is the call graph for this function:

◆ LoadChannelMap()

Int_t LRBCorrector::LoadChannelMap ( const std::string & mapfile)
overridevirtual

Reimplemented from VQwDataHandler.

Definition at line 46 of file LRBCorrector.cc.

47{
48 // Open the file
49 QwParameterFile map(mapfile);
50
51 // Read the sections of dependent variables
52 std::pair<EQwHandleType,std::string> type_name;
53
54 // Add independent variables and sensitivities
55 while (map.ReadNextLine()) {
56 // Throw away comments, whitespace, empty lines
57 map.TrimComment();
58 map.TrimWhitespace();
59 if (map.LineIsEmpty()) continue;
60 // Get first token: label (dv or iv), second token is the name like "asym_blah"
61 string primary_token = map.GetNextToken(" ");
62 string current_token = map.GetNextToken(" ");
63 // Parse current token into independent variable type and name
64 type_name = ParseHandledVariable(current_token);
65
66 if (primary_token == "iv") {
67 fIndependentType.push_back(type_name.first);
68 fIndependentName.push_back(type_name.second);
69 fIndependentFull.push_back(current_token);
70 }
71 else if (primary_token == "dv") {
72 fDependentType.push_back(type_name.first);
73 fDependentName.push_back(type_name.second);
74 fDependentFull.push_back(current_token);
75 }
76 else if (primary_token == "treetype") {
77 QwMessage << "Tree Type read, ignoring." << QwLog::endl;
78 }
79 else {
80 QwError << "Function LoadChannelMap in LRBCorrector.cc read in invalid primary_token." << QwLog::endl;
81 }
82 }
83
84 // Construct slope file name
85 std::string SlopeFileName = fAlphaFileBase + run_label.Data() + fAlphaFileSuff;
86 std::string SlopeFilePath = fAlphaFilePath + "/";
87 std::string SlopeFile = SlopeFilePath + SlopeFileName;
88 TString corFileName(SlopeFile);
89
90 QwMessage << "Trying to open " << corFileName << QwLog::endl;
91 TFile* corFile = new TFile(corFileName);
92 if (! corFile->IsOpen()) {
93 QwWarning << "Failed to open " << corFileName << ", slopes NOT found" << QwLog::endl;
94 return 0;
95 }
96
97 // DV names
98 TH1 *dvnames = (TH1 *) corFile->Get("DVname");
99 if (dvnames == 0) {
100 QwWarning << "DV names matrix is null" << QwLog::endl;
101 corFile->Close();
102 return 0;
103 }
104 for (size_t i = 0; i != fDependentName.size(); i++) {
105 TString name = dvnames->GetXaxis()->GetBinLabel(i+1);
106 name.Remove(0, name.First("_") + 1);
107 if (fDependentName[i] != name) {
108 QwWarning << "DV name expected differs from found: "
109 << fDependentName[i] << " != " << name
110 << QwLog::endl;
111 corFile->Close();
112 return 0;
113 }
114 }
115
116 // IV names
117 TH1 *ivnames = (TH1 *) corFile->Get("IVname");
118 if (ivnames == 0) {
119 QwWarning << "IV names matrix is null" << QwLog::endl;
120 corFile->Close();
121 return 0;
122 }
123 for (size_t i = 0; i != fIndependentName.size(); i++) {
124 TString name = ivnames->GetXaxis()->GetBinLabel(i+1);
125 name.Remove(0, name.First("_") + 1);
126 if (fIndependentName[i] != name) {
127 QwWarning << "IV name expected differs from found: "
128 << fIndependentName[i] << " != " << name
129 << QwLog::endl;
130 corFile->Close();
131 return 0;
132 }
133 }
134
135 // Slope matrix
136 TKey* key = corFile->GetKey("slopes");
137 if (key == 0) {
138 QwWarning << "No slope matrix found" << QwLog::endl;
139 corFile->Close();
140 return 0;
141 }
142 QwMessage << "Slope matrix has " << key->GetCycle() << " cycle(s)." << QwLog::endl;
143 fLastCycle = key->GetCycle(); // last cycle
144 for (Short_t cycle=1; cycle<=fLastCycle; cycle++){
145 TKey* key_cycle = corFile->GetKey("slopes", cycle);
146 TMatrixD *alphasM = (TMatrixD *) key_cycle->ReadObj();
147 if (alphasM == 0) {
148 QwWarning << "Slope matrix is null" << QwLog::endl;
149 corFile->Close();
150 return 0;
151 }
152 if (alphasM->GetNrows() != Int_t(fIndependentType.size())) {
153 QwWarning << "Slope matrix has wrong number of rows: "
154 << alphasM->GetNrows() << " != " << fIndependentType.size()
155 << QwLog::endl;
156 corFile->Close();
157 return 0;
158 }
159 if (alphasM->GetNcols() != Int_t(fDependentType.size())) {
160 QwWarning << "Slope matrix has wrong number of cols: "
161 << alphasM->GetNcols() << " != " << fDependentType.size()
162 << QwLog::endl;
163 corFile->Close();
164 return 0;
165 }
166
167 // Assign sensitivities
168 fSensitivity[cycle].resize(fDependentType.size());
169 for (size_t i = 0; i != fDependentType.size(); ++i) {
170 fSensitivity[cycle].at(i).resize(fIndependentType.size());
171 for (size_t j = 0; j != fIndependentType.size(); ++j) {
172 fSensitivity[cycle].at(i).at(j) = -1.0*(*alphasM)(j,i);
173 }
174 }
175 }
176
177 corFile->Close();
178 return 0;
179}
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
std::map< Short_t, std::vector< std::vector< Double_t > > > fSensitivity
std::string fAlphaFilePath
std::string fAlphaFileBase
std::string fAlphaFileSuff
std::vector< std::string > fDependentFull
std::pair< EQwHandleType, std::string > ParseHandledVariable(const std::string &variable)
std::vector< EQwHandleType > fDependentType
std::vector< std::string > fDependentName

References QwLog::endl(), fAlphaFileBase, fAlphaFilePath, fAlphaFileSuff, VQwDataHandler::fDependentFull, VQwDataHandler::fDependentName, VQwDataHandler::fDependentType, fIndependentFull, fIndependentName, fIndependentType, fLastCycle, fSensitivity, QwParameterFile::GetNextToken(), QwParameterFile::LineIsEmpty(), VQwDataHandler::ParseHandledVariable(), QwError, QwMessage, QwWarning, QwParameterFile::ReadNextLine(), VQwDataHandler::run_label, QwParameterFile::TrimComment(), and QwParameterFile::TrimWhitespace().

+ Here is the call graph for this function:

◆ ParseConfigFile()

void LRBCorrector::ParseConfigFile ( QwParameterFile & file)
overridevirtual

Parse configuration file for map file, priority, tree settings.

Reimplemented from VQwDataHandler.

Definition at line 38 of file LRBCorrector.cc.

39{
41 file.PopValue("slope-file-base", fAlphaFileBase);
42 file.PopValue("slope-file-suff", fAlphaFileSuff);
43 file.PopValue("slope-path", fAlphaFilePath);
44}
Bool_t PopValue(const std::string keyname, T &retvalue)
virtual void ParseConfigFile(QwParameterFile &file)

References fAlphaFileBase, fAlphaFilePath, fAlphaFileSuff, VQwDataHandler::ParseConfigFile(), and QwParameterFile::PopValue().

+ Here is the call graph for this function:

◆ ProcessData()

void LRBCorrector::ProcessData ( )
overridevirtual

Copy dependent variables to output variables (default processing).

Reimplemented from VQwDataHandler.

Definition at line 222 of file LRBCorrector.cc.

222 {
223 Short_t cycle = fBurstCounter+1;
224 if (fSensitivity.count(cycle) == 0) return;
225 for (size_t i = 0; i < fDependentVar.size(); ++i) {
227 }
228}
void CalcOneOutput(const VQwHardwareChannel *dv, VQwHardwareChannel *output, std::vector< const VQwHardwareChannel * > &ivs, std::vector< Double_t > &sens)
Short_t fBurstCounter
When a datahandler array is processed, handlers with lower priority will be processed before handlers...
std::vector< VQwHardwareChannel * > fOutputVar

References VQwDataHandler::CalcOneOutput(), VQwDataHandler::fBurstCounter, VQwDataHandler::fDependentVar, fIndependentVar, VQwDataHandler::fOutputVar, and fSensitivity.

+ Here is the call graph for this function:

◆ UpdateBurstCounter()

void LRBCorrector::UpdateBurstCounter ( Short_t burstcounter)
inlineoverridevirtual

Reimplemented from VQwDataHandler.

Definition at line 36 of file LRBCorrector.h.

36 {
37 if (burstcounter<fLastCycle){
38 fBurstCounter=burstcounter;
39 } else if (fLastCycle==1){
40 // If we have a single cycle of slopes, don't throw a warning.
41 fBurstCounter = 0;
42 } else {
44 QwWarning << "LRBCorrector, " << GetName()
45 << ": Burst counter, " << burstcounter
46 << ", is greater than the stored number of sets of slopes. "
47 << "Using the last set of slopes (cycle=" << fLastCycle
48 << ")" << QwLog::endl;
49 }
50 };
TString GetName()

References QwLog::endl(), VQwDataHandler::fBurstCounter, fLastCycle, VQwDataHandler::GetName(), and QwWarning.

+ Here is the call graph for this function:

Field Documentation

◆ fAlphaFileBase

std::string LRBCorrector::fAlphaFileBase
protected

Definition at line 55 of file LRBCorrector.h.

Referenced by LoadChannelMap(), and ParseConfigFile().

◆ fAlphaFilePath

std::string LRBCorrector::fAlphaFilePath
protected

Definition at line 57 of file LRBCorrector.h.

Referenced by LoadChannelMap(), and ParseConfigFile().

◆ fAlphaFileSuff

std::string LRBCorrector::fAlphaFileSuff
protected

Definition at line 56 of file LRBCorrector.h.

Referenced by LoadChannelMap(), and ParseConfigFile().

◆ fIndependentFull

std::vector< std::string > LRBCorrector::fIndependentFull
protected

Definition at line 61 of file LRBCorrector.h.

Referenced by ConnectChannels(), and LoadChannelMap().

◆ fIndependentName

std::vector< std::string > LRBCorrector::fIndependentName
protected

Definition at line 60 of file LRBCorrector.h.

Referenced by ConnectChannels(), and LoadChannelMap().

◆ fIndependentType

std::vector< EQwHandleType > LRBCorrector::fIndependentType
protected

Definition at line 59 of file LRBCorrector.h.

Referenced by ConnectChannels(), and LoadChannelMap().

◆ fIndependentValues

std::vector< Double_t > LRBCorrector::fIndependentValues
protected

Definition at line 64 of file LRBCorrector.h.

◆ fIndependentVar

std::vector< const VQwHardwareChannel* > LRBCorrector::fIndependentVar
protected

Definition at line 63 of file LRBCorrector.h.

Referenced by ConnectChannels(), and ProcessData().

◆ fLastCycle

Short_t LRBCorrector::fLastCycle
protected

Definition at line 66 of file LRBCorrector.h.

Referenced by LoadChannelMap(), LRBCorrector(), and UpdateBurstCounter().

◆ fSensitivity

std::map<Short_t,std::vector<std::vector<Double_t> > > LRBCorrector::fSensitivity
protected

Definition at line 67 of file LRBCorrector.h.

Referenced by LoadChannelMap(), and ProcessData().


The documentation for this class was generated from the following files: