JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
VQwDataHandler.cc
Go to the documentation of this file.
1/*!
2 * \file VQwDataHandler.cc
3 * \brief Virtual base class implementation for data handlers accessing multiple subsystems
4 *
5 * Base implementation for data handlers that process inputs from multiple
6 * subsystems, creating corrected output channels and managing tree/RNTuple
7 * construction, running sums, and variable publication. Used by concrete
8 * handlers like combiners and correctors. Documentation-only edits; runtime
9 * behavior unchanged.
10 */
11
12#include <iostream>
13
14using namespace std;
15
16//header file
17#include "VQwDataHandler.h"
18
19#ifdef HAS_RNTUPLE_SUPPORT
20// ROOT headers for RNTuple support
21#include <ROOT/RNTupleModel.hxx>
22#include <ROOT/RNTupleWriter.hxx>
23#endif
24
25#include "QwParameterFile.h"
26#include "QwRootFile.h"
27#include "QwVQWK_Channel.h"
28#include "QwPromptSummary.h"
29
30#ifdef __USE_DATABASE__
31#include "QwParityDB.h"
32#endif // __USE_DATABASE__
33
34
35/** Constructor: initialize base data handler with name and defaults. */
37: fPriority(0),
39 fName(name),
40 fMapFile(""),
41 fTreeName(""),
42 fTreeComment(""),
43 fPrefix(""),
47 fKeepRunningSum(kFALSE)
48{ fRunningsum=NULL;}
49
50/** Copy constructor: deep-copy output variables and running sums. */
52: fPriority(source.fPriority),
54 fName(source.fName),
55 fMapFile(source.fMapFile),
56 fTreeName(source.fTreeName),
58 fPrefix(source.fPrefix),
63{
68 // Create new objects for the the outputs.
69 fOutputVar.resize(source.fOutputVar.size());
70 for (size_t i = 0; i < this->fDependentVar.size(); i++) {
71 //const QwVQWK_Channel* vqwk = dynamic_cast<const QwVQWK_Channel*>(source.fOutputVar[i]);
72 this->fOutputVar[i] = source.fOutputVar[i]->Clone(VQwDataElement::kDerived);
73 //this->fOutputVar[i] = new QwVQWK_Channel(*vqwk, VQwDataElement::kDerived);
74 }
75 if (source.fRunningsum!=NULL){
76 fRunningsum = source.fRunningsum->Clone();
77 } else {
78 fRunningsum = NULL;
79 }
80}
81
82
83/** Destructor: clean up owned output variable clones. */
85
86 for (size_t i = 0; i < fOutputVar.size(); ++i) {
87 if (fOutputVar.at(i) != NULL){
88 delete fOutputVar.at(i);
89 }
90 }
91 fOutputVar.clear();
92
93}
94
95/** Parse configuration file for map file, priority, tree settings. */
97 file.RewindToFileStart();
98 file.EnableGreediness();
99 while (file.ReadNextLine()) {
100 QwMessage << file.GetLine() << QwLog::endl;
101 }
102 // Check for and process key-value pairs
103 file.PopValue("map",fMapFile);
104 file.PopValue("priority",fPriority);
105 file.PopValue("tree-name",fTreeName);
106 file.PopValue("tree-comment",fTreeComment);
107 file.PopValue("prefix",fPrefix);
108}
109
110
111/**
112 * Calculate one corrected output by copying the dependent variable and
113 * applying sensitivity-weighted corrections from independent variables.
114 */
116 vector< const VQwHardwareChannel* > &ivs,
117 vector< Double_t > &sens) {
118
119 // if second is NULL, can't do corrector
120 if (output == NULL){
121 QwError<<"Second is value is NULL, unable to calculate corrector."<<QwLog::endl;
122 return;
123 }
124 // For correct type (asym, diff, mps)
125 // if (fDependentType.at(dv) != type) continue;
126
127 // Clear data in second, if first is NULL
128 if (dv == NULL){
129 output->ClearEventData();
130 }else{
131 // Update second value
132 output->AssignValueFrom(dv);
133 }
134
135 // Add corrections
136 for (size_t iv = 0; iv < ivs.size(); iv++) {
137 output->ScaledAdd(sens.at(iv), ivs.at(iv));
138 }
139
140}
141
142/** Copy dependent variables to output variables (default processing). */
144
145 for (size_t i = 0; i < fDependentVar.size(); ++i) {
146 *(fOutputVar.at(i)) = *(fDependentVar[i]);
147 }
148 for (size_t i = 0; i < fDependentValues.size(); ++i) {
150 }
151
152}
153
154
155/**
156 * Connect to external pointers for dependent variables from asym/diff
157 * subsystem arrays, creating corresponding output variables as clones.
158 */
161
162 /// Fill vector of pointers to the relevant data elements
163 for (size_t dv = 0; dv < fDependentName.size(); dv++) {
164 // Get the dependent variables
165 const VQwHardwareChannel* dv_ptr = 0;
166 VQwHardwareChannel* new_ptr = NULL;
167 string name = "";
168 string cor = "cor_";
169
170 if (fDependentType.at(dv)==kHandleTypeMps) {
171 // Quietly ignore the MPS type when we're connecting the asym & diff
172 continue;
173 } else if(fDependentName.at(dv).at(0) == '@' ) {
174 name = fDependentName.at(dv).substr(1,fDependentName.at(dv).length());
175 }else{
176 dv_ptr = this->RequestExternalPointer(fDependentFull.at(dv));
177 if (dv_ptr == NULL){
178 switch (fDependentType.at(dv)) {
179 case kHandleTypeAsym:
180 dv_ptr = asym.RequestExternalPointer(fDependentName.at(dv));
181 break;
182 case kHandleTypeDiff:
183 dv_ptr = diff.RequestExternalPointer(fDependentName.at(dv));
184 break;
185 default:
186 QwWarning << "QwCombiner::ConnectChannels(QwSubsystemArrayParity& asym, QwSubsystemArrayParity& diff): Dependent variable, "
187 << fDependentName.at(dv)
188 << ", for asym/diff processing does not have proper type, type=="
189 << fDependentType.at(dv) << "."<< QwLog::endl;
190 break;
191 }
192 }
193 if (dv_ptr == NULL){
194 QwWarning << "VQwDataHandler::ConnectChannels(QwSubsystemArrayParity& asym, QwSubsystemArrayParity& diff): Dependent variable, "
195 << fDependentName.at(dv)
196 << ", was not found (fullname=="
197 << fDependentFull.at(dv)<< ")." << QwLog::endl;
198 continue;
199 }
200
201 name = dv_ptr->GetElementName().Data();
202 name.insert(0, cor);
203 new_ptr = dv_ptr->Clone(VQwDataElement::kDerived);
204 new_ptr->SetElementName(name);
205 new_ptr->SetSubsystemName(fName);
206 }
207
208 // alias
209 if(fDependentName.at(dv).at(0) == '@') {
210 //QwMessage << "dv: " << name << QwLog::endl;
211 new_ptr = new QwVQWK_Channel(name, VQwDataElement::kDerived);
212 new_ptr->SetSubsystemName(fName);
213 }
214 // defined type
215 else if(dv_ptr!=NULL) {
216 //QwMessage << "dv: " << fDependentName.at(dv) << QwLog::endl;
217 }else {
218 QwWarning << "Dependent variable " << fDependentName.at(dv) << " could not be found, "
219 << "or is not a VQWK channel." << QwLog::endl;
220 continue;
221 }
222
223 // pair creation
224 if(new_ptr != NULL) {
225 fDependentVar.push_back(dv_ptr);
226 fOutputVar.push_back(new_ptr);
227 }
228 }
229 return 0;
230}
231
232
233/**
234 * Parse a variable string to extract type (asym/diff/yield/mps) and name.
235 *
236 * @return Pair of handle type and variable name.
237 */
238pair<VQwDataHandler::EQwHandleType,string> VQwDataHandler::ParseHandledVariable(const string& variable) {
239
240 pair<EQwHandleType,string> type_name;
241 size_t len = variable.length();
242 size_t pos1 = variable.find_first_of(ParseSeparator);
243 size_t pos2 = variable.find_first_not_of(ParseSeparator,pos1);
244 if (pos1 == string::npos) {
245 type_name.first = kHandleTypeUnknown;
246 type_name.second = variable;
247 } else {
248 string type = variable.substr(0,pos1);
249 string name = variable.substr(pos2,len-pos2);
250 if (type == "asym")
251 {type_name.first = kHandleTypeAsym;}
252 else if (type == "diff")
253 {type_name.first = kHandleTypeDiff;}
254 else if (type == "yield")
255 {type_name.first = kHandleTypeYield;}
256 else if (type == "mps")
257 {type_name.first = kHandleTypeMps;}
258 else
259 {type_name.first = kHandleTypeUnknown;}
260 type_name.second = name;
261 }
262 return type_name;
263
264}
265
266/**
267 * Construct TTree branches for output variables, using running sum if
268 * configured for statistics trees.
269 */
271 QwRootFile *treerootfile,
272 const std::string& treeprefix,
273 const std::string& branchprefix)
274{
275 if (fTreeName.size() > 0) {
276 if (fOutputVar.size() == 0) {
277 QwWarning << "No data handler output; not creating tree "
278 << treeprefix + fTreeName
279 << QwLog::endl;
280 } else {
281 TString tmp_branchprefix(branchprefix.c_str());
282 if (tmp_branchprefix.Contains("stat") && fKeepRunningSum
283 && fRunningsum!=NULL){
284 fRunningsumFillsTree = kTRUE;
285 } else {
286 fRunningsumFillsTree = kFALSE;
287 }
288 fTreeName = treeprefix+fTreeName;
290 treerootfile->ConstructTreeBranches(fTreeName, fTreeComment, *fRunningsum, fPrefix+branchprefix);
291 }else {
292 treerootfile->ConstructTreeBranches(fTreeName, fTreeComment, *this, fPrefix+branchprefix);
293 }
294 }
295 }
296}
297
299 TTree *tree,
300 TString& prefix,
302{
303 for (size_t i = 0; i < fOutputVar.size(); ++i) {
304 fOutputVar.at(i)->ConstructBranchAndVector(tree, prefix, values);
305 }
306}
307
309{
310 if (fTreeName.size()>0){
312 treerootfile->FillTreeBranches(*fRunningsum);
313 } else {
314 treerootfile->FillTreeBranches(*this);
315 }
316 treerootfile->FillTree(fTreeName);
317 }
318}
319
321 QwRootFile *treerootfile,
322 const std::string& treeprefix,
323 const std::string& branchprefix)
324{
325 if (fTreeName.size() > 0) {
326 if (fOutputVar.size() == 0) {
327 QwWarning << "No data handler output; not creating RNTuple "
328 << treeprefix + fTreeName
329 << QwLog::endl;
330 } else {
331 TString tmp_branchprefix(branchprefix.c_str());
332 if (tmp_branchprefix.Contains("stat") && fKeepRunningSum
333 && fRunningsum!=NULL){
334 fRunningsumFillsTree = kTRUE;
335 } else {
336 fRunningsumFillsTree = kFALSE;
337 }
338 fTreeName = treeprefix+fTreeName;
339#ifdef HAS_RNTUPLE_SUPPORT
341 treerootfile->ConstructNTupleFields(fTreeName, fTreeComment, *fRunningsum, fPrefix+branchprefix);
342 }else {
343 treerootfile->ConstructNTupleFields(fTreeName, fTreeComment, *this, fPrefix+branchprefix);
344 }
345#endif
346 }
347 }
348}
349
351{
352 if (fTreeName.size()>0){
353#ifdef HAS_RNTUPLE_SUPPORT
355 treerootfile->FillNTupleFields(*fRunningsum);
356 } else {
357 treerootfile->FillNTupleFields(*this);
358 }
359#endif
360 }
361}
362
363
364/**
365 * Fill tree vector with current output variable values.
366 */
368{
369 // Fill the data element
370 for (size_t i = 0; i < fOutputVar.size(); ++i) {
371 if (fOutputVar.at(i) == NULL) {continue;}
372 fOutputVar.at(i)->FillTreeVector(values);
373 }
374}
375
376#ifdef HAS_RNTUPLE_SUPPORT
377void VQwDataHandler::ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs)
378{
379 for (size_t i = 0; i < fOutputVar.size(); ++i) {
380 fOutputVar.at(i)->ConstructNTupleAndVector(model, prefix, values, fieldPtrs);
381 }
382}
383
384void VQwDataHandler::FillNTupleVector(std::vector<Double_t>& values) const
385{
386 // Fill the data element
387 for (size_t i = 0; i < fOutputVar.size(); ++i) {
388 if (fOutputVar.at(i) == NULL) {continue;}
389 fOutputVar.at(i)->FillNTupleVector(values);
390 }
391}
392#endif // HAS_RNTUPLE_SUPPORT
393
394/** Initialize running sum accumulator as a clone of this handler. */
396{
397 if (fKeepRunningSum && fRunningsum == NULL){
398 fRunningsum = this->Clone();
399 fRunningsum->fKeepRunningSum = kFALSE;
400 fRunningsum->ClearEventData();
401 }
402}
403
404/** Accumulate current event into running sum if no error flags are set. */
406{
407 if (fKeepRunningSum && fErrorFlagPtr!=NULL && (*fErrorFlagPtr)==0){
408 fRunningsum->AccumulateRunningSum(*this);
409 }
410}
411
412void VQwDataHandler::AccumulateRunningSum(VQwDataHandler &value, Int_t count, Int_t ErrorMask)
413{
414 for (size_t i = 0; i < fOutputVar.size(); i++){
415 this->fOutputVar[i]->AccumulateRunningSum(value.fOutputVar[i], count, ErrorMask);
416 }
417}
418
419
421{
422 if (fKeepRunningSum && (fRunningsum != NULL)){
423 for(size_t i = 0; i < fRunningsum->fOutputVar.size(); i++) {
424 // calling CalculateRunningAverage in scope of VQwHardwareChannel
425 fRunningsum->fOutputVar[i]->CalculateRunningAverage();
426 }
427 }
428}
429
431{
432 QwMessage<<"=== "<< fName << " ==="<<QwLog::endl<<QwLog::endl;
433 for(size_t i = 0; i < fOutputVar.size(); i++) {
434 fOutputVar[i]->PrintValue();
435 }
436}
437
438
440{
441 for(size_t i = 0; i < fOutputVar.size(); i++) {
442 fOutputVar[i]->ClearEventData();
443 }
444 if (fKeepRunningSum && fRunningsum!=NULL){
445 fRunningsum->ClearEventData();
446 }
447}
448
449
450//*****************************************************************//
452{
453 // Publish variables
454 Bool_t status = kTRUE;
455 VQwHardwareChannel* tmp_channel = nullptr;
456
457 // Publish variables through map file
458 for (size_t pp = 0; pp < fPublishList.size(); pp++) {
459 TString publish_name = fPublishList.at(pp).at(0);
460 TString device_type = fPublishList.at(pp).at(1);
461 TString device_name = fPublishList.at(pp).at(2);
462 TString device_prop = fPublishList.at(pp).at(3);
463 device_type.ToLower();
464 device_prop.ToLower();
465
466 tmp_channel = NULL;
467
468 for(size_t i=0;i<fOutputVar.size(); ++i) {
469 if(device_name.CompareTo(fOutputVar.at(i)->GetElementName())==0){
470 tmp_channel = fOutputVar.at(i);
471 break;
472 }
473 }
474 if (tmp_channel == NULL) {
475 QwError << "VQwDataHandler::PublishInternalValues(): " << publish_name
476 << " not found" << QwLog::endl;
477 status &= kFALSE;
478 } else {
479 QwDebug << "VQwDataHandler::PublishInternalValues(): " << publish_name
480 << " found" << QwLog::endl;
481 status &= PublishInternalValue(publish_name, "published-value", tmp_channel);
482 }
483 }
484 return status;
485}
486
487Bool_t VQwDataHandler::PublishByRequest(TString device_name)
488{
489 Bool_t status = kFALSE;
490 for(size_t i=0;i<fOutputVar.size(); ++i) {
491 if(device_name.CompareTo(fOutputVar.at(i)->GetElementName())==0){
492 status = PublishInternalValue(device_name, "published-by-request",
493 fOutputVar.at(i));
494 break;
495 }
496 }
497 if (!status && fOutputVar.size()>0)
498 QwDebug << "VQwDataHandler::PublishByRequest: Failed to publish channel name: " << device_name << QwLog::endl;
499 return status;
500}
501
503{
504 Bool_t local_print_flag = false;
505 Bool_t local_add_element= type.Contains("asy");
506
507
508 if(local_print_flag){
509 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
510 QwMessage << " QwDataHandlerArray::WritePromptSummary() " << QwLog::endl;
511 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
512 }
513
514 const VQwHardwareChannel* tmp_channel = 0;
515 TString element_name = "";
516 Double_t element_value = 0.0;
517 Double_t element_value_err = 0.0;
518 Double_t element_value_width = 0.0;
519
520 PromptSummaryElement *local_ps_element = NULL;
521 Bool_t local_add_these_elements= false;
522
523 for (size_t i = 0; i < fOutputVar.size(); i++)
524 {
525 element_name = fOutputVar[i]->GetElementName();
526 tmp_channel = fOutputVar[i];
527 element_value = 0.0;
528 element_value_err = 0.0;
529 element_value_width = 0.0;
530
531
532 local_add_these_elements=element_name.Contains("dd")||element_name.Contains("da"); // Need to change this to add other detectors in summary
533
534 if(local_add_these_elements && local_add_element){
535 ps->AddElement(new PromptSummaryElement(element_name));
536 }
537
538 local_ps_element=ps->GetElementByName(element_name);
539
540 if(local_ps_element) {
541 element_value = tmp_channel->GetValue();
542 element_value_err = tmp_channel->GetValueError();
543 element_value_width = tmp_channel->GetValueWidth();
544
545 local_ps_element->Set(type, element_value, element_value_err, element_value_width);
546 }
547
548 if( local_print_flag && local_ps_element) {
549 printf("Type %12s, Element %32s, value %12.4e error %8.4e width %12.4e\n", type.Data(), element_name.Data(), element_value, element_value_err, element_value_width);
550 }
551 }
552
553}
554
555
556#ifdef __USE_DATABASE__
557void VQwDataHandler::FillDB(QwParityDB *db, TString datatype)
558{
559
560 Bool_t local_print_flag = kTRUE;
561
562 UInt_t analysis_id = db->GetAnalysisID();
563
564 TString measurement_type;
565 measurement_type = QwDBInterface::DetermineMeasurementTypeID(datatype);
566
567 std::vector<QwDBInterface> interface;
568
569 std::vector<QwParitySchema::beam_row> beamlist;
570 std::vector<QwParitySchema::md_data_row> mdlist;
571 std::vector<QwParitySchema::lumi_data_row> lumilist;
572
574
575
576 for(size_t i = 0; i < fOutputVar.size(); i++) {
577 interface.clear();
578 fOutputVar[i]->AddEntriesToList(interface);
579 for(size_t j=0; j<interface.size(); j++) {
580 interface.at(j).SetAnalysisID( analysis_id ) ;
581 interface.at(j).SetMeasurementTypeID( measurement_type );
582 tabletype = interface.at(j).SetDetectorID( db );
583 if (tabletype==QwDBInterface::kQwDBI_OtherTable){
584 TString tmp_name = interface.at(j).GetDeviceName();
585 tmp_name.Remove(0,5);
586 interface.at(j).SetDetectorName(tmp_name);
587 tabletype = interface.at(j).SetDetectorID( db );
588 }
589 if (tabletype==QwDBInterface::kQwDBI_BeamTable){
590 interface.at(j).AddThisEntryToList( beamlist );
591 } else if (tabletype==QwDBInterface::kQwDBI_MDTable){
592 interface.at(j).AddThisEntryToList( mdlist );
593 } else if (tabletype==QwDBInterface::kQwDBI_LumiTable){
594 interface.at(j).AddThisEntryToList( lumilist );
595 } else {
596 QwError << "QwCombiner::FillDB: Unrecognized detector name: "
597 << interface.at(j).GetDeviceName() << QwLog::endl;
598 }
599 interface.at(j).PrintStatus( local_print_flag);
600 }
601 }
602
603 // Database operations with scoped connection
604 {
605 auto c = db->GetScopedConnection();
606
607 // Check the entrylist size, if it isn't zero, start to query..
608 if( beamlist.size() ) {
609 for (const auto& entry: beamlist) {
610 c->QueryExecute(entry.insert_into());
611 }
612 } else {
613 QwMessage << "QwCombiner::FillDB :: This is the case when the beamlist contains nothing for type="<< measurement_type.Data()
614 << QwLog::endl;
615 }
616 if( mdlist.size() ) {
617 for (const auto& entry: mdlist) {
618 c->QueryExecute(entry.insert_into());
619 }
620 } else {
621 QwMessage << "QwCombiner::FillDB :: This is the case when the mdlist contains nothing for type="<< measurement_type.Data()
622 << QwLog::endl;
623 }
624 if( lumilist.size() ) {
625 for (const auto& entry: lumilist) {
626 c->QueryExecute(entry.insert_into());
627 }
628 } else {
629 QwMessage << "QwCombiner::FillDB :: This is the case when the lumilist contains nothing for type="<< measurement_type.Data()
630 << QwLog::endl;
631 }
632 }
633 return;
634}
635#endif // __USE_DATABASE__
636
#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
#define QwDebug
Predefined log drain for debugging output.
Definition QwLog.h:59
ROOT file and tree management wrapper classes.
Parameter file parsing and management.
Prompt summary data management.
Decoding and management for VQWK ADC channels (6x32-bit datawords)
Virtual base class for data handlers accessing multiple subsystems.
Bool_t PublishInternalValue(const TString name, const TString desc, const VQwHardwareChannel *element) const
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 TString DetermineMeasurementTypeID(TString type, TString suffix="", Bool_t forcediffs=kFALSE)
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)
std::string GetLine()
void Set(TString type, const Double_t a, const Double_t a_err, const Double_t a_width)
void AddElement(PromptSummaryElement *in)
PromptSummaryElement * GetElementByName(TString name)
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:53
A wrapper class for a ROOT file or memory mapped file.
Definition QwRootFile.h:827
Int_t FillTree(const std::string &name)
Fill the tree with name.
Definition QwRootFile.h:949
void ConstructTreeBranches(const std::string &name, const std::string &desc, T &object, const std::string &prefix="")
Construct the tree branches of a generic object.
void FillTreeBranches(const std::string &name, const T &object)
Fill the tree branches of a generic object by tree name.
Concrete hardware channel for VQWK ADC modules (6x32-bit words)
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.
Abstract base for concrete hardware channels implementing dual-operator pattern.
Double_t GetValueError() const
virtual void ScaledAdd(Double_t scale, const VQwHardwareChannel *value)=0
virtual VQwHardwareChannel * Clone() const
void ClearEventData() override
Clear the event data in this element.
Double_t GetValueWidth() const
Double_t GetValue() const
void AssignValueFrom(const VQwDataElement *valueptr) override=0
Subsystem array container specialized for parity analysis with asymmetry calculations.
const UInt_t * GetEventcutErrorFlagPointer() const
void CalcOneOutput(const VQwHardwareChannel *dv, VQwHardwareChannel *output, std::vector< const VQwHardwareChannel * > &ivs, std::vector< Double_t > &sens)
virtual void ProcessData()
QwSubsystemArrayParity * fSubsystemArray
Single event pointer.
void FillTreeVector(QwRootTreeBranchVector &values) const
std::vector< std::string > fDependentFull
void AccumulateRunningSum()
void SetEventcutErrorFlagPointer(const UInt_t *errorflagptr)
void WritePromptSummary(QwPromptSummary *ps, TString type)
std::vector< std::vector< TString > > fPublishList
void PrintValue() const
std::string fPrefix
Bool_t PublishInternalValues() const override
Publish all variables of the subsystem.
virtual void ConstructTreeBranches(QwRootFile *treerootfile, const std::string &treeprefix="", const std::string &branchprefix="")
std::vector< Double_t > fDependentValues
virtual void ParseConfigFile(QwParameterFile &file)
virtual void ClearEventData()
std::vector< const VQwHardwareChannel * > fDependentVar
~VQwDataHandler() override
void CalculateRunningAverage()
virtual void FillTreeBranches(QwRootFile *treerootfile)
VQwDataHandler(const TString &name)
VQwDataHandler * fRunningsum
const UInt_t * fErrorFlagPtr
Error flag pointer.
virtual void FillNTupleFields(QwRootFile *treerootfile)
std::string ParseSeparator
virtual void ConstructNTupleFields(QwRootFile *treerootfile, const std::string &treeprefix="", const std::string &branchprefix="")
RNTuple methods.
Short_t fBurstCounter
When a datahandler array is processed, handlers with lower priority will be processed before handlers...
std::string fTreeComment
std::vector< VQwHardwareChannel * > fOutputVar
std::string fTreeName
Bool_t PublishByRequest(TString device_name) override
Try to publish an internal variable matching the submitted name.
Bool_t fRunningsumFillsTree
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values)
std::pair< EQwHandleType, std::string > ParseHandledVariable(const std::string &variable)
std::vector< EQwHandleType > fDependentType
QwHelicityPattern * fHelicityPattern
Helicity pattern pointer.
std::vector< Double_t > fOutputValues
std::vector< std::string > fDependentName
virtual Int_t ConnectChannels(QwSubsystemArrayParity &, QwSubsystemArrayParity &asym, QwSubsystemArrayParity &diff)
std::string fMapFile