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

Data handler that forms linear combinations of channels. More...

#include <QwCombiner.h>

+ Inheritance diagram for QwCombiner:
+ Collaboration diagram for QwCombiner:

Public Types

typedef std::vector< VQwHardwareChannel * >::iterator Iterator_HdwChan
 
typedef std::vector< VQwHardwareChannel * >::const_iterator ConstIterator_HdwChan
 
- 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
 

Public Member Functions

 QwCombiner (const TString &name)
 Constructor with name.
 
 QwCombiner (const QwCombiner &source)
 Copy constructor.
 
 ~QwCombiner () override
 Virtual destructor.
 
Int_t LoadChannelMap (const std::string &mapfile) override
 Load the channels and sensitivities.
 
Int_t ConnectChannels (QwSubsystemArrayParity &event) override
 Connect to Channels (event only)
 
Int_t ConnectChannels (QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff) override
 Connect to Channels (asymmetry/difference only)
 
void ProcessData () override
 
- Public Member Functions inherited from VQwDataHandler
 VQwDataHandler (const TString &name)
 
 VQwDataHandler (const VQwDataHandler &source)
 
virtual void ParseConfigFile (QwParameterFile &file)
 
void SetPointer (QwHelicityPattern *ptr)
 
void SetPointer (QwSubsystemArrayParity *ptr)
 
virtual Int_t ConnectChannels (QwSubsystemArrayParity &, QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff)
 
Int_t ConnectChannels (QwHelicityPattern &helicitypattern)
 
virtual void UpdateBurstCounter (Short_t burstcounter)
 
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, QwCombiner >
 ~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

 QwCombiner ()
 Default constructor (Protected for child class access)
 
- 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

UInt_t fErrorFlagMask
 Error flag mask.
 
const UInt_t * fErrorFlagPointer
 
std::vector< std::vector< EQwHandleType > > fIndependentType
 List of channels to use in the combiner.
 
std::vector< std::vector< std::string > > fIndependentName
 
std::vector< std::vector< const VQwHardwareChannel * > > fIndependentVar
 
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

- Static Public Member Functions inherited from MQwCloneable< VQwDataHandler, QwCombiner >
static VQwDataHandlerCreate (const std::string &name)
 Object creation.
 
static QwCombinerCast (QwCombiner *type)
 Object dynamic cast.
 

Detailed Description

Data handler that forms linear combinations of channels.

QwCombiner connects to one or more input variables (yields/asymmetries/ differences) and computes configured linear combinations using provided sensitivities. It is typically used to derive composite diagnostics or to prepare inputs for higher-level analyses.

Definition at line 23 of file QwCombiner.h.

Member Typedef Documentation

◆ ConstIterator_HdwChan

typedef std::vector<VQwHardwareChannel*>::const_iterator QwCombiner::ConstIterator_HdwChan

Definition at line 27 of file QwCombiner.h.

◆ Iterator_HdwChan

typedef std::vector<VQwHardwareChannel*>::iterator QwCombiner::Iterator_HdwChan

Definition at line 26 of file QwCombiner.h.

Constructor & Destructor Documentation

◆ QwCombiner() [1/3]

QwCombiner::QwCombiner ( const TString & name)

Constructor with name.

Definition at line 28 of file QwCombiner.cc.

29: VQwDataHandler(name)
30{
31 ParseSeparator = ":";
32 fKeepRunningSum = kTRUE;
35}
const UInt_t * fErrorFlagPointer
Definition QwCombiner.h:61
UInt_t fErrorFlagMask
Error flag mask.
Definition QwCombiner.h:60
std::string ParseSeparator

References fErrorFlagMask, fErrorFlagPointer, VQwDataHandler::fKeepRunningSum, VQwDataHandler::ParseSeparator, and VQwDataHandler::VQwDataHandler().

Referenced by QwCombiner(), QwCombinerSubsystem::QwCombinerSubsystem(), and QwCombinerSubsystem::QwCombinerSubsystem().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ QwCombiner() [2/3]

QwCombiner::QwCombiner ( const QwCombiner & source)

Copy constructor.

Definition at line 37 of file QwCombiner.cc.

38: VQwDataHandler(source)
39{
42}

References fErrorFlagMask, fErrorFlagPointer, QwCombiner(), and VQwDataHandler::VQwDataHandler().

+ Here is the call graph for this function:

◆ ~QwCombiner()

QwCombiner::~QwCombiner ( )
override

Virtual destructor.

Destructor.

Definition at line 45 of file QwCombiner.cc.

46{
47 Iterator_HdwChan element;
48 for (element = fOutputVar.begin(); element != fOutputVar.end(); element++) {
49 if (*element != NULL){
50 delete *element;
51 }
52 }
53 fOutputVar.clear();
54}
std::vector< VQwHardwareChannel * >::iterator Iterator_HdwChan
Definition QwCombiner.h:26
std::vector< VQwHardwareChannel * > fOutputVar

References VQwDataHandler::fOutputVar.

◆ QwCombiner() [3/3]

QwCombiner::QwCombiner ( )
inlineprotected

Default constructor (Protected for child class access)

Definition at line 57 of file QwCombiner.h.

57{ };

Member Function Documentation

◆ ConnectChannels() [1/2]

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

Connect to Channels (asymmetry/difference only)

Parameters
asymSubsystem array with asymmetries
diffSubsystem array with differences
Returns
0 on success, non-zero on failure

Connect to the dependent and independent channels

Fill vector of pointers to the relevant data elements

Reimplemented from VQwDataHandler.

Definition at line 174 of file QwCombiner.cc.

177{
178 // Return if correction is not enabled
179
180 /// Fill vector of pointers to the relevant data elements
181 fIndependentVar.resize(fDependentName.size());
182 fDependentVar.resize(fDependentName.size());
183 fOutputVar.resize(fDependentName.size());
184
185 for (size_t dv = 0; dv < fDependentName.size(); dv++) {
186 // Add independent variables
187 for (size_t iv = 0; iv < fIndependentName.at(dv).size(); iv++) {
188 // Get the independent variables
189 const VQwHardwareChannel* iv_ptr = 0;
190 iv_ptr = RequestExternalPointer(fIndependentName.at(dv).at(iv));
191 if (iv_ptr == NULL){
192 switch (fIndependentType.at(dv).at(iv)) {
193 case kHandleTypeAsym:
194 iv_ptr = asym.RequestExternalPointer(fIndependentName.at(dv).at(iv));
195 break;
196 case kHandleTypeDiff:
197 iv_ptr = diff.RequestExternalPointer(fIndependentName.at(dv).at(iv));
198 break;
199 default:
200 QwWarning << "Independent variable for combiner has unknown type."
201 << QwLog::endl;
202 break;
203 }
204 }
205 if (iv_ptr) {
206 fIndependentVar[dv].push_back(iv_ptr);
207 } else {
208 QwWarning << "Independent variable " << fIndependentName.at(dv).at(iv) << " for combiner of "
209 << "dependent variable " << fDependentName.at(dv) << " could not be found."
210 << QwLog::endl;
211 }
212 }
213
214 // Get the dependent variables
215 const VQwHardwareChannel* dv_ptr = 0;
216 VQwHardwareChannel* new_chan = NULL;
217 const VQwHardwareChannel* chan = NULL;
218 string name = "";
219 string calc = "calc_";
220
221 if (fDependentType.at(dv)==kHandleTypeMps){
222 // Quietly ignore the MPS type when we're connecting the asym & diff
223 continue;
224 } else if(fDependentName.at(dv).at(0) == '@' ){
225 name = fDependentName.at(dv).substr(1,fDependentName.at(dv).length());
226 }else{
227 dv_ptr = this->RequestExternalPointer(fDependentFull.at(dv));
228 if (dv_ptr==NULL){
229 switch (fDependentType.at(dv)) {
230 case kHandleTypeAsym:
231 dv_ptr = asym.RequestExternalPointer(fDependentName.at(dv));
232 break;
233 case kHandleTypeDiff:
234 dv_ptr = diff.RequestExternalPointer(fDependentName.at(dv));
235 break;
236 default:
237 QwWarning << "QwCombiner::ConnectChannels(QwSubsystemArrayParity& asym, QwSubsystemArrayParity& diff): Dependent variable, "
238 << fDependentName.at(dv)
239 << ", for asym/diff combiner does not have proper type, type=="
240 << fDependentType.at(dv) << "."<< QwLog::endl;
241 break;
242 }
243 }
244
245 name = dv_ptr->GetElementName().Data();
246 name.insert(0, calc);
247
248 new_chan = dv_ptr->Clone(VQwDataElement::kDerived);
249 new_chan->SetElementName(name);
250 new_chan->SetSubsystemName(fName);
251 }
252
253 // alias
254 if(fDependentName.at(dv).at(0) == '@'){
255 //QwMessage << "dv: " << name << QwLog::endl;
256 if (fIndependentVar.at(dv).empty()) {
257 // Throw exception: alias cannot be created without independent variables
258 throw std::runtime_error("Cannot create alias '" + name +
259 "' for dependent variable '" + fDependentName.at(dv) +
260 "': no independent variables found to determine channel type");
261 } else {
262 // Preferred: use Clone() from first independent variable to preserve channel type
263 new_chan = fIndependentVar.at(dv).front()->Clone(VQwDataElement::kDerived);
264 }
265 new_chan->SetElementName(name);
266 new_chan->SetSubsystemName(fName);
267 }
268 // defined type
269 else if(dv_ptr!=NULL){
270 //QwMessage << "dv: " << fDependentName.at(dv) << QwLog::endl;
271 } else {
272 QwWarning << "Dependent variable " << fDependentName.at(dv) << " could not be found, "
273 << "or is not a known channel type." << QwLog::endl;
274 continue;
275 }
276
277 // pair creation
278 if(new_chan != NULL){
279 fDependentVar[dv] = chan;
280 fOutputVar[dv] = new_chan;
281 }
282 }
283
284 // Store error flag pointer
285 QwMessage << "Using asymmetry error flag" << QwLog::endl;
287
288 return 0;
289}
#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
void SetSubsystemName(TString sysname)
Set the name of the inheriting subsystem name.
virtual const TString & GetElementName() const
Get the name of this element.
void SetElementName(const TString &name)
Set the name of this element.
virtual VQwHardwareChannel * Clone() const
std::vector< std::vector< std::string > > fIndependentName
Definition QwCombiner.h:65
std::vector< std::vector< const VQwHardwareChannel * > > fIndependentVar
Definition QwCombiner.h:66
std::vector< std::vector< EQwHandleType > > fIndependentType
List of channels to use in the combiner.
Definition QwCombiner.h:64
const UInt_t * GetEventcutErrorFlagPointer() const
std::vector< std::string > fDependentFull
std::vector< const VQwHardwareChannel * > fDependentVar
std::vector< EQwHandleType > fDependentType
std::vector< std::string > fDependentName

References VQwHardwareChannel::Clone(), QwLog::endl(), VQwDataHandler::fDependentFull, VQwDataHandler::fDependentName, VQwDataHandler::fDependentType, VQwDataHandler::fDependentVar, fErrorFlagPointer, fIndependentName, fIndependentType, fIndependentVar, VQwDataHandler::fName, VQwDataHandler::fOutputVar, VQwDataElement::GetElementName(), QwSubsystemArrayParity::GetEventcutErrorFlagPointer(), VQwDataElement::kDerived, VQwDataHandler::kHandleTypeAsym, VQwDataHandler::kHandleTypeDiff, VQwDataHandler::kHandleTypeMps, QwMessage, QwWarning, MQwPublishable< U, T >::RequestExternalPointer(), MQwPublishable_child< QwDataHandlerArray, VQwDataHandler >::RequestExternalPointer(), VQwDataElement::SetElementName(), and VQwDataElement::SetSubsystemName().

+ Here is the call graph for this function:

◆ ConnectChannels() [2/2]

Int_t QwCombiner::ConnectChannels ( QwSubsystemArrayParity & event)
overridevirtual

Connect to Channels (event only)

Parameters
eventSubsystem array with per-MPS yields
Returns
0 on success, non-zero on failure

Connect to the dependent and independent channels

Parameters: event Helicity event structure Returns: Zero on success

Fill vector of pointers to the relevant data elements

Reimplemented from VQwDataHandler.

Definition at line 296 of file QwCombiner.cc.

297{
298 // Return if correction is not enabled
299
300 /// Fill vector of pointers to the relevant data elements
301 fIndependentVar.resize(fDependentName.size());
302 fDependentVar.resize(fDependentName.size());
303 fOutputVar.resize(fDependentName.size());
304
305 for (size_t dv = 0; dv < fDependentName.size(); dv++) {
306
307 // Add independent variables
308 for (size_t iv = 0; iv < fIndependentName.at(dv).size(); iv++) {
309 // Get the independent variables
310 const VQwHardwareChannel* iv_ptr = 0;
311 if(fIndependentType.at(dv).at(iv) == kHandleTypeMps){
312 iv_ptr = event.RequestExternalPointer(fIndependentName.at(dv).at(iv));
313 } else {
314 QwWarning << "Independent variable for MPS combiner has unknown type."
315 << QwLog::endl;
316 }
317 if (iv_ptr) {
318 fIndependentVar[dv].push_back(iv_ptr);
319 } else {
320 QwWarning << "Independent variable " << fIndependentName.at(dv).at(iv) << " for combiner of "
321 << "dependent variable " << fDependentName.at(dv) << " could not be found."
322 << QwLog::endl;
323 }
324 }
325
326 // Get the dependent variables
327 const VQwHardwareChannel* dv_ptr = 0;
328 VQwHardwareChannel* new_chan = NULL;
329 const VQwHardwareChannel* chan = NULL;
330 string name = " s";
331 string calc = "calc_";
332
334 // Quietly skip the asymmetry or difference types.
335 continue;
336 } else if(fDependentType.at(dv) != kHandleTypeMps){
337 QwWarning << "QwCombiner::ConnectChannels(QwSubsystemArrayParity& event): Dependent variable, "
338 << fDependentName.at(dv)
339 << ", for MPS combiner does not have MPS type, type=="
340 << fDependentType.at(dv) << "."<< QwLog::endl;
341 continue;
342 } else {
343 if(fDependentName.at(dv).at(0) == '@' ){
344 name = fDependentName.at(dv).substr(1,fDependentName.at(dv).length());
345 new_chan = fIndependentVar.at(dv).front()->Clone(VQwDataElement::kDerived);
346 new_chan->SetElementName(name);
347 new_chan->SetSubsystemName(fName);
348 } else {
349 dv_ptr = event.RequestExternalPointer(fDependentName.at(dv));
350
351 name = dv_ptr->GetElementName().Data();
352 name.insert(0,calc);
353
354 new_chan = dv_ptr->Clone(VQwDataElement::kDerived);
355 new_chan->SetElementName(name);
356 new_chan->SetSubsystemName(fName);
357 }
358 }
359
360 // alias
361 if(new_chan==NULL){
362 QwWarning << "Dependent variable " << fDependentName.at(dv) << " could not be found, "
363 << "or is not a known channel type." << QwLog::endl;
364 continue;
365 } else {
366 // pair creation
367 fDependentVar[dv] = chan;
368 fOutputVar[dv] = new_chan;
369 }
370 }
371
372 // Store error flag pointer
373 QwMessage << "Using event error flag" << QwLog::endl;
374 fErrorFlagPointer = event.GetEventcutErrorFlagPointer();
375
376 return 0;
377}

References VQwHardwareChannel::Clone(), QwLog::endl(), VQwDataHandler::fDependentName, VQwDataHandler::fDependentType, VQwDataHandler::fDependentVar, fErrorFlagPointer, fIndependentName, fIndependentType, fIndependentVar, VQwDataHandler::fName, VQwDataHandler::fOutputVar, VQwDataElement::GetElementName(), VQwDataElement::kDerived, VQwDataHandler::kHandleTypeAsym, VQwDataHandler::kHandleTypeDiff, VQwDataHandler::kHandleTypeMps, QwMessage, QwWarning, VQwDataElement::SetElementName(), and VQwDataElement::SetSubsystemName().

+ Here is the call graph for this function:

◆ LoadChannelMap()

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

Load the channels and sensitivities.

Load the channel map

Parameters
mapfileFilename of map file
Returns
Zero when success

Reimplemented from VQwDataHandler.

Definition at line 71 of file QwCombiner.cc.

72{
73 // Open the file
74 QwParameterFile map(mapfile);
75
76 // Read the preamble
77 std::unique_ptr<QwParameterFile> preamble = map.ReadSectionPreamble();
78 TString mask;
79 if (preamble->FileHasVariablePair("=", "mask", mask)) {
81 }
82
83
84 // Read the sections of dependent variables
85 bool keep_header = true;
86 std::string section_name;
87 std::unique_ptr<QwParameterFile> section = nullptr;
88 std::pair<EQwHandleType,std::string> type_name;
89 while ((section = map.ReadNextSection(section_name,keep_header))) {
90 if(section_name=="PUBLISH") continue;
91
92 // Store index to the current position in the dv vector
93 size_t current_dv_start = fDependentName.size();
94
95 // Add dependent variables from the section header
96 section->ReadNextLine();
97 if (section->LineHasSectionHeader()) {
98 section->TrimSectionHeader();
99 section->TrimWhitespace();
100 // Parse section header into tokens separated by a comma
101 std::string current_token;
102 std::string previous_token;
103 do {
104 previous_token = current_token;
105 current_token = section->GetNextToken(",");
106 if (current_token.size() == 0) continue;
107 // Parse current token into dependent variable type and name
108 type_name = ParseHandledVariable(current_token);
109 fDependentType.push_back(type_name.first);
110 fDependentName.push_back(type_name.second);
111 fDependentFull.push_back(current_token);
112 // Resize the vectors of sensitivities and independent variables
113 fSensitivity.resize(fDependentName.size());
114 fIndependentType.resize(fDependentName.size());
115 fIndependentName.resize(fDependentName.size());
116 } while (current_token.size() != 0);
117 } else QwError << "Section does not start with header." << QwLog::endl;
118
119 // Add independent variables and sensitivities
120 while (section->ReadNextLine()) {
121 // Throw away comments, whitespace, empty lines
122 section->TrimComment();
123 section->TrimWhitespace();
124 if (section->LineIsEmpty()) continue;
125 // Get first token: independent variable
126 std::string current_token = section->GetNextToken(",");
127 // Parse current token into independent variable type and name
128 type_name = ParseHandledVariable(current_token);
129 // Loop over dependent variables to set sensitivities
130 for (size_t dv = current_dv_start; dv < fDependentName.size(); dv++) {
131 Double_t sensitivity = atof(section->GetNextToken(",").c_str());
132 fSensitivity.at(dv).push_back(sensitivity);
133 fIndependentType.at(dv).push_back(type_name.first);
134 fIndependentName.at(dv).push_back(type_name.second);
135 }
136 }
137 }
138
139 TString varvalue;
140 // Now load the variables to publish
141 std::vector<std::vector<TString> > fPublishList;
142 map.RewindToFileStart();
143 std::unique_ptr<QwParameterFile> section2;
144 std::vector<TString> publishinfo;
145 while ((section2=map.ReadNextSection(varvalue))) {
146 if (varvalue == "PUBLISH") {
147 fPublishList.clear();
148 while (section2->ReadNextLine()) {
149 section2->TrimComment(); // Remove everything after a comment character
150 section2->TrimWhitespace(); // Get rid of leading and trailing spaces
151 for (int ii = 0; ii < 4; ii++) {
152 varvalue = section2->GetTypedNextToken<TString>();
153 if (varvalue.Length()) {
154 publishinfo.push_back(varvalue);
155 }
156 }
157 if (publishinfo.size() == 4)
158 fPublishList.push_back(publishinfo);
159 publishinfo.clear();
160 }
161 }
162 }
163 // Print list of variables to publish
164 if (fPublishList.size()>0){
165 QwMessage << "Variables to publish:" << QwLog::endl;
166 for (size_t jj = 0; jj < fPublishList.size(); jj++){
167 QwMessage << fPublishList.at(jj).at(0) << " " << fPublishList.at(jj).at(1) << " " << fPublishList.at(jj).at(2) << " " << fPublishList.at(jj).at(3) << QwLog::endl;
168 }
169 }
170 return 0;
171}
#define QwError
Predefined log drain for errors.
Definition QwLog.h:39
static UInt_t GetUInt(const TString &varvalue)
std::vector< std::vector< Double_t > > fSensitivity
Definition QwCombiner.h:67
std::vector< std::vector< TString > > fPublishList
std::pair< EQwHandleType, std::string > ParseHandledVariable(const std::string &variable)

References QwLog::endl(), VQwDataHandler::fDependentFull, VQwDataHandler::fDependentName, VQwDataHandler::fDependentType, fErrorFlagMask, fIndependentName, fIndependentType, VQwDataHandler::fPublishList, fSensitivity, QwParameterFile::GetUInt(), VQwDataHandler::ParseHandledVariable(), QwError, QwMessage, QwParameterFile::ReadNextSection(), QwParameterFile::ReadSectionPreamble(), and QwParameterFile::RewindToFileStart().

+ Here is the call graph for this function:

◆ ProcessData()

void QwCombiner::ProcessData ( )
overridevirtual

Copy dependent variables to output variables (default processing).

Reimplemented from VQwDataHandler.

Reimplemented in QwCombinerSubsystem.

Definition at line 379 of file QwCombiner.cc.

380{
381 if (fErrorFlagMask!=0 && fErrorFlagPointer!=NULL) {
382 if ((*fErrorFlagPointer & fErrorFlagMask)!=0) {
383 //QwMessage << "0x" << std::hex << *fErrorFlagPointer << " passed mask " << "0x" << fErrorFlagMask << std::dec << QwLog::endl;
384 for (size_t i = 0; i < fDependentVar.size(); ++i) {
386 }
387 //} else {
388 //QwMessage << "0x" << std::hex << *fErrorFlagPointer << " failed mask " << "0x" << fErrorFlagMask << std::dec << QwLog::endl;
389 }
390 }
391 else{
392 for (size_t i = 0; i < fDependentVar.size(); ++i) {
394 }
395 }
396}
void CalcOneOutput(const VQwHardwareChannel *dv, VQwHardwareChannel *output, std::vector< const VQwHardwareChannel * > &ivs, std::vector< Double_t > &sens)

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

Referenced by QwCombinerSubsystem::ProcessData().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Field Documentation

◆ fErrorFlagMask

UInt_t QwCombiner::fErrorFlagMask
protected

Error flag mask.

Definition at line 60 of file QwCombiner.h.

Referenced by LoadChannelMap(), ProcessData(), QwCombiner(), and QwCombiner().

◆ fErrorFlagPointer

const UInt_t* QwCombiner::fErrorFlagPointer
protected

Definition at line 61 of file QwCombiner.h.

Referenced by ConnectChannels(), ConnectChannels(), ProcessData(), QwCombiner(), and QwCombiner().

◆ fIndependentName

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

Definition at line 65 of file QwCombiner.h.

Referenced by ConnectChannels(), ConnectChannels(), and LoadChannelMap().

◆ fIndependentType

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

List of channels to use in the combiner.

Definition at line 64 of file QwCombiner.h.

Referenced by ConnectChannels(), ConnectChannels(), and LoadChannelMap().

◆ fIndependentVar

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

Definition at line 66 of file QwCombiner.h.

Referenced by ConnectChannels(), ConnectChannels(), and ProcessData().

◆ fSensitivity

std::vector< std::vector< Double_t > > QwCombiner::fSensitivity
protected

Definition at line 67 of file QwCombiner.h.

Referenced by LoadChannelMap(), and ProcessData().


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