JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
VQwDetectorArray.cc
Go to the documentation of this file.
1/*!
2 * \file VQwDetectorArray.cc
3 * \brief Virtual base class implementation for detector arrays managing PMT collections
4 *
5 * Base detector array implementation managing PMT collections (integration
6 * and combined), including channel mapping, event cuts, normalization options,
7 * publishing, tree construction, and running sums. Derived classes implement
8 * specific detector systems. Documentation-only edits; runtime behavior unchanged.
9 */
10
11#include "VQwDetectorArray.h"
12
13// System headers
14#include <sstream>
15
16#ifdef HAS_RNTUPLE_SUPPORT
17// ROOT headers for RNTuple support
18#include <ROOT/RNTupleModel.hxx>
19#include <ROOT/RNTupleWriter.hxx>
20#endif
21
22// Qweak headers
23#include "QwSubsystemArray.h"
24#include "QwLog.h"
25#ifdef __USE_DATABASE__
26#include "QwParitySchema.h"
27#include "QwParityDB.h"
28#endif
29#include "QwPromptSummary.h"
30
31/**
32 * Define command-line options for detector array normalization.
33 *
34 * @param options Options object to configure.
35 */
37 // Define the execution options
38 options.AddOptions()
39 ("QwDetectorArray.normalize",
40 po::value<bool>()->default_bool_value(true),
41 "Normalize the detectors by beam current");
42
43 options.AddOptions()
44 ("QwDetectorArray.norm_threshold",
45 po::value<double>()->default_value(2.5),
46 "Normalize the detectors for currents above this value");
47}
48
49/**
50 * Load detector array configuration from parsed command-line options.
51 *
52 * @param options Options object.
53 */
55
56 bNormalization = options.GetValue<bool>("QwDetectorArray.normalize");
57
58 if (! bNormalization) {
59
60 QwWarning << "QwDetectorArray::ProcessOptions: "
61 << "Detector yields WILL NOT be normalized."
62 << QwLog::endl;
63
64 }
65
66 fNormThreshold = options.GetValue<double>("QwDetectorArray.norm_threshold");
67
68}
69
70
71//*****************************************************************//
72/**
73 * Publish internal detector channels according to the configured
74 * publish list (integration and combined PMTs).
75 *
76 * @return true if all requested channels are successfully published.
77 */
79
80 // Publish variables
81
82 Bool_t status = kTRUE;
83
84/*
85 status = status && PublishInternalValue("qwk_md1neg", "qwk_md1neg", GetIntegrationPMT("qwk_md1neg")->GetChannel("qwk_md1neg"));
86 status = status && PublishInternalValue("qwk_md1pos", "qwk_md1pos", GetIntegrationPMT("qwk_md1pos")->GetChannel("qwk_md1pos"));
87 status = status && PublishInternalValue("qwk_md2neg", "qwk_md2neg", GetIntegrationPMT("qwk_md2neg")->GetChannel("qwk_md2neg"));
88 status = status && PublishInternalValue("qwk_md2pos", "qwk_md2pos", GetIntegrationPMT("qwk_md2pos")->GetChannel("qwk_md2pos"));
89 status = status && PublishInternalValue("qwk_md3neg", "qwk_md3neg", GetIntegrationPMT("qwk_md3neg")->GetChannel("qwk_md3neg"));
90 status = status && PublishInternalValue("qwk_md3pos", "qwk_md3pos", GetIntegrationPMT("qwk_md3pos")->GetChannel("qwk_md3pos"));
91 status = status && PublishInternalValue("qwk_md4neg", "qwk_md4neg", GetIntegrationPMT("qwk_md4neg")->GetChannel("qwk_md4neg"));
92 status = status && PublishInternalValue("qwk_md4pos", "qwk_md4pos", GetIntegrationPMT("qwk_md4pos")->GetChannel("qwk_md4pos"));
93 status = status && PublishInternalValue("qwk_md5neg", "qwk_md5neg", GetIntegrationPMT("qwk_md5neg")->GetChannel("qwk_md5neg"));
94 status = status && PublishInternalValue("qwk_md5pos", "qwk_md5pos", GetIntegrationPMT("qwk_md5pos")->GetChannel("qwk_md5pos"));
95 status = status && PublishInternalValue("qwk_md6neg", "qwk_md6neg", GetIntegrationPMT("qwk_md6neg")->GetChannel("qwk_md6neg"));
96 status = status && PublishInternalValue("qwk_md6pos", "qwk_md6pos", GetIntegrationPMT("qwk_md6pos")->GetChannel("qwk_md6pos"));
97 status = status && PublishInternalValue("qwk_md7neg", "qwk_md7neg", GetIntegrationPMT("qwk_md7neg")->GetChannel("qwk_md7neg"));
98 status = status && PublishInternalValue("qwk_md7pos", "qwk_md7pos", GetIntegrationPMT("qwk_md7pos")->GetChannel("qwk_md7pos"));
99 status = status && PublishInternalValue("qwk_md8neg", "qwk_md8neg", GetIntegrationPMT("qwk_md8neg")->GetChannel("qwk_md8neg"));
100 status = status && PublishInternalValue("qwk_md8pos", "qwk_md8pos", GetIntegrationPMT("qwk_md8pos")->GetChannel("qwk_md8pos"));
101 status = status && PublishInternalValue("qwk_md9neg", "qwk_md9neg", GetIntegrationPMT("qwk_md9neg")->GetChannel("qwk_md9neg"));
102 status = status && PublishInternalValue("qwk_md9pos", "qwk_md9pos", GetIntegrationPMT("qwk_md9pos")->GetChannel("qwk_md9pos"));
103*/
104
105/*
106 status = status && PublishInternalValue("qwk_md1barsum","qwk_md1barsum", GetCombinedPMT("qwk_md1barsum")->GetChannel("qwk_md1barsum"));
107 status = status && PublishInternalValue("qwk_md2barsum","qwk_md2barsum", GetCombinedPMT("qwk_md2barsum")->GetChannel("qwk_md2barsum"));
108 status = status && PublishInternalValue("qwk_md3barsum","qwk_md3barsum", GetCombinedPMT("qwk_md3barsum")->GetChannel("qwk_md3barsum"));
109 status = status && PublishInternalValue("qwk_md4barsum","qwk_md4barsum", GetCombinedPMT("qwk_md4barsum")->GetChannel("qwk_md4barsum"));
110 status = status && PublishInternalValue("qwk_md5barsum","qwk_md5barsum", GetCombinedPMT("qwk_md5barsum")->GetChannel("qwk_md5barsum"));
111 status = status && PublishInternalValue("qwk_md6barsum","qwk_md6barsum", GetCombinedPMT("qwk_md6barsum")->GetChannel("qwk_md6barsum"));
112 status = status && PublishInternalValue("qwk_md7barsum","qwk_md7barsum", GetCombinedPMT("qwk_md7barsum")->GetChannel("qwk_md7barsum"));
113 status = status && PublishInternalValue("qwk_md8barsum","qwk_md8barsum", GetCombinedPMT("qwk_md8barsum")->GetChannel("qwk_md8barsum"));
114
115 status = status && PublishInternalValue("qwk_mdallbars","qwk_mdallbars", GetCombinedPMT("qwk_mdallbars")->GetChannel("qwk_mdallbars"));
116*/
117
118 //return status;
119
120
121 // TODO:
122 // The variables should be published based on the parameter file.
123 // See QwBeamLine class for an implementation.
124
125 // Publish variables through map file
126 // This should work with bcm, bpmstripline, bpmcavity, combo bpm and combo bcm
127
128 for (size_t pp = 0; pp < fPublishList.size(); pp++) {
129
130 TString publish_name = fPublishList.at(pp).at(0);
131 TString device_type = fPublishList.at(pp).at(1);
132 TString device_name = fPublishList.at(pp).at(2);
133 TString device_prop = fPublishList.at(pp).at(3);
134 device_type.ToLower();
135 device_prop.ToLower();
136
137 const VQwHardwareChannel* tmp_channel = NULL;
138
139 if (device_type == "integrationpmt") {
140
141 tmp_channel = GetIntegrationPMT(device_name)->GetChannel(device_name);
142
143 } else if (device_type == "combinedpmt") {
144
145 tmp_channel = GetCombinedPMT(device_name)->GetChannel(device_name);
146
147 } else {
148
149 QwError << "QwBeamLine::PublishInternalValues() error "<< QwLog::endl;
150
151 }
152
153 if (tmp_channel == NULL) {
154
155 QwError << "QwBeamLine::PublishInternalValues(): " << publish_name << " not found" << QwLog::endl;
156 status |= kFALSE;
157
158 } else {
159
160 QwDebug << "QwBeamLine::PublishInternalValues(): " << publish_name << " found" << QwLog::endl;
161
162 }
163
164 status = status && PublishInternalValue(publish_name, publish_name, tmp_channel);
165
166 }
167
168 return status;
169
170}
171
172
173/**
174 * Publish a specific device channel on-demand by name lookup.
175 *
176 * @param device_name Name of the detector channel to publish.
177 * @return true if the channel is found and published successfully.
178 */
179Bool_t VQwDetectorArray::PublishByRequest(TString device_name) {
180
181 Bool_t status = kFALSE;
182 //std::cerr << "##### device_name==\"" << device_name << "\"" << std::endl;
183
184 for(size_t i=0;i<fMainDetID.size();i++) {
185
186 //std::cerr << "fMainDetID[i].fdetectorname==\"" << fMainDetID[i].fdetectorname << "\"" << std::endl;
187
188 if(device_name.CompareTo(fMainDetID[i].fdetectorname)!=0) continue;
189
190 if (fMainDetID[i].fTypeID == kQwCombinedPMT){
191
192 status = PublishInternalValue(device_name, "published-by-request",
193 fCombinedPMT[fMainDetID[i].fIndex].GetChannel(device_name));
194
195 } else if (fMainDetID[i].fTypeID == kQwIntegrationPMT) {
196
197 status = PublishInternalValue(device_name, "published-by-request",
198 fIntegrationPMT[fMainDetID[i].fIndex].GetChannel(device_name));
199
200 } else {
201
202 QwError << "Unknown channel name: " << device_name << QwLog::endl;
203
204 }
205
206 break;
207
208 }
209
210 if (!status)
211 QwDebug << "VQwDetectorArray::PublishByRequest: Failed to publish channel name: " << device_name << QwLog::endl;
212
213 return status;
214
215}
216
217
218//*****************************************************************//
219/**
220 * Load detector channel map file, creating integration and combined PMTs
221 * and configuring buffer layout, saturation limits, and sample sizes.
222 *
223 * @param mapfile Path to the channel map file.
224 * @return 0 on success.
225 */
226Int_t VQwDetectorArray::LoadChannelMap(TString mapfile) {
227
228 Bool_t ldebug=kFALSE;
229
230 std::vector<TString> combinedchannelnames;
231 std::vector<Double_t> weight;
232 Int_t wordsofar=0;
233 Int_t currentsubbankindex=-1;
234 Int_t sample_size=0;
235 Double_t abs_saturation_limit = 8.5; // default saturation limit(volt)
236 Bool_t bAssignedLimit = kFALSE;
237
238 // Open the file
239 QwParameterFile mapstr(mapfile.Data());
240 TString varname, varvalue;
241
243 mapstr.EnableGreediness();
244 mapstr.SetCommentChars("!");
245
246 UInt_t value;
247 size_t vqwk_buffer_offset = 0;
248
249 while (mapstr.ReadNextLine()) {
250
251 RegisterRocBankMarker(mapstr);
252 if (mapstr.PopValue("abs_saturation_limit",value)) {
253 abs_saturation_limit=value;
254 bAssignedLimit = kTRUE;
255 }
256
257 if (mapstr.PopValue("sample_size",value)) {
258 sample_size=value;
259 }
260
261 if (mapstr.PopValue("vqwk_buffer_offset",value)) {
262 vqwk_buffer_offset=value;
263 }
264
265 mapstr.TrimComment('!'); // Remove everything after a '!' character.
266 mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
267
268 if (mapstr.LineIsEmpty()) continue;
269
270 Bool_t lineok = kTRUE;
271 TString keyword = "";
272 TString keyword2 = "";
273 TString modtype = "";
274 TString dettype = "";
275 TString namech = "";
276 Int_t modnum = 0;
277 Int_t channum = 0;
278
279 modtype = mapstr.GetTypedNextToken<TString>(); // module type
280
281 modtype.ToUpper();
282
283 if (modtype == "VPMT") {
284
285 channum = mapstr.GetTypedNextToken<Int_t>(); //channel number
286 Int_t combinedchans = mapstr.GetTypedNextToken<Int_t>(); //number of combined channels
287 dettype = mapstr.GetTypedNextToken<TString>(); //type-purpose of the detector
288 dettype.ToLower();
289 namech = mapstr.GetTypedNextToken<TString>(); //name of the detector
290 namech.ToLower();
291 combinedchannelnames.clear();
292
293 for (int i=0; i<combinedchans; i++){
294
295 TString nameofcombinedchan = mapstr.GetTypedNextToken<TString>();
296 nameofcombinedchan.ToLower();
297 combinedchannelnames.push_back(nameofcombinedchan);
298 }
299
300 weight.clear();
301
302 for (int i=0; i<combinedchans; i++) {
303
304 weight.push_back( mapstr.GetTypedNextToken<Double_t>());
305 }
306
307 keyword = mapstr.GetTypedNextToken<TString>();
308 keyword.ToLower();
309 keyword2 = mapstr.GetTypedNextToken<TString>();
310 keyword2.ToLower();
311
312 } else {
313
314 modnum = mapstr.GetTypedNextToken<Int_t>(); //slot number
315 channum = mapstr.GetTypedNextToken<Int_t>(); //channel number
316 dettype = mapstr.GetTypedNextToken<TString>(); //type-purpose of the detector
317 dettype.ToLower();
318 namech = mapstr.GetTypedNextToken<TString>(); //name of the detector
319 namech.ToLower();
320
321 keyword = mapstr.GetTypedNextToken<TString>();
322 keyword.ToLower();
323 keyword2 = mapstr.GetTypedNextToken<TString>();
324 keyword2.ToLower();
325 }
326
327
328 if (currentsubbankindex!=GetSubbankIndex(fCurrentROC_ID,fCurrentBank_ID)) {
329
331 wordsofar=0;
332 }
333
334 QwDetectorArrayID localMainDetID;
335 localMainDetID.fdetectorname=namech;
336 localMainDetID.fmoduletype=modtype;
337 localMainDetID.fSubbankIndex=currentsubbankindex;
338 localMainDetID.fdetectortype=dettype;
339
340 //localMainDetID.fWordInSubbank=wordsofar;
341
342 if (modtype=="MOLLERADC") {
343
344 Int_t offset = QwMollerADC_Channel::GetBufferOffset(modnum, channum)+vqwk_buffer_offset;
345
346 if (offset>=0){
347
348 localMainDetID.fWordInSubbank = wordsofar + offset;
349 }
350
351 } else if (modtype=="VPMT") {
352
353 localMainDetID.fCombinedChannelNames = combinedchannelnames;
354 localMainDetID.fWeight = weight;
355
356 //std::cout<<"Add in a combined channel"<<std::endl;
357 } else {
358
359 QwError << "VQwDetectorArray::LoadChannelMap: Unknown module type: "
360 << modtype <<", the detector "<<namech<<" will not be decoded "
361 << QwLog::endl;
362 lineok=kFALSE;
363 continue;
364 }
365
366 localMainDetID.fTypeID=GetDetectorTypeID(dettype);
367
368 if (localMainDetID.fTypeID==kQwUnknownPMT) {
369
370 QwError << "VQwDetectorArray::LoadChannelMap: Unknown detector type: "
371 << dettype <<", the detector "<<namech<<" will not be decoded "
372 << QwLog::endl;
373 lineok=kFALSE;
374 continue;
375 }
376
377 localMainDetID.fIndex= GetDetectorIndex(localMainDetID.fTypeID,
378 localMainDetID.fdetectorname);
379
380 if (localMainDetID.fIndex==-1){
381
382 if (localMainDetID.fTypeID==kQwIntegrationPMT){
383
384 QwIntegrationPMT localIntegrationPMT(GetName(),localMainDetID.fdetectorname);
385
386 if (keyword=="not_blindable" || keyword2=="not_blindable")
387 localIntegrationPMT.SetBlindability(kFALSE);
388
389 else
390 localIntegrationPMT.SetBlindability(kTRUE);
391
392 if (keyword=="not_normalizable" || keyword2=="not_normalizable")
393 localIntegrationPMT.SetNormalizability(kFALSE);
394
395 else
396 localIntegrationPMT.SetNormalizability(kTRUE);
397
398 fIntegrationPMT.push_back(localIntegrationPMT);
399 fIntegrationPMT[fIntegrationPMT.size()-1].SetDefaultSampleSize(sample_size);
400
401 if(bAssignedLimit)
402 fIntegrationPMT[fIntegrationPMT.size()-1].SetSaturationLimit(abs_saturation_limit);
403
404 localMainDetID.fIndex=fIntegrationPMT.size()-1;
405
406 } else if (localMainDetID.fTypeID==kQwCombinedPMT) {
407
408 QwCombinedPMT localcombinedPMT(GetName(),localMainDetID.fdetectorname);
409
410 if (keyword=="not_normalizable" || keyword2=="not_normalizable")
411 localcombinedPMT.SetNormalizability(kFALSE);
412
413 else
414 localcombinedPMT.SetNormalizability(kTRUE);
415
416 if (keyword=="not_blindable" || keyword2 =="not_blindable")
417 localcombinedPMT.SetBlindability(kFALSE);
418
419 else
420 localcombinedPMT.SetBlindability(kTRUE);
421
422 fCombinedPMT.push_back(localcombinedPMT);
423 fCombinedPMT[fCombinedPMT.size()-1].SetDefaultSampleSize(sample_size);
424 localMainDetID.fIndex=fCombinedPMT.size()-1;
425 }
426 }
427
428 if (ldebug) {
429
430 localMainDetID.Print();
431 std::cout<<"line ok=";
432
433 if (lineok)
434 std::cout<<"TRUE"<<std::endl;
435
436 else
437 std::cout<<"FALSE"<<std::endl;
438 }
439
440 if (lineok)
441 fMainDetID.push_back(localMainDetID);
442
443 } // End of "while (mapstr.ReadNextLine())"
444
445 for (size_t i=0; i<fMainDetID.size(); i++) {
446
447 if (fMainDetID[i].fTypeID==kQwCombinedPMT) {
448
449 Int_t ind = fMainDetID[i].fIndex;
450
451 //check to see if all required channels are available
452 if (ldebug) {
453
454 std::cout<<"fMainDetID[i].fCombinedChannelNames.size()="
455 <<fMainDetID[i].fCombinedChannelNames.size()<<std::endl<<"name list: ";
456
457 for (size_t n=0; n<fMainDetID[i].fCombinedChannelNames.size(); n++)
458 std::cout<<" "<<fMainDetID[i].fCombinedChannelNames[n];
459
460 std::cout<<std::endl;
461 }
462
463 Int_t chanmatched=0;
464
465 for (size_t j=0; j<fMainDetID[i].fCombinedChannelNames.size(); j++) {
466
467 for (size_t k=0; k<fMainDetID.size(); k++) {
468
469 if (fMainDetID[i].fCombinedChannelNames[j]==fMainDetID[k].fdetectorname) {
470
471 if (ldebug)
472 std::cout<<"found a to-be-combined channel candidate"<<std::endl;
473
474 chanmatched ++;
475 break;
476 }
477 }
478 }
479
480 if ((Int_t) fMainDetID[i].fCombinedChannelNames.size()==chanmatched) {
481
482 for (size_t l=0; l<fMainDetID[i].fCombinedChannelNames.size(); l++) {
483
484 Int_t ind_pmt = GetDetectorIndex(GetDetectorTypeID("integrationpmt"),
485 fMainDetID[i].fCombinedChannelNames[l]);
486
487 fCombinedPMT[ind].Add(&fIntegrationPMT[ind_pmt],fMainDetID[i].fWeight[l]);
488 }
489
490 fCombinedPMT[ind].LinkChannel(fMainDetID[i].fdetectorname);
491
492 if (ldebug)
493 std::cout<<"linked a combined channel"<<std::endl;
494 } else {
495
496 std::cerr<<"cannot combine void channels for "<<fMainDetID[i].fdetectorname<<std::endl;
497 fMainDetID[i].fIndex = -1;
498 continue;
499 }
500 }
501 }
502
503
504 // Now load the variables to publish
505 mapstr.RewindToFileStart();
506 std::unique_ptr<QwParameterFile> section;
507 std::vector<TString> publishinfo;
508 while ((section = mapstr.ReadNextSection(varvalue))) {
509
510 if (varvalue == "PUBLISH") {
511
512 fPublishList.clear();
513
514 while (section->ReadNextLine()) {
515
516 section->TrimComment(); // Remove everything after a comment character
517 section->TrimWhitespace(); // Get rid of leading and trailing spaces
518
519 for (int ii = 0; ii < 4; ii++) {
520
521 varvalue = section->GetNextToken().c_str();
522
523 if (varvalue.Length()) {
524
525 publishinfo.push_back(varvalue);
526 }
527 }
528
529 if (publishinfo.size() == 4)
530 fPublishList.push_back(publishinfo);
531
532 publishinfo.clear();
533 }
534 }
535 }
536
537 // Print list of variables to publish
538 if (fPublishList.size()>0){
539
540 QwMessage << "Variables to publish:" << QwLog::endl;
541
542 for (size_t jj = 0; jj < fPublishList.size(); jj++)
543 QwMessage << fPublishList.at(jj).at(0) << " " << fPublishList.at(jj).at(1) << " "
544 << fPublishList.at(jj).at(2) << " " << fPublishList.at(jj).at(3) << QwLog::endl;
545 }
546
547 if (ldebug) {
548
549 std::cout<<"Done with Load channel map\n";
550
551 for (size_t i=0;i<fMainDetID.size();i++)
552 if (fMainDetID[i].fIndex>=0)
553 fMainDetID[i].Print();
554 }
555
556 ldebug=kFALSE;
557 mapstr.Close(); // Close the file (ifstream)
558 return 0;
559}
560
561void VQwDetectorArray::LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag) {
562 TString device_type = mapstr.GetTypedNextToken<TString>();
563 device_type.ToLower();
564 TString device_name = mapstr.GetTypedNextToken<TString>();
565 device_name.ToLower();
566
567 Int_t det_index = GetDetectorIndex(GetDetectorTypeID(device_type),device_name);
568 if (det_index == -1) {
569 QwWarning << " Device not found " << device_name << " of type " << device_type << QwLog::endl;
570 //continue;
571 }
572
573 Double_t LLX = mapstr.GetTypedNextToken<Double_t>(); //lower limit for IntegrationPMT value
574 Double_t ULX = mapstr.GetTypedNextToken<Double_t>(); //upper limit for IntegrationPMT value
575 varvalue = mapstr.GetTypedNextToken<TString>();//global/local
576 varvalue.ToLower();
577
578 Double_t burplevel = mapstr.GetTypedNextToken<Double_t>();
579 Double_t stabilitycut = mapstr.GetTypedNextToken<Double_t>();
580
582 QwMessage << "VQwDetectorArray Error Code passing to QwIntegrationPMT " << GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut) << QwLog::endl;
583 fIntegrationPMT[det_index].SetSingleEventCuts(GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut),LLX,ULX,stabilitycut,burplevel);
584 } else if (device_type == GetQwPMTInstrumentTypeName(kQwCombinedPMT)){
585 QwMessage << "VQwDetectorArray Error Code passing to QwCombinedPMT " << GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut) << QwLog::endl;
586 fCombinedPMT[det_index].SetSingleEventCuts(GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut),LLX,ULX,stabilitycut,burplevel);
587 }
588}
589
590void VQwDetectorArray::LoadEventCuts_Fin(Int_t &eventcut_flag) {
591 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
592 fIntegrationPMT[i].SetEventCutMode(eventcut_flag);
593 for (size_t i = 0; i < fCombinedPMT.size(); i++)
594 fCombinedPMT[i].SetEventCutMode(eventcut_flag);
595
596 fMainDetErrorCount = 0;//set the error counter to zero
597}
598
599
600
601Int_t VQwDetectorArray::LoadInputParameters(TString pedestalfile) {
602
603 Bool_t ldebug=kFALSE;
604 TString varname;
605 Double_t varped;
606 Double_t varcal;
607
608 // Double_t varbaserate;
609 Double_t varnormrate;
610 Double_t varvoltperhz;
611 Double_t varasym;
612 Double_t varcx;
613 Double_t varcy;
614 Double_t varcxp;
615 Double_t varcyp;
616 Double_t varce;
617
618 TString localname;
619
620 Int_t lineread=0;
621
622 QwParameterFile mapstr(pedestalfile.Data()); //Open the file
624
625 while (mapstr.ReadNextLine()) {
626
627 lineread+=1;
628 if (ldebug)std::cout<<" line read so far ="<<lineread<<"\n";
629 mapstr.TrimComment('!'); // Remove everything after a '!' character.
630
631 mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
632
633 if (mapstr.LineIsEmpty()) continue;
634
635 else {
636 varname = mapstr.GetTypedNextToken<TString>(); //name of the channel
637 varname.ToLower();
638 varname.Remove(TString::kBoth,' ');
639 varped = mapstr.GetTypedNextToken<Double_t>(); // value of the pedestal
640 varcal = mapstr.GetTypedNextToken<Double_t>(); // value of the calibration factor
641
642
643
644
645 if (ldebug)
646 std::cout << "Inputs for channel " << varname << ": ped=" << varped << ": cal=" << varcal << "\n"
647 << ": varnormrate=" << varnormrate << "\n"
648 << ": varvoltperhz=" << varvoltperhz << "\n"
649 << ": asym=" << varasym << "\n"
650 << ": C_x=" << varcx << ": C_y=" << varcy << "\n"
651 << ": C_xp=" << varcxp << ": C_yp=" << varcyp << "\n"
652 << ": C_e=" << varce << "\n";
653
654 // Bool_t notfound=kTRUE;
655
656 // if (notfound)
657 for (size_t i=0;i<fIntegrationPMT.size();i++){
658 if (fIntegrationPMT[i].GetElementName()==varname) {
659
660 fIntegrationPMT[i].SetPedestal(varped);
661 fIntegrationPMT[i].SetCalibrationFactor(varcal);
662
663
664
665 // i=fIntegrationPMT.size()+1;
666 // notfound=kFALSE;
667 // i=fIntegrationPMT.size()+1;
668
669 break;
670
671 }
672 }
673
674 }
675
676 }
677
678 if (ldebug)
679 std::cout<<" line read in the pedestal + cal file ="<<lineread<<" \n";
680
681 ldebug=kFALSE;
682 mapstr.Close(); // Close the file (ifstream)
683 return 0;
684
685}
686
687
688
689
690
691
692
693
694
696
697 Bool_t ldebug=kFALSE;
698 TString varname;
699 Double_t varped;
700 Double_t varcal;
701
702 // Double_t varbaserate;
703 Double_t varnormrate;
704 Double_t varvoltperhz;
705 Double_t varasym;
706 Double_t varcx;
707 Double_t varcy;
708 Double_t varcxp;
709 Double_t varcyp;
710 Double_t varce;
711
712 TString localname;
713
714 Int_t lineread=0;
715
716 QwParameterFile mapstr(pedestalfile.Data()); //Open the file
718
719 while (mapstr.ReadNextLine()) {
720
721 lineread+=1;
722 if (ldebug)std::cout<<" line read so far ="<<lineread<<"\n";
723 mapstr.TrimComment('!'); // Remove everything after a '!' character.
724
725 mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
726
727 if (mapstr.LineIsEmpty()) continue;
728
729 else {
730 varname = mapstr.GetTypedNextToken<TString>(); //name of the channel
731 varname.ToLower();
732 varname.Remove(TString::kBoth,' ');
733
734
735 varnormrate = mapstr.GetTypedNextToken<Double_t>(); // value of the NormRate
736 varvoltperhz = mapstr.GetTypedNextToken<Double_t>(); // value of the VoltPerHz
737 varasym = mapstr.GetTypedNextToken<Double_t>(); // value of the asymmetry
738 varcx = mapstr.GetTypedNextToken<Double_t>(); // value of the coefficient C_x
739 varcy = mapstr.GetTypedNextToken<Double_t>(); // value of the coefficient C_y
740 varcxp = mapstr.GetTypedNextToken<Double_t>(); // value of the coefficient C_xp
741 varcyp = mapstr.GetTypedNextToken<Double_t>(); // value of the coefficient C_yp
742 varce = mapstr.GetTypedNextToken<Double_t>(); // value of the coefficient C_e
743
744
745 if (ldebug)
746 std::cout << "Inputs for channel " << varname << ": ped=" << varped << ": cal=" << varcal << "\n"
747 << ": varnormrate=" << varnormrate << "\n"
748 << ": varvoltperhz=" << varvoltperhz << "\n"
749 << ": asym=" << varasym << "\n"
750 << ": C_x=" << varcx << ": C_y=" << varcy << "\n"
751 << ": C_xp=" << varcxp << ": C_yp=" << varcyp << "\n"
752 << ": C_e=" << varce << "\n";
753
754 // Bool_t notfound=kTRUE;
755
756 // if (notfound)
757 for (size_t i=0;i<fIntegrationPMT.size();i++){
758 if (fIntegrationPMT[i].GetElementName()==varname) {
759
760
761 fIntegrationPMT[i].SetNormRate(varnormrate);
762 fIntegrationPMT[i].SetVoltPerHz(varvoltperhz);
763 fIntegrationPMT[i].SetAsymmetry(varasym);
764 fIntegrationPMT[i].SetCoefficientCx(varcx);
765 fIntegrationPMT[i].SetCoefficientCy(varcy);
766 fIntegrationPMT[i].SetCoefficientCxp(varcxp);
767 fIntegrationPMT[i].SetCoefficientCyp(varcyp);
768 fIntegrationPMT[i].SetCoefficientCe(varce);
769
770 // i=fIntegrationPMT.size()+1;
771 // notfound=kFALSE;
772 // i=fIntegrationPMT.size()+1;
773
774 break;
775
776 }
777 }
778
779 }
780
781 }
782
783 if (ldebug)
784 std::cout<<" line read in the pedestal + cal file ="<<lineread<<" \n";
785
786 ldebug=kFALSE;
787 mapstr.Close(); // Close the file (ifstream)
788
789}
790
791
792
793
794
795
796
797
798
799
800
801
803
804 Bool_t test=kTRUE;
805 return test;
806
807}
808
810
811 for (size_t i=0;i<fIntegrationPMT.size();i++){
812
813 fIntegrationPMT[i].ClearEventData();
814
815 }
816
817 for (size_t i=0;i<fCombinedPMT.size();i++)
819
820 return;
821
822}
823
824
825/********************************************************/
826
827void VQwDetectorArray::SetRandomEventParameters(Double_t mean, Double_t sigma) {
828
829 for (size_t i = 0; i < fMainDetID.size(); i++) {
830
831 // This is a QwIntegrationPMT
832 if (fMainDetID.at(i).fTypeID == kQwIntegrationPMT)
833 fIntegrationPMT[fMainDetID.at(i).fIndex].SetRandomEventParameters(mean, sigma);
834
835 }
836
837}
838
840
841 for (size_t i = 0; i < fMainDetID.size(); i++) {
842
843 // This is a QwIntegrationPMT
844 if (fMainDetID.at(i).fTypeID == kQwIntegrationPMT)
845 fIntegrationPMT[fMainDetID.at(i).fIndex].SetRandomEventAsymmetry(asymmetry);
846
847 }
848
849}
850
851void VQwDetectorArray::RandomizeEventData(int helicity, double time) {
852
853 for (size_t i = 0; i < fMainDetID.size(); i++) {
854
855 // This is a QwIntegrationPMT
856 if (fMainDetID.at(i).fTypeID == kQwIntegrationPMT)
857 fIntegrationPMT[fMainDetID.at(i).fIndex].RandomizeEventData(helicity, time);
858
859 }
860
861}
862
863void VQwDetectorArray::EncodeEventData(std::vector<UInt_t> &buffer) {
864
865 std::vector<UInt_t> elements;
866 elements.clear();
867
868 // Get all buffers in the order they are defined in the map file
869 for (size_t i = 0; i < fMainDetID.size(); i++) {
870
871 // This is a QwIntegrationPMT
872 if (fMainDetID.at(i).fTypeID == kQwIntegrationPMT)
873 fIntegrationPMT[fMainDetID.at(i).fIndex].EncodeEventData(elements);
874
875 }
876
877
878 // If there is element data, generate the subbank header
879 std::vector<UInt_t> subbankheader;
880 std::vector<UInt_t> rocheader;
881
882 if (elements.size() > 0) {
883
884 // Form CODA subbank header
885 subbankheader.clear();
886 subbankheader.push_back(elements.size() + 1); // subbank size
887 subbankheader.push_back((fCurrentBank_ID << 16) | (0x01 << 8) | (1 & 0xff));
888 // subbank tag | subbank type | event number
889
890 // Form CODA bank/roc header
891 rocheader.clear();
892 rocheader.push_back(subbankheader.size() + elements.size() + 1); // bank/roc size
893 rocheader.push_back((fCurrentROC_ID << 16) | (0x10 << 8) | (1 & 0xff));
894 // bank tag == ROC | bank type | event number
895
896 // Add bank header, subbank header and element data to output buffer
897 buffer.insert(buffer.end(), rocheader.begin(), rocheader.end());
898 buffer.insert(buffer.end(), subbankheader.begin(), subbankheader.end());
899 buffer.insert(buffer.end(), elements.begin(), elements.end());
900
901 }
902
903}
904
905void VQwDetectorArray::RandomizeMollerEvent(int helicity /*, const QwBeamCharge& charge, const QwBeamPosition& xpos, const QwBeamPosition& ypos, const QwBeamAngle& xprime, const QwBeamAngle& yprime, const QwBeamEnergy& energy*/) {
906
907 /* fTargetCharge.PrintInfo();
908 fTargetX.PrintInfo();
909 fTargetY.PrintInfo();
910 fTargetXprime.PrintInfo();
911 fTargetYprime.PrintInfo();
912 fTargetEnergy.PrintInfo();*/
913
914 if(RequestExternalValue("x_targ", &fTargetX)){
915
916 if (bDEBUG){
917
918 dynamic_cast<QwMollerADC_Channel*>(&fTargetX)->PrintInfo();
919 QwWarning << "VQwDetectorArray::RandomizeMollerEvent Found "<<fTargetX.GetElementName()<< QwLog::endl;
920 }
921
922 }else{
923
924 bIsExchangedDataValid = kFALSE;
925 QwError << GetName() << " could not get external value for "
926 << fTargetX.GetElementName() << QwLog::endl;
927
928 }
929
930 if(RequestExternalValue("y_targ", &fTargetY)){
931
932 if (bDEBUG){
933 dynamic_cast<QwMollerADC_Channel*>(&fTargetY)->PrintInfo();
934 QwWarning << "VQwDetectorArray::RandomizeMollerEvent Found "<<fTargetY.GetElementName()<< QwLog::endl;
935 }
936
937 }else{
938
939 bIsExchangedDataValid = kFALSE;
940 QwError << GetName() << " could not get external value for "
941 << fTargetY.GetElementName() << QwLog::endl;
942 }
943
944 if(RequestExternalValue("xp_targ", &fTargetXprime)){
945
946 if (bDEBUG){
947
948 dynamic_cast<QwMollerADC_Channel*>(&fTargetXprime)->PrintInfo();
949 QwWarning << "VQwDetectorArray::RandomizeMollerEvent Found "<<fTargetXprime.GetElementName()<< QwLog::endl;
950 }
951
952 }else{
953
954 bIsExchangedDataValid = kFALSE;
955 QwError << GetName() << " could not get external value for "
956 << fTargetXprime.GetElementName() << QwLog::endl;
957
958 }
959
960 if(RequestExternalValue("yp_targ", &fTargetYprime)){
961
962 if (bDEBUG){
963
964 dynamic_cast<QwMollerADC_Channel*>(&fTargetYprime)->PrintInfo();
965 QwWarning << "VQwDetectorArray::RandomizeMollerEvent Found "<<fTargetYprime.GetElementName()<< QwLog::endl;
966
967 }
968
969 }else{
970
971 bIsExchangedDataValid = kFALSE;
972 QwError << GetName() << " could not get external value for "
973 << fTargetYprime.GetElementName() << QwLog::endl;
974
975 }
976
977 if(RequestExternalValue("e_targ", &fTargetEnergy)){
978
979 if (bDEBUG){
980
981 dynamic_cast<QwMollerADC_Channel*>(&fTargetEnergy)->PrintInfo();
982 QwWarning << "VQwDetectorArray::RandomizeMollerEvent Found "<<fTargetEnergy.GetElementName()<< QwLog::endl;
983
984 }
985
986 }else{
987
988 bIsExchangedDataValid = kFALSE;
989 QwError << GetName() << " could not get external value for "
990 << fTargetEnergy.GetElementName() << QwLog::endl;
991
992 }
993
994 for (size_t i = 0; i < fMainDetID.size(); i++) {
995
997 //fIntegrationPMT[i].PrintInfo();
998
999 }
1000
1001}
1002
1003Int_t VQwDetectorArray::ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words) {
1004
1005 /* Int_t index = GetSubbankIndex(roc_id,bank_id);
1006 if (index>=0 && num_words>0){
1007 // We want to process the configuration data for this ROC.
1008 UInt_t words_read = 0;
1009 for (size_t i = 0; i < fMainDetID.size(); i++) {
1010 words_read += fIntegrationPMT[i].ProcessConfigurationBuffer(&(buffer[words_read]),
1011 num_words-words_read);
1012 }
1013 }*/
1014 return 0;
1015
1016}
1017
1018
1019Int_t VQwDetectorArray::ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words) {
1020
1021 Bool_t lkDEBUG=kFALSE;
1022
1023 Int_t index = GetSubbankIndex(roc_id,bank_id);
1024
1025 if (index>=0 && num_words>0) {
1026
1027 // We want to process this ROC. Begin looping through the data.
1028 if (lkDEBUG)
1029 std::cout << "VQwDetectorArray::ProcessEvBuffer: "
1030 << "Begin processing ROC" << roc_id
1031 << " and subbank "<<bank_id
1032 << " number of words="<<num_words<<std::endl;
1033
1034 for (size_t i=0;i<fMainDetID.size();i++) {
1035
1036 if (fMainDetID[i].fSubbankIndex==index) {
1037
1038 if (fMainDetID[i].fTypeID == kQwIntegrationPMT) {
1039
1040 if (lkDEBUG) {
1041
1042 std::cout<<"found IntegrationPMT data for "<<fMainDetID[i].fdetectorname<<std::endl;
1043 std::cout<<"word left to read in this buffer:"<<num_words-fMainDetID[i].fWordInSubbank<<std::endl;
1044
1045 }
1046
1047 fIntegrationPMT[fMainDetID[i].fIndex].ProcessEvBuffer(&(buffer[fMainDetID[i].fWordInSubbank]),
1048 num_words-fMainDetID[i].fWordInSubbank);
1049
1050 }
1051
1052 }
1053
1054 }
1055
1056 }
1057
1058 return 0;
1059
1060}
1061
1062
1063
1065
1066 Bool_t status=kTRUE;
1067 for(size_t i=0;i<fIntegrationPMT.size();i++){
1068
1069 status &= fIntegrationPMT[i].ApplySingleEventCuts();
1070 if(!status && bDEBUG)
1071 std::cout<<"******* VQwDetectorArray::SingleEventCuts()->IntegrationPMT[ "<<i<<" , "<<fIntegrationPMT[i].GetElementName()<<" ] ******\n";
1072 }
1073
1074 for(size_t i=0;i<fCombinedPMT.size();i++){
1075
1076 status &= fCombinedPMT[i].ApplySingleEventCuts();
1077 if(!status && bDEBUG)
1078 std::cout<<"******* VQwDetectorArray::SingleEventCuts()->CombinedPMT[ "<<i<<" , "<<fCombinedPMT[i].GetElementName()<<" ] ******\n";
1079 }
1080
1081 if (!status)
1082 fMainDetErrorCount++; //failed event counter for VQwDetectorArray
1083
1084 return status;
1085}
1086
1087
1088UInt_t VQwDetectorArray::GetEventcutErrorFlag() { //return the error flag
1089
1090 UInt_t ErrorFlag;
1091
1092 ErrorFlag=0;
1093
1094 for(size_t i=0;i<fIntegrationPMT.size();i++){
1095
1096 ErrorFlag |= fIntegrationPMT[i].GetEventcutErrorFlag();
1097
1098 }
1099
1100 for(size_t i=0;i<fCombinedPMT.size();i++){
1101
1102 ErrorFlag |= fCombinedPMT[i].GetEventcutErrorFlag();
1103
1104 }
1105
1106 return ErrorFlag;
1107
1108}
1109
1111
1112 for(size_t i=0;i<fIntegrationPMT.size();i++){
1113
1114 fIntegrationPMT[i].IncrementErrorCounters();
1115
1116 }
1117
1118 for(size_t i=0;i<fCombinedPMT.size();i++){
1119
1120 fCombinedPMT[i].IncrementErrorCounters();
1121
1122 }
1123
1124}
1125
1126//inherited from the VQwSubsystemParity; this will display the error summary
1128
1129 QwMessage<<"*********VQwDetectorArray Error Summary****************"<<QwLog::endl;
1131
1132 for(size_t i=0;i<fIntegrationPMT.size();i++){
1133
1134 //std::cout<<" IntegrationPMT ["<<i<<"] "<<std::endl;
1135 fIntegrationPMT[i].PrintErrorCounters();
1136
1137 }
1138
1139 for(size_t i=0;i<fCombinedPMT.size();i++){
1140
1141 //std::cout<<" CombinedPMT ["<<i<<"] "<<std::endl;
1142 fCombinedPMT[i].PrintErrorCounters();
1143
1144 }
1145
1147
1148}
1149
1151
1152 Bool_t burpstatus = kFALSE;
1153 VQwSubsystem* tmp = const_cast<VQwSubsystem *>(subsys);
1154
1155 if(Compare(tmp)) {
1156
1157 const VQwDetectorArray* input = dynamic_cast<const VQwDetectorArray*>(subsys);
1158
1159 for(size_t i=0;i<input->fIntegrationPMT.size();i++)
1160 burpstatus |= (this->fIntegrationPMT[i]).CheckForBurpFail(&(input->fIntegrationPMT[i]));
1161
1162 for(size_t i=0;i<input->fCombinedPMT.size();i++)
1163 burpstatus |= (this->fCombinedPMT[i]).CheckForBurpFail(&(input->fCombinedPMT[i]));
1164
1165 }
1166
1167 return burpstatus;
1168
1169}
1170
1171
1173
1174 VQwSubsystem* tmp = const_cast<VQwSubsystem*>(ev_error);
1175
1176 if(Compare(tmp)){
1177
1178 const VQwDetectorArray* input = dynamic_cast<const VQwDetectorArray*> (ev_error);
1179
1180 for (size_t i=0;i<input->fIntegrationPMT.size();i++)
1181 this->fIntegrationPMT[i].UpdateErrorFlag(&(input->fIntegrationPMT[i]));
1182
1183 for (size_t i=0;i<input->fCombinedPMT.size();i++)
1184 this->fCombinedPMT[i].UpdateErrorFlag(&(input->fCombinedPMT[i]));
1185
1186 }
1187
1188}
1189
1190
1192
1193 for (size_t i=0;i<fIntegrationPMT.size();i++)
1195
1196 for (size_t i=0;i<fCombinedPMT.size();i++) {
1197
1198 //std::cout<<"Process combination "<<i<<std::endl;
1199 fCombinedPMT[i].ProcessEvent();
1200
1201 }
1202
1203 return;
1204
1205}
1206
1207/**
1208 * Exchange data between subsystems
1209 */
1210
1212
1213 //QwWarning << "VQwDetectorArray::ExchangeProcessedData "<< QwLog::endl;
1214 bIsExchangedDataValid = kTRUE;
1215
1216 if (1==1 || bNormalization) {
1217
1218 if(RequestExternalValue("q_targ", &fTargetCharge)) {
1219
1220 if (bDEBUG) {
1221
1222 QwWarning << "VQwDetectorArray::ExchangeProcessedData Found "<<fTargetCharge.GetElementName()<< QwLog::endl;
1223 //QwWarning <<"****VQwDetectorArray****"<< QwLog::endl;
1224 (dynamic_cast<QwMollerADC_Channel*>(&fTargetCharge))->PrintInfo();
1225
1226 }
1227
1228 } else {
1229
1230 bIsExchangedDataValid = kFALSE;
1231 QwError << GetName() << " could not get external value for "
1232 << fTargetCharge.GetElementName() << QwLog::endl;
1233
1234 }
1235
1236 }
1237
1238}
1239
1240
1242
1244
1245 //data is valid, process it
1246
1247 if (bDEBUG) {
1248
1249 Double_t pedestal = fTargetCharge.GetPedestal();
1250 Double_t calfactor = fTargetCharge.GetCalibrationFactor();
1251 Double_t volts = fTargetCharge.GetAverageVolts();
1252
1253 std::cout<<"VQwDetectorArray::ProcessEvent_2(): processing with exchanged data"<<std::endl;
1254 std::cout<<"pedestal, calfactor, average volts = "<<pedestal<<", "<<calfactor<<", "<<volts<<std::endl;
1255
1256 }
1257
1259 this->DoNormalization();
1260
1261 } else {
1262
1263 QwWarning<<"VQwDetectorArray::ProcessEvent_2(): could not get all external values."<<QwLog::endl;
1264
1265 }
1266
1267}
1268
1269
1270
1271
1272void VQwDetectorArray::ConstructHistograms(TDirectory *folder, TString &prefix) {
1273
1274 for (size_t i=0;i<fIntegrationPMT.size();i++)
1275 fIntegrationPMT[i].ConstructHistograms(folder,prefix);
1276
1277 for (size_t i=0;i<fCombinedPMT.size();i++)
1278 fCombinedPMT[i].ConstructHistograms(folder,prefix);
1279
1280 return;
1281
1282}
1283
1284
1286
1287 for (size_t i=0;i<fIntegrationPMT.size();i++)
1289
1290 for (size_t i=0;i<fCombinedPMT.size();i++)
1292
1293 return;
1294
1295}
1296
1297
1298void VQwDetectorArray::ConstructBranchAndVector(TTree *tree, TString & prefix, QwRootTreeBranchVector &values) {
1299
1300 for (size_t i=0;i<fIntegrationPMT.size();i++)
1301 fIntegrationPMT[i].ConstructBranchAndVector(tree, prefix, values);
1302
1303 for (size_t i=0;i<fCombinedPMT.size();i++)
1304 fCombinedPMT[i].ConstructBranchAndVector(tree, prefix, values);
1305
1306 return;
1307
1308}
1309
1310void VQwDetectorArray::ConstructBranch(TTree *tree, TString & prefix) {
1311
1312 for (size_t i=0;i<fIntegrationPMT.size();i++)
1313 fIntegrationPMT[i].ConstructBranch(tree, prefix);
1314
1315 for (size_t i=0;i<fCombinedPMT.size();i++)
1316 fCombinedPMT[i].ConstructBranch(tree, prefix);
1317
1318 return;
1319
1320}
1321
1322void VQwDetectorArray::ConstructBranch(TTree *tree, TString & prefix, QwParameterFile& trim_file) {
1323
1324 TString tmp;
1325 std::unique_ptr<QwParameterFile> nextmodule;
1326 trim_file.RewindToFileStart();
1327 tmp="QwIntegrationPMT";
1328 trim_file.RewindToFileStart();
1329
1330 if (trim_file.FileHasModuleHeader(tmp)){
1331
1332 nextmodule=trim_file.ReadUntilNextModule();//This section contains sub modules and or channels to be included in the tree
1333
1334 for (size_t i=0;i<fIntegrationPMT.size();i++)
1335 fIntegrationPMT[i].ConstructBranch(tree, prefix, *nextmodule);
1336
1337 }
1338
1339 tmp="QwCombinedPMT";
1340 trim_file.RewindToFileStart();
1341
1342 if (trim_file.FileHasModuleHeader(tmp)){
1343
1344 nextmodule=trim_file.ReadUntilNextModule();//This section contains sub modules and or channels to be included in the tree
1345
1346 for (size_t i=0;i<fCombinedPMT.size();i++)
1347 fCombinedPMT[i].ConstructBranch(tree, prefix, *nextmodule );
1348 }
1349
1350 return;
1351}
1352
1354
1355 for (size_t i=0;i<fIntegrationPMT.size();i++)
1357
1358 for (size_t i=0;i<fCombinedPMT.size();i++)
1359 fCombinedPMT[i].FillTreeVector(values);
1360
1361 return;
1362
1363}
1364
1365#ifdef HAS_RNTUPLE_SUPPORT
1366void VQwDetectorArray::ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs) {
1367
1368 for (size_t i=0;i<fIntegrationPMT.size();i++)
1369 fIntegrationPMT[i].ConstructNTupleAndVector(model, prefix, values, fieldPtrs);
1370
1371 for (size_t i=0;i<fCombinedPMT.size();i++)
1372 fCombinedPMT[i].ConstructNTupleAndVector(model, prefix, values, fieldPtrs);
1373
1374 return;
1375
1376}
1377
1378void VQwDetectorArray::FillNTupleVector(std::vector<Double_t>& values) const {
1379
1380 for (size_t i=0;i<fIntegrationPMT.size();i++)
1381 fIntegrationPMT[i].FillNTupleVector(values);
1382
1383 for (size_t i=0;i<fCombinedPMT.size();i++)
1384 fCombinedPMT[i].FillNTupleVector(values);
1385
1386 return;
1387
1388}
1389#endif
1390
1391
1392const QwIntegrationPMT* VQwDetectorArray::GetChannel(const TString name) const {
1393
1394 return GetIntegrationPMT(name);
1395
1396}
1397
1398
1400
1401 // std::cout<<" Here in VQwDetectorArray::Compare \n";
1402
1403 Bool_t res=kTRUE;
1404
1405 if (typeid(*value)!=typeid(*this)) {
1406
1407 res=kFALSE;
1408 //std::cout<<" types are not ok \n";
1409 //std::cout<<" this is bypassed just for now but should be fixed eventually \n";
1410
1411 } else {
1412
1413 VQwDetectorArray* input = dynamic_cast<VQwDetectorArray*>(value);
1414
1415 if (input->fIntegrationPMT.size()!=fIntegrationPMT.size()
1416 || input->fCombinedPMT.size()!=fCombinedPMT.size() ) {
1417
1418 res=kFALSE;
1419 //std::cout<<" not the same number of channels \n";
1420
1421 }
1422
1423 }
1424
1425 return res;
1426
1427}
1428
1429
1431
1432 // std::cout<<" here in VQwDetectorArray::operator= \n";
1433
1434 if (this != value && Compare(value)) {
1435
1436 //VQwSubsystem::operator=(value);
1437 VQwDetectorArray* input = dynamic_cast<VQwDetectorArray*> (value);
1438
1439 for (size_t i=0;i<input->fIntegrationPMT.size();i++)
1440 this->fIntegrationPMT[i]=input->fIntegrationPMT[i];
1441
1442 for (size_t i=0;i<input->fCombinedPMT.size();i++)
1443 (this->fCombinedPMT[i])=(input->fCombinedPMT[i]);
1444
1445 }
1446
1447 return *this;
1448
1449}
1450
1451
1453
1454 if (Compare(value)) {
1455
1456 VQwDetectorArray* input= dynamic_cast<VQwDetectorArray*>(value) ;
1457
1458 for (size_t i=0;i<input->fIntegrationPMT.size();i++)
1459 this->fIntegrationPMT[i]+=input->fIntegrationPMT[i];
1460
1461 for (size_t i=0;i<input->fCombinedPMT.size();i++)
1462 this->fCombinedPMT[i]+=input->fCombinedPMT[i];
1463
1464 }
1465
1466 return *this;
1467
1468}
1469
1470
1472
1473 if (Compare(value)) {
1474
1475 VQwDetectorArray* input= dynamic_cast<VQwDetectorArray*>(value);
1476
1477 for (size_t i=0;i<input->fIntegrationPMT.size();i++)
1478 this->fIntegrationPMT[i]-=input->fIntegrationPMT[i];
1479
1480 for (size_t i=0;i<input->fCombinedPMT.size();i++)
1481 this->fCombinedPMT[i]-=input->fCombinedPMT[i];
1482
1483 }
1484
1485 return *this;
1486
1487}
1488
1489
1491
1492 if (Compare(numer)&&Compare(denom)) {
1493
1494 VQwDetectorArray* innumer= dynamic_cast<VQwDetectorArray*>(numer) ;
1495 VQwDetectorArray* indenom= dynamic_cast<VQwDetectorArray*>(denom) ;
1496
1497 for (size_t i=0;i<innumer->fIntegrationPMT.size();i++)
1498 this->fIntegrationPMT[i].Ratio(innumer->fIntegrationPMT[i],indenom->fIntegrationPMT[i]);
1499
1500 for (size_t i=0;i<innumer->fCombinedPMT.size();i++)
1501 this->fCombinedPMT[i].Ratio(innumer->fCombinedPMT[i],indenom->fCombinedPMT[i]);
1502
1503 }
1504
1505 return;
1506
1507}
1508
1509
1510void VQwDetectorArray::Scale(Double_t factor) {
1511
1512 for (size_t i=0;i<fIntegrationPMT.size();i++)
1513 fIntegrationPMT[i].Scale(factor);
1514
1515 for (size_t i=0;i<fCombinedPMT.size();i++)
1516 fCombinedPMT[i].Scale(factor);
1517
1518 return;
1519
1520}
1521
1522//*****************************************************************//
1523
1525
1526 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
1527 fIntegrationPMT[i].Normalize(denom);
1528
1529 for (size_t i = 0; i < fCombinedPMT.size(); i++)
1530 fCombinedPMT[i].Normalize(denom);
1531
1532}
1533
1534
1535
1537
1538 for (size_t i=0;i<fIntegrationPMT.size();i++)
1540
1541 for (size_t i=0;i<fCombinedPMT.size();i++)
1543
1544 return;
1545
1546}
1547
1548
1549void VQwDetectorArray::AccumulateRunningSum(VQwSubsystem* value1, Int_t count, Int_t ErrorMask) {
1550
1551 if (Compare(value1)) {
1552
1553 VQwDetectorArray* value = dynamic_cast<VQwDetectorArray*>(value1);
1554
1555 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
1556 fIntegrationPMT[i].AccumulateRunningSum(value->fIntegrationPMT[i], count, ErrorMask);
1557
1558 for (size_t i = 0; i < fCombinedPMT.size(); i++)
1559 fCombinedPMT[i].AccumulateRunningSum(value->fCombinedPMT[i], count, ErrorMask);
1560
1561 }
1562
1563}
1564
1566
1567 if (Compare(value1)) {
1568
1569 VQwDetectorArray* value = dynamic_cast<VQwDetectorArray*>(value1);
1570
1571 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
1573
1574 for (size_t i = 0; i < fCombinedPMT.size(); i++)
1575 fCombinedPMT[i].DeaccumulateRunningSum(value->fCombinedPMT[i], ErrorMask);
1576
1577 }
1578
1579}
1580
1581
1587
1588//*****************************************************************//
1589
1591
1592 Bool_t ldebug=kFALSE;
1593
1594 if (ldebug) {
1595
1596 std::cout<<"VQwDetectorArray::GetDetectorIndex\n";
1597 std::cout<<"type_id=="<<type_id<<" name="<<name<<"\n";
1598 std::cout<<fMainDetID.size()<<" already registered detector\n";
1599 }
1600
1601 Int_t result=-1;
1602 for (size_t i=0;i<fMainDetID.size();i++) {
1603
1604 if (fMainDetID[i].fTypeID==type_id)
1605 if (fMainDetID[i].fdetectorname==name) {
1606
1607 result=fMainDetID[i].fIndex;
1608
1609 if (ldebug)
1610 std::cout<<"testing against ("<<fMainDetID[i].fTypeID
1611 <<","<<fMainDetID[i].fdetectorname<<")=>"<<result<<"\n";
1612
1613 }
1614 }
1615
1616 return result;
1617
1618}
1619
1621
1622 TString tmpname = name;
1623 tmpname.ToLower();
1624 if (! fIntegrationPMT.empty()) {
1625
1626 for (size_t i=0;i<fIntegrationPMT.size();i++) {
1627
1628 if (fIntegrationPMT.at(i).GetElementName() == tmpname) {
1629
1630 //std::cout<<"Get IntegrationPMT "<<tmpname<<std::endl;
1631 return &(fIntegrationPMT.at(i));
1632
1633 }
1634
1635 }
1636
1637 }
1638
1639 QwMessage << "VQwDetectorArray::GetIntegrationPMT: cannot find channel " << tmpname << QwLog::endl;
1640
1641 return NULL;
1642
1643}
1644
1645
1646const QwCombinedPMT* VQwDetectorArray::GetCombinedPMT(const TString name) const {
1647 TString tmpname = name;
1648 tmpname.ToLower();
1649 if (! fCombinedPMT.empty())
1650 {
1651 for (size_t i=0;i<fCombinedPMT.size();i++)
1652 {
1653 if (fCombinedPMT.at(i).GetElementName() == tmpname)
1654 {
1655 //std::cout<<"Get CombinedPMT "<<tmpname<<std::endl;
1656 return &(fCombinedPMT.at(i));
1657 }
1658 }
1659 }
1660 QwMessage << "VQwDetectorArray::GetCombinedPMT: cannot find channel " << tmpname << QwLog::endl;
1661 return NULL;
1662}
1663
1665
1667
1668 try {
1669
1670 this->Normalize(&fTargetCharge);
1671
1672 }
1673
1674 catch (std::exception& e) {
1675
1676 std::cerr << e.what() << std::endl;
1677
1678 }
1679
1680 }
1681
1682}
1683
1684#ifdef __USE_DATABASE__
1685void VQwDetectorArray::FillDB(QwParityDB *db, TString datatype) {
1686
1687 Bool_t local_print_flag = false;
1688
1689 if(local_print_flag) {
1690
1691 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1692 QwMessage << " VQwDetectorArray::FillDB " << QwLog::endl;
1693 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1694
1695 }
1696
1697 std::vector<QwDBInterface> interface;
1698 std::vector<QwParitySchema::md_data_row> entrylist;
1699
1700 UInt_t analysis_id = db->GetAnalysisID();
1701
1702 TString measurement_type;
1703 measurement_type = QwDBInterface::DetermineMeasurementTypeID(datatype);
1704
1705 UInt_t i,j;
1706 i = j = 0;
1707
1708 if(local_print_flag) QwMessage << QwColor(Qw::kGreen) << "IntegrationPMT" <<QwLog::endl;
1709
1710 for(i=0; i<fIntegrationPMT.size(); i++) {
1711
1712 interface.clear();
1713 interface = fIntegrationPMT[i].GetDBEntry();
1714
1715 for(j=0; j<interface.size(); j++) {
1716
1717 interface.at(j).SetAnalysisID( analysis_id );
1718 interface.at(j).SetMainDetectorID( db );
1719 interface.at(j).SetMeasurementTypeID( measurement_type );
1720 interface.at(j).PrintStatus( local_print_flag );
1721 interface.at(j).AddThisEntryToList( entrylist );
1722 }
1723 }
1724
1725 if(local_print_flag) QwMessage << QwColor(Qw::kGreen) << "Combined PMT" <<QwLog::endl;
1726
1727 for(i=0; i< fCombinedPMT.size(); i++) {
1728
1729 interface.clear();
1730 interface = fCombinedPMT[i].GetDBEntry();
1731
1732 for(j=0; j<interface.size(); j++) {
1733
1734 interface.at(j).SetAnalysisID( analysis_id );
1735 interface.at(j).SetMainDetectorID( db );
1736 interface.at(j).SetMeasurementTypeID( measurement_type );
1737 interface.at(j).PrintStatus( local_print_flag );
1738 interface.at(j).AddThisEntryToList( entrylist );
1739 }
1740
1741 }
1742
1743 if(local_print_flag) {
1744
1745 QwMessage << QwColor(Qw::kGreen) << "Entrylist Size : "
1746 << QwColor(Qw::kBoldRed) << entrylist.size()
1747 << QwColor(Qw::kNormal) << QwLog::endl;
1748
1749 }
1750
1751 // Check the entrylist size, if it isn't zero, start to query..
1752 if( entrylist.size() ) {
1753 auto c = db->GetScopedConnection();
1754 for (const auto& entry : entrylist) {
1755 c->QueryExecute(entry.insert_into());
1756 }
1757 } else {
1758 QwMessage << "VQwDetectorArray::FillDB :: This is the case when the entrylist contains nothing in "<< datatype.Data() << QwLog::endl;
1759 }
1760}
1761#endif
1762
1764
1765 QwMessage << "=== VQwDetectorArray: " << GetName() << " ===" << QwLog::endl;
1766
1767 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
1769
1770 for (size_t i = 0; i < fCombinedPMT.size(); i++)
1772
1773}
1774
1776
1777 std::cout<<"Name of the subsystem ="<<fSystemName<<"\n";
1778
1779 std::cout<<"there are "<<fIntegrationPMT.size()<<" IntegrationPMT \n";
1780 std::cout<<" "<<fCombinedPMT.size()<<" CombinedPMT \n";
1781
1782 std::cout<<" Printing Running AVG and other channel info"<<std::endl;
1783
1784 for (size_t i = 0; i < fIntegrationPMT.size(); i++)
1786 for (size_t i = 0; i < fCombinedPMT.size(); i++)
1788
1789}
1790
1792
1793 for (size_t i=0;i<fMainDetID.size();i++) {
1794
1795 std::cout<<"============================="<<std::endl;
1796 std::cout<<" Detector ID="<<i<<std::endl;
1797 fMainDetID[i].Print();
1798
1799 }
1800
1801 return;
1802
1803}
1804
1805#ifdef __USE_DATABASE__
1806void VQwDetectorArray::FillErrDB(QwParityDB *db, TString datatype) {
1807
1808 Bool_t local_print_flag = false;
1809 if(local_print_flag){
1810
1811 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1812 QwMessage << " QwDetectorArrayID::FillErrDB " << QwLog::endl;
1813 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1814
1815 }
1816
1817
1818 std::vector<QwErrDBInterface> interface;
1819 std::vector<QwParitySchema::md_errors_row> entrylist;
1820
1821 UInt_t analysis_id = db->GetAnalysisID();
1822
1823 UInt_t i,j;
1824 i = j = 0;
1825 if(local_print_flag) QwMessage << QwColor(Qw::kGreen) << "IntegrationPMT" <<QwLog::endl;
1826
1827 for(i=0; i<fIntegrationPMT.size(); i++) {
1828
1829 interface.clear();
1830 interface = fIntegrationPMT[i].GetErrDBEntry();
1831
1832 for(j=0; j<interface.size(); j++) {
1833
1834 interface.at(j).SetAnalysisID ( analysis_id );
1835 interface.at(j).SetMainDetectorID ( db );
1836 interface.at(j).PrintStatus ( local_print_flag );
1837 interface.at(j).AddThisEntryToList( entrylist );
1838
1839 }
1840
1841 }
1842
1843 if(local_print_flag) QwMessage << QwColor(Qw::kGreen) << "Combined PMT" <<QwLog::endl;
1844
1845 for(i=0; i< fCombinedPMT.size(); i++) {
1846
1847 interface.clear();
1848 interface = fCombinedPMT[i].GetErrDBEntry();
1849
1850 for(j=0; j<interface.size(); j++) {
1851
1852 interface.at(j).SetAnalysisID ( analysis_id );
1853 interface.at(j).SetMainDetectorID ( db );
1854 interface.at(j).PrintStatus ( local_print_flag );
1855 interface.at(j).AddThisEntryToList( entrylist );
1856
1857 }
1858
1859 }
1860
1861 if(local_print_flag) {
1862
1863 QwMessage << QwColor(Qw::kGreen) << "Entrylist Size : "
1864 << QwColor(Qw::kBoldRed) << entrylist.size()
1865 << QwColor(Qw::kNormal) << QwLog::endl;
1866
1867 }
1868
1869 // Check the entrylist size, if it isn't zero, start to query..
1870 if( entrylist.size() ) {
1871 auto c = db->GetScopedConnection();
1872 for (const auto& entry : entrylist) {
1873 c->QueryExecute(entry.insert_into());
1874 }
1875 } else {
1876 QwMessage << "VQwDetectorArray::FillErrDB :: This is the case when the entrylist contains nothing in "<< datatype.Data() << QwLog::endl;
1877 }
1878}
1879#endif
1880
1881
1883
1884 Bool_t local_print_flag = false;
1885 Bool_t local_add_element= type.Contains("yield");
1886
1887 if(local_print_flag){
1888
1889 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1890 QwMessage << " QwDetectorArrayID::WritePromptSummary() " << QwLog::endl;
1891 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
1892 }
1893
1894 const VQwHardwareChannel* tmp_channel = 0;
1895 TString element_name = "";
1896 Double_t element_value = 0.0;
1897 Double_t element_value_err = 0.0;
1898 Double_t element_value_width = 0.0;
1899
1900 PromptSummaryElement *local_ps_element = NULL;
1901 Bool_t local_add_these_elements= false;
1902
1903 for (size_t i = 0; i < fMainDetID.size(); i++) {
1904
1905 element_name = fMainDetID[i].fdetectorname;
1906 tmp_channel=GetIntegrationPMT(element_name)->GetChannel(element_name);
1907 element_value = 0.0;
1908 element_value_err = 0.0;
1909 element_value_width = 0.0;
1910
1911
1912 local_add_these_elements=element_name.Contains("sam"); // Need to change this to add other detectorss in summary
1913
1914 if(local_add_these_elements&&local_add_element){
1915
1916 ps->AddElement(new PromptSummaryElement(element_name));
1917
1918 }
1919
1920
1921 local_ps_element=ps->GetElementByName(element_name);
1922
1923
1924 if(local_ps_element) {
1925
1926 element_value = tmp_channel->GetValue();
1927 element_value_err = tmp_channel->GetValueError();
1928 element_value_width = tmp_channel->GetValueWidth();
1929
1930 local_ps_element->Set(type, element_value, element_value_err, element_value_width);
1931
1932 }
1933
1934 if( local_print_flag && local_ps_element) {
1935
1936 printf("Type %12s, Element %32s, value %12.4e error %8.4e width %12.4e\n",
1937 type.Data(), element_name.Data(), element_value, element_value_err, element_value_width);
1938
1939 }
1940
1941 }
1942
1943 return;
1944
1945}
1946
1947
1948
1950
1951 std::cout<<std::endl<<"Detector name= "<<fdetectorname<<std::endl;
1952 std::cout<<"SubbankkIndex= "<<fSubbankIndex<<std::endl;
1953 std::cout<<"word index in subbank= "<<fWordInSubbank<<std::endl;
1954 std::cout<<"module type= "<<fmoduletype<<std::endl;
1955 std::cout<<"detector type= "<<fdetectortype<<" index= "<<fTypeID<<std::endl;
1956 std::cout<<"Index of this detector in the vector of similar detector= "<<fIndex<<std::endl;
1957 std::cout<<"Subelement index= "<<fSubelement<<std::endl;
1958 std::cout<<"==========================================\n";
1959
1960 return;
1961
1962}
1963
1964
1965
1966
A logfile class, based on an identical class in the Hermes analyzer.
#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
Array container for managing multiple subsystems.
Prompt summary data management.
EQwPMTInstrumentType GetQwPMTInstrumentType(TString name)
Definition QwTypes.cc:16
TString GetQwPMTInstrumentTypeName(EQwPMTInstrumentType type)
Definition QwTypes.cc:81
ULong64_t BankID_t
Definition QwTypes.h:21
UInt_t GetGlobalErrorFlag(TString evtype, Int_t evMode, Double_t stabilitycut)
Definition QwTypes.cc:132
EQwPMTInstrumentType
Definition QwTypes.h:135
@ kQwUnknownPMT
Definition QwTypes.h:136
@ kQwCombinedPMT
Definition QwTypes.h:139
@ kQwIntegrationPMT
Definition QwTypes.h:137
UInt_t ROCID_t
Definition QwTypes.h:20
Virtual base class for detector arrays (PMTs, etc.)
@ kBoldRed
Definition QwColor.h:79
@ kGreen
Definition QwColor.h:77
@ kNormal
Definition QwColor.h:81
Bool_t RequestExternalValue(const TString &name, VQwHardwareChannel *value) const
Bool_t PublishInternalValue(const TString name, const TString desc, const VQwHardwareChannel *element) const
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
Concrete hardware channel for Moller ADC modules (6x32-bit words)
static void PrintErrorCounterHead()
static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex)
static void PrintErrorCounterTail()
Command-line and configuration file options processor.
Definition QwOptions.h:141
T GetValue(const std::string &key)
Get a templated value.
Definition QwOptions.h:236
po::options_description_easy_init AddOptions(const std::string &blockname="Specialized options")
Add an option to a named block or create new block.
Definition QwOptions.h:170
Configuration file parser with flexible tokenization and search capabilities.
T GetTypedNextToken()
Get next token into specific type.
std::unique_ptr< QwParameterFile > ReadUntilNextModule(const bool add_current_line=false)
Bool_t PopValue(const std::string keyname, T &retvalue)
void TrimWhitespace(TString::EStripType head_tail=TString::kBoth)
void SetCommentChars(const std::string value)
Set various sets of special characters.
void TrimComment(const char commentchar)
const std::pair< TString, TString > GetParamFileNameContents()
Bool_t FileHasModuleHeader(const std::string &secname)
std::unique_ptr< QwParameterFile > ReadNextSection(std::string &secname, const bool keep_header=false)
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
The pure virtual base class of all data elements.
Abstract base for concrete hardware channels implementing dual-operator pattern.
Double_t GetValueError() const
Double_t GetValueWidth() const
Double_t GetValue() const
std::vector< std::vector< TString > > fPublishList
BankID_t fCurrentBank_ID
Bank ID (and Marker word) that is currently being processed;.
TString fSystemName
Name of this subsystem.
Int_t GetSubbankIndex() const
void RegisterRocBankMarker(QwParameterFile &mapstr)
static void DefineOptions()
Define options function (note: no virtual static functions in C++)
TString GetName() const
std::map< TString, TString > fDetectorMaps
Map of file name to full path or content.
VQwSubsystem(const TString &name)
Constructor with name.
ROCID_t fCurrentROC_ID
ROC ID that is currently being processed.
Combines multiple integration PMTs into weighted sum/average.
void SetNormalizability(Bool_t isnormalizable)
void SetBlindability(Bool_t isblindable)
const QwMollerADC_Channel * GetChannel(const TString name) const
Integration PMT providing yield/diff/asym readout from Moller ADC.
const QwMollerADC_Channel * GetChannel(const TString name) const
void SetBlindability(Bool_t isblindable)
void SetNormalizability(Bool_t isnormalizable)
Identifier and mapping information for detector-array channels.
std::vector< Double_t > fWeight
std::vector< TString > fCombinedChannelNames
EQwPMTInstrumentType fTypeID
Int_t LoadInputParameters(TString pedestalfile) override
Mandatory parameter file definition.
Bool_t ApplySingleEventCuts() override
Apply the single event cuts.
void LoadEventCuts_Fin(Int_t &eventcut_flag) override
const QwCombinedPMT * GetCombinedPMT(const TString name) const
void Ratio(VQwSubsystem *numer, VQwSubsystem *denom) override
Int_t GetDetectorIndex(EQwPMTInstrumentType TypeID, TString name)
void ProcessOptions(QwOptions &options) override
void AccumulateRunningSum(VQwSubsystem *value, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
Update the running sums for devices.
QwBeamPosition fTargetX
const QwIntegrationPMT * GetChannel(const TString name) const
EQwPMTInstrumentType GetDetectorTypeID(TString name)
void DoNormalization(Double_t factor=1.0)
void LoadMockDataParameters(TString pedestalfile) override
VQwDetectorArray()
Private default constructor (not implemented, will throw linker error on use)
QwBeamAngle fTargetXprime
QwBeamPosition fTargetY
void FillTreeVector(QwRootTreeBranchVector &values) const override
Fill the tree vector.
Int_t LoadChannelMap(TString mapfile) override
void IncrementErrorCounters() override
Increment the error counters.
void DeaccumulateRunningSum(VQwSubsystem *value, Int_t ErrorMask=0xFFFFFFF) override
remove one entry from the running sums for devices
void ConstructBranch(TTree *tree, TString &prefix) override
Construct the branch and tree vector.
void WritePromptSummary(QwPromptSummary *ps, TString type) override
Int_t ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
TODO: The non-event-type-aware ProcessEvBuffer routine should be replaced with the event-type-aware v...
std::vector< QwIntegrationPMT > fIntegrationPMT
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override
Report the number of events failed due to HW and event cut failures.
void LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag) override
Bool_t Compare(VQwSubsystem *source)
VQwSubsystem & operator-=(VQwSubsystem *value) override
void PrintDetectorID() const
void PrintErrorCounters() const override
QwBeamAngle fTargetYprime
Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
static const Bool_t bDEBUG
void ExchangeProcessedData() override
Bool_t PublishInternalValues() const override
std::vector< QwDetectorArrayID > fMainDetID
void EncodeEventData(std::vector< UInt_t > &buffer) override
void Normalize(VQwDataElement *denom)
VQwSubsystem & operator+=(VQwSubsystem *value) override
QwBeamCharge fTargetCharge
std::vector< QwCombinedPMT > fCombinedPMT
VQwSubsystem & operator=(VQwSubsystem *value) override
Assignment Note: Must be called at the beginning of all subsystems routine call to operator=(VQwSubsy...
void PrintValue() const override
Print values of all channels.
void SetRandomEventAsymmetry(Double_t asymmetry)
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
Construct the branch and tree vector.
UInt_t GetEventcutErrorFlag() override
Return the error flag to the top level routines related to stability checks and ErrorFlag updates.
Bool_t PublishByRequest(TString device_name) override
void SetRandomEventParameters(Double_t mean, Double_t sigma)
void PrintInfo() const override
Print some information about the subsystem.
void CalculateRunningAverage() override
Calculate the average for all good events.
void Scale(Double_t factor) override
void ProcessEvent() override
void RandomizeMollerEvent(int helicity)
void ClearEventData() override
void ProcessEvent_2() override
Process the event data again, including data from other subsystems. Not all derived classes will requ...
const QwIntegrationPMT * GetIntegrationPMT(const TString name) const
void RandomizeEventData(int helicity=0, Double_t time=0.0) override
QwBeamEnergy fTargetEnergy
void FillHistograms() override
Fill the histograms for this subsystem.
virtual void FillDB(QwParityDB *, TString)
Fill the database.
virtual UInt_t UpdateErrorFlag()
Uses the error flags of contained data elements to update Returns the error flag to the top level rou...
virtual void FillErrDB(QwParityDB *, TString)