JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
LRBCorrector.cc
Go to the documentation of this file.
1/*!
2 * \file LRBCorrector.cc
3 * \brief Implementation of linear regression blue corrector data handler
4 * \author Michael Vallee
5 * \date 2018-08-01
6 */
7
8#include <iostream>
9
10#include "QwHelicityPattern.h"
11
12#include "LRBCorrector.h"
13
14// Qweak headers
15#include "VQwDataElement.h"
16#include "QwVQWK_Channel.h"
17#include "QwParameterFile.h"
18
19#include <TFile.h>
20#include <TH2.h>
21#include <TKey.h>
22#include <TTree.h>
23#include <TChain.h>
24#include <TObjArray.h>
25#include <TEventList.h>
26
27#include <TMatrixD.h>
28
29
30/// \brief Constructor with name
31LRBCorrector::LRBCorrector(const TString& name):VQwDataHandler(name),
32 fLastCycle(0)
33{
34 ParseSeparator = "_";
35 fKeepRunningSum = kTRUE;
36}
37
39{
41 file.PopValue("slope-file-base", fAlphaFileBase);
42 file.PopValue("slope-file-suff", fAlphaFileSuff);
43 file.PopValue("slope-path", fAlphaFilePath);
44}
45
46Int_t LRBCorrector::LoadChannelMap(const std::string& mapfile)
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}
180
181
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}
220
221
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}
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
#define QwWarning
Predefined log drain for warnings.
Definition QwLog.h:44
#define QwMessage
Predefined log drain for regular messages.
Definition QwLog.h:49
Parameter file parsing and management.
Definition of the pure virtual base class of all data elements.
Decoding and management for VQWK ADC channels (6x32-bit datawords)
Linear regression blue corrector data handler class.
Helicity pattern analysis and management.
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
Configuration file parser with flexible tokenization and search capabilities.
Bool_t PopValue(const std::string keyname, T &retvalue)
void TrimWhitespace(TString::EStripType head_tail=TString::kBoth)
void TrimComment(const char commentchar)
std::string GetNextToken(const std::string &separatorchars)
Get next token as a string.
Abstract base for concrete hardware channels implementing dual-operator pattern.
std::map< Short_t, std::vector< std::vector< Double_t > > > fSensitivity
std::string fAlphaFilePath
Short_t fLastCycle
std::string fAlphaFileBase
void ParseConfigFile(QwParameterFile &file) override
std::vector< std::string > fIndependentName
void ProcessData() override
Int_t ConnectChannels(QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff) override
std::vector< EQwHandleType > fIndependentType
std::vector< const VQwHardwareChannel * > fIndependentVar
std::string fAlphaFileSuff
std::vector< std::string > fIndependentFull
Subsystem array container specialized for parity analysis with asymmetry calculations.
void CalcOneOutput(const VQwHardwareChannel *dv, VQwHardwareChannel *output, std::vector< const VQwHardwareChannel * > &ivs, std::vector< Double_t > &sens)
std::vector< std::string > fDependentFull
virtual void ParseConfigFile(QwParameterFile &file)
std::vector< const VQwHardwareChannel * > fDependentVar
VQwDataHandler(const TString &name)
std::string ParseSeparator
Short_t fBurstCounter
When a datahandler array is processed, handlers with lower priority will be processed before handlers...
std::vector< VQwHardwareChannel * > fOutputVar
std::pair< EQwHandleType, std::string > ParseHandledVariable(const std::string &variable)
std::vector< EQwHandleType > fDependentType
std::vector< std::string > fDependentName
virtual Int_t ConnectChannels(QwSubsystemArrayParity &, QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff)