JAPAn
Just Another Parity Analyzer
Loading...
Searching...
No Matches
QwBeamMod.cc
Go to the documentation of this file.
1/*!
2 * \file QwBeamMod.cc
3 * \brief Beam modulation subsystem implementation
4 * \author Joshua Hoskins, Ezekiel Wertz, Paul King
5 * \date 2010-05-25
6 */
7
8#include "QwBeamMod.h"
9
10// System headers
11#include <iostream>
12#include <stdexcept>
13
14#ifdef HAS_RNTUPLE_SUPPORT
15// ROOT headers for RNTuple support
16#include <ROOT/RNTupleModel.hxx>
17#include <ROOT/RNTupleWriter.hxx>
18#endif // HAS_RNTUPLE_SUPPORT
19
20// Qweak headers
21#include "QwLog.h"
22#include "QwParameterFile.h"
23#include "QwHistogramHelper.h"
24#include "QwRootFile.h"
25#ifdef __USE_DATABASE__
26#include "QwParitySchema.h"
27#include "QwParityDB.h"
28#endif // __USE_DATABASE__
29
30#include "QwScaler_Channel.h"
31
32// Root plotting headers
33#include "TCanvas.h"
34#include "TF1.h"
35#include "TMath.h"
36
37
38//*****************************************************************
40 //Handle command line options
41}
42
43Int_t QwBeamMod::LoadChannelMap(TString mapfile){
44 // std::cout <<"Here in LoadChannelMap" << std::endl;
45 Bool_t ldebug=kFALSE;
46
47 Int_t wordsofar=0;
48 Int_t bankindex=-1;
49
50
51 TString varname, varvalue;
52
53 // Open the file
54 QwParameterFile mapstr(mapfile.Data());
56 mapstr.EnableGreediness();
57 mapstr.SetCommentChars("!");
58 while (mapstr.ReadNextLine()) {
60
63 if ((bankindex+1)>0){
64 UInt_t numbanks = UInt_t(bankindex+1);
65 if (fWordsPerSubbank.size()<numbanks){
66 fWordsPerSubbank.resize(numbanks,
67 std::pair<Int_t, Int_t>(fWord.size(),fWord.size()));
68 }
69 }
70 wordsofar=0;
71 }
72 mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
73 if (mapstr.LineIsEmpty()) continue;
74
75
76 if (mapstr.HasVariablePair("=",varname,varvalue))
77 {
78 // This is a declaration line. Decode it.
79 varname.ToLower();
80 QwWarning << "QwBeamMod::LoadChannelMap: Unrecognized declaration "
81 << varname << " = " << varvalue << QwLog::endl;
82
83 } else {
84 Bool_t lineok=kTRUE;
85 // QwModChannelID localModChannelID(bankindex, wordsofar,namech, modtype, this);
86 QwModChannelID localModChannelID(bankindex, mapstr);
87 TString namech = localModChannelID.fmodulename;
88 if(localModChannelID.fmoduletype=="VQWK") {
89 wordsofar+=6;
90
91 if (lineok){
92 VQwHardwareChannel* localchan1 = new QwVQWK_Channel();
93
94 localchan1->InitializeChannel(GetName(),"QwBeamMod",localModChannelID.fmodulename,"raw");
95 fModChannel.push_back(localchan1);
96 fModChannel[fModChannel.size()-1]->LoadChannelParameters(mapstr);
97 localModChannelID.fIndex=fModChannel.size()-1;
98 fModChannelID.push_back(localModChannelID);
99 /* localchan1->PrintInfo();
100 localchan1->PrintValue();
101 localchan1->PrintErrorCounters();*/
102 // Store reference to ramp channel
103 if (localModChannelID.fmodulename == "ramp") {
104 fRampChannelIndex = fModChannel.size() - 1;
105 }
106
107 }
108
109 if(ldebug)
110 {
111 localModChannelID.Print();
112 std::cout<<"line ok=";
113 if(lineok) std::cout<<"TRUE"<<std::endl;
114 else std::cout<<"FALSE"<<std::endl;
115 }
116 }
117
118 if(localModChannelID.fmoduletype=="SCALER" || localModChannelID.fmoduletype == "SIS3801D24") {
119 wordsofar+=1;
120
121 if (lineok){
122 VQwHardwareChannel* localchan2 = new QwSIS3801D24_Channel();
123
124 localchan2->InitializeChannel(GetName(),"QwBeamMod",localModChannelID.fmodulename,"raw");
125 fModChannel.push_back(localchan2);
126 fModChannel[fModChannel.size()-1]->LoadChannelParameters(mapstr);
127 localModChannelID.fIndex=fModChannel.size()-1;
128 fModChannelID.push_back(localModChannelID);
129 /* localchan2->PrintInfo();
130 localchan2->PrintValue();
131 localchan2->PrintErrorCounters();*/
132 // Store reference to ramp channel
133 if (localModChannelID.fmodulename == "ramp") {
134 fRampChannelIndex = fModChannel.size() - 1;
135 }
136 }
137
138 if(ldebug)
139 {
140 localModChannelID.Print();
141 std::cout<<"line ok=";
142 if(lineok) std::cout<<"TRUE"<<std::endl;
143 else std::cout<<"FALSE"<<std::endl;
144 }
145
146 }
147
148 if (localModChannelID.fmoduletype == "SIS3801D32") {
149 wordsofar+=1;
150
151 if (lineok){
152 VQwHardwareChannel* localchan3 = new QwSIS3801D32_Channel();
153
154 localchan3->InitializeChannel(GetName(),"QwBeamMod",localModChannelID.fmodulename,"raw");
155 fModChannel.push_back(localchan3);
156 fModChannel[fModChannel.size()-1]->LoadChannelParameters(mapstr);
157 localModChannelID.fIndex=fModChannel.size()-1;
158 fModChannelID.push_back(localModChannelID);
159 /* localchan3->PrintInfo();
160 localchan3->PrintValue();
161 localchan3->PrintErrorCounters();*/
162 // Store reference to ramp channel
163 if (localModChannelID.fmodulename == "ramp") {
164 fRampChannelIndex = fModChannel.size() - 1;
165 }
166 }
167
168 if(ldebug)
169 {
170 localModChannelID.Print();
171 std::cout<<"line ok=";
172 if(lineok) std::cout<<"TRUE"<<std::endl;
173 else std::cout<<"FALSE"<<std::endl;
174 }
175 }
176
177 if(localModChannelID.fmoduletype=="SKIP"){
178 if (localModChannelID.modnum<=0) wordsofar+=1;
179 else wordsofar+= localModChannelID.modnum;
180
181 }
182
183 if(localModChannelID.fmoduletype == "WORD") {
184 // std::cout << "Decoding QwWord :: " << namech << std::endl;
185 QwWord localword;
186 localword.fSubbankIndex=bankindex;
187 localword.fWordInSubbank=wordsofar;
188 wordsofar+=1;
189 // I assume that one data = one word here. But it is not always the case, for
190 // example the triumf adc gives 6 words per channel
191 localword.fModuleType=localModChannelID.fmoduletype;
192 localword.fWordName=namech;
193 //localword.fWordType=dettype; // FIXME dettype is undefined so commented this out
194 fWord.push_back(localword);
195
196 // Store reference to pattern number
197 if (localword.fWordName == "bm_pattern_number") {
198 fPatternWordIndex = fWord.size() - 1;
199 }
200
201 if(namech=="ffb_status")//save the location of this word to access this later
202 fFFB_Index=fWord.size()-1;
203 //store reference to bmwobj
204 if(namech == "bmwobj"){
205 fBmwObj_Index = fWord.size() - 1;
206 }
207
208 fWordsPerSubbank[bankindex].second = fWord.size();
209 QwDebug << "--" << namech << "--" << fWord.size()-1 << QwLog::endl;
210 }
211
212 }
213 }
214
215 if(ldebug)
216 {
217 std::cout<<"Done with Load map channel \n";
218 for(size_t i=0;i<fModChannelID.size();i++)
219 fModChannelID[i].Print();
220
221 for(size_t i=0;i<fWord.size();i++)
222 fWord[i].PrintID();
223 }
224 return 0;
225}
226
227
228
229//QwModChannelID::QwModChannelID(Int_t subbankid, Int_t wordssofar,
230// TString name, TString modtype, QwBeamMod * obj):
231// fSubbankIndex(subbankid),fWordInSubbank(wordssofar),
232// fmoduletype(modtype),fmodulename(name),kUnknownDeviceType(-1)
233//{
234// fTypeID = kQwUnknownDeviceType;
235//}
236
238 QwParameterFile &paramfile):
239 fSubbankIndex(subbankid), fIndex(-1)
240{
241
242 fmoduletype = paramfile.GetTypedNextToken<TString>();
243 modnum = paramfile.GetTypedNextToken<Int_t>(); //slot number
244 channum = paramfile.GetTypedNextToken<Int_t>(); //channel number
245 // fdetectortype = paramfile.GetTypedNextToken<TString>(); //type-purpose of the detector
246 fmodulename = paramfile.GetTypedNextToken<TString>(); //name of the detector
247 fmoduletype.ToUpper();
248 // fdetectortype.ToLower();
249 fmodulename.ToLower();
250
251 Int_t offset;
252 if (fmoduletype == "VQWK") {
254 if (paramfile.ReturnValue("vqwk_buffer_offset",offset)) {
255 fWordInSubbank += offset;
256 }
257 // } else if (fmoduletype == "ADC18") {
258 //fWordInSubbank = QwADC18_Channel::GetBufferOffset(modnum, channum);
259 //if (paramfile.ReturnValue("adc18_buffer_offset",offset)) {
260 // fWordInSubbank += offset;
261 // }
262 } else if (fmoduletype == "SCALER"|| fmoduletype == "SIS3801D24" || fmoduletype == "SIS3801D32" ) {
264 if (paramfile.ReturnValue("scaler_buffer_offset",offset)) {
265 fWordInSubbank += offset;
266 }
267 } else {
268 fWordInSubbank = -1;
269 }
270}
271
272
273
274
275
276//*****************************************************************
277
279 fFFB_holdoff=0;//Default holdoff for the FFB pause
280}
281
282void QwBeamMod::LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag){
283 TString device_type = mapstr.GetTypedNextToken<TString>();
284 device_type.ToUpper();
285 TString device_name = mapstr.GetTypedNextToken<TString>();
286 device_name.ToUpper();
287
288 if (device_type == "VQWK"||device_type=="SCALER" ||device_type=="SIS3801D24" ||device_type=="SIS3801D32"){
289 device_name.ToLower();
290 Double_t LLX = mapstr.GetTypedNextToken<Double_t>(); //lower limit for BCM value
291 Double_t ULX = mapstr.GetTypedNextToken<Double_t>(); //upper limit for BCM value
292 varvalue = mapstr.GetTypedNextToken<TString>();//global/local
293 Double_t burplevel = mapstr.GetTypedNextToken<Double_t>();
294 varvalue.ToLower();
295 Double_t stabilitycut = mapstr.GetTypedNextToken<Double_t>();
296 QwMessage<<"QwBeamMod Error Code "<<GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut)<<QwLog::endl;
297 Int_t det_index=GetDetectorIndex(device_name);
298 QwMessage << "*****************************" << QwLog::endl;
299 QwMessage << " Type " << device_type << " Name " << device_name << " Index [" << det_index << "] "
300 << " device flag " << eventcut_flag << QwLog::endl;
301 fModChannel[det_index]->SetSingleEventCuts((GetGlobalErrorFlag(varvalue,eventcut_flag,stabilitycut)|kBModErrorFlag),LLX,ULX,stabilitycut,burplevel);
302 QwMessage << "*****************************" << QwLog::endl;
303
304 } else if (device_type == "WORD" && device_name== "FFB_STATUS"){
305 fFFB_holdoff=mapstr.GetTypedNextToken<UInt_t>();//Read the FFB OFF interval
306 } else if (device_type == "WORD" && device_name== "BMWOBJ"){
307 fBMWObj_LL=mapstr.GetTypedNextToken<Int_t>();
308 fBMWObj_UL=mapstr.GetTypedNextToken<Int_t>();
309 QwMessage << "bmwobj error cuts"
310 << "LowerLimit=" << fBMWObj_LL
311 << "UpperLimit=" << fBMWObj_UL <<QwLog::endl;
312 }
313}
314
315void QwBeamMod::LoadEventCuts_Fin(Int_t &eventcut_flag){
316 //update the event cut ON/OFF for all the devices
317 for (size_t i=0;i<fModChannel.size();i++)
318 fModChannel[i]->SetEventCutMode(eventcut_flag);
319}
320
321//*****************************************************************
323 Bool_t burpstatus = kFALSE;
324 VQwSubsystem* tmp = const_cast<VQwSubsystem *>(subsys);
325 if(Compare(tmp)) {
326 const QwBeamMod* input = dynamic_cast<const QwBeamMod*>(subsys);
327 for(size_t i=0;i<input->fModChannel.size();i++){
328 //QwError << "************* test Clock *****************" << QwLog::endl;
329 burpstatus |= (this->fModChannel[i])->CheckForBurpFail(input->fModChannel[i]);
330 }
331 }
332 return burpstatus;
333}
334
335//*****************************************************************
336Int_t QwBeamMod::LoadInputParameters(TString pedestalfile)
337{
338 Bool_t ldebug=kFALSE;
339 TString varname;
340
341 Int_t lineread=1;
342
343 if(ldebug)std::cout<<"QwBeamMod::LoadInputParameters("<< pedestalfile<<")\n";
344
345 QwParameterFile mapstr(pedestalfile.Data()); //Open the file
346 while (mapstr.ReadNextLine())
347 {
348 lineread+=1;
349 if(ldebug)std::cout<<" line read so far ="<<lineread<<"\n";
350 mapstr.TrimComment('!'); // Remove everything after a '!' character.
351 mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
352 if (mapstr.LineIsEmpty()) continue;
353 else
354 {
355 TString varname = mapstr.GetTypedNextToken<TString>(); //name of the channel
356 varname.ToLower();
357 varname.Remove(TString::kBoth,' ');
358 Double_t varped = mapstr.GetTypedNextToken<Double_t>(); // value of the pedestal
359 Double_t varcal = mapstr.GetTypedNextToken<Double_t>(); // value of the calibration factor
360 /* Double_t varweight = */ mapstr.GetTypedNextToken<Double_t>(); // value of the statistical weight
361
362 //if(ldebug) std::cout<<"inputs for channel "<<varname
363 // <<": ped="<<varped<<": cal="<<varcal<<": weight="<<varweight<<"\n";
364 for(size_t i=0;i<fModChannel.size();i++)
365 if(fModChannel[i]->GetElementName()==varname)
366 {
367 fModChannel[i]->SetPedestal(varped);
368 fModChannel[i]->SetCalibrationFactor(varcal);
369 break;
370 }
371 }
372
373 }
374 if(ldebug) std::cout<<" line read in the pedestal + cal file ="<<lineread<<" \n";
375
376 ldebug=kFALSE;
377 return 0;
378}
379
380
381Int_t QwBeamMod::ProcessEvBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words)
382{
383
384
385 Bool_t lkDEBUG=kFALSE;
386
387 Int_t index = GetSubbankIndex(roc_id,bank_id);
388
389
390 if (index>=0 && num_words>0){
391 // We want to process this ROC. Begin looping through the data.
392 if (lkDEBUG)
393 QwMessage << "QwBeamMod::ProcessEvBuffer: "
394 << "Begin processing ROC" << roc_id
395 << " and subbank "<<bank_id
396 << " number of words="<<num_words<<QwLog::endl;
397 if (buffer[0]==0xf0f0f0f0 && num_words%2==1){
398 buffer++;
399
400 if (lkDEBUG)
401 QwMessage << "QwBeamMod::ProcessEvBuffer: "
402 << "Skipped padding word 0xf0f0f0f0 at beginning of buffer."
403 << QwLog::endl;
404 }
405
406 for(size_t i=0;i<fModChannelID.size();i++) {
407 if(fModChannelID[i].fSubbankIndex==index) {
408 if (lkDEBUG) {
409 std::cout<<"found modulation data for "<<fModChannelID[i].fmodulename<<std::endl;
410 std::cout<<"word left to read in this buffer:"<<num_words-fModChannelID[i].fWordInSubbank<<std::endl;
411 }
412 fModChannel[fModChannelID[i].fIndex]->ProcessEvBuffer(&(buffer[fModChannelID[i].fWordInSubbank]),
413 num_words-fModChannelID[i].fWordInSubbank);
414 //fModChannel[fModChannelID[i].fIndex]->PrintInfo();
415 //fModChannel[fModChannelID[i].fIndex]->PrintValue();
416 //fModChannel[fModChannelID[i].fIndex]->PrintErrorCounters();
417 }
418 }
419
420 for(Int_t i=fWordsPerSubbank[index].first; i<fWordsPerSubbank[index].second; i++) {
421 if(fWord[i].fWordInSubbank+1<= (Int_t) num_words) {
422 fWord[i].fValue=buffer[fWord[i].fWordInSubbank];
423 //fWord[i].Print();
424 } else {
425 QwWarning << "QwBeamMod::ProcessEvBuffer: There is not enough word in the buffer to read data for "
426 << fWord[i].fWordName << QwLog::endl;
427 QwWarning << "QwBeamMod::ProcessEvBuffer: Words in this buffer:" << num_words
428 << " trying to read word number =" << fWord[i].fWordInSubbank << QwLog::endl;
429 }
430 }
431 if(lkDEBUG) {
432 QwDebug << "QwBeamMod::ProcessEvBuffer: Done with Processing this event" << QwLog::endl;
433 for(size_t i=0;i<fWord.size();i++) {
434 std::cout << " word number = " << i << " ";
435 fWord[i].Print();
436 }
437 }
438
439 }
440 lkDEBUG=kFALSE;
441
442 SetDataLoaded(kTRUE);
443
444 return 0;
445}
446
448 //std::cout << "Here in ApplySingleEventCuts() " << std::endl;
449 // std::cout << "ErrorFlag: " << fFFB_ErrorFlag << std::endl;
450 Bool_t test_Mod=kTRUE;
451 Bool_t test_BCM1=kTRUE;
452
453 for(size_t i=0;i<fModChannel.size();i++){
454 //std::cout<<" BCM ["<<i<<"] "<<std::endl;
455 //fModChannel[i]->PrintValue();
456 //if(fModChannel[i]->GetEventcutErrorFlag() >0){
457 // fModChannel[i]->PrintValue();
458 // }
459 test_BCM1 = fModChannel[i]->ApplySingleEventCuts();
460 test_Mod &= test_BCM1;
461
462 // fModChannel[i]->PrintValue();
463 //std::cout << "ErrorFlag: "<< fModChannel[i]->GetEventcutErrorFlag()<< std::endl;
464
465 //if(fModChannel[i]->GetEventcutErrorFlag() >0){
466 // fModChannel[i]->PrintValue();
467 //}
468
469
470 if(!test_BCM1 && bDEBUG) std::cout<<"******* QwBeamMod::SingleEventCuts()->BCM[ "<<i<<" , "<<fModChannel[i]->GetElementName()<<" ] ******" << std::endl;
471 }
472
473
474 /* if (fWord[fFFB_Index].fValue==8 && fFFB_Flag && fFFB_holdoff_Counter==0){
475 fFFB_holdoff_Counter=fFFB_holdoff;
476 fFFB_Flag=kFALSE;
477
478 }
479 */
480 // std::cout << "Upper Limit: " << fBMWObj_UL << " , " << "Lower Limit: " << fBMWObj_LL << std::endl;
481 //copy for fBmwObj
482 if (fBMWObj_LL < fBMWObj_UL) {
483 //std::cout << "Upper Limit: " << fBMWObj_UL << " , " << "Lower Limit: " << fBMWObj_LL << std::endl;
484
485 // Cuts are valid
486 if (fWord[fBmwObj_Index].GetValue()>fBMWObj_UL
487 || fWord[fBmwObj_Index].GetValue()<fBMWObj_LL){
488 // Value is outside of cuts range.
489
490 //std::cout << "bmwobj value: "<< fWord[fBmwObj_Index].GetValue() << " , " << "ErrorFlag: "<< fFFB_ErrorFlag << std::endl;
492 //std::cout << "ErrorFlag: " << fFFB_ErrorFlag << std::endl;
493 }
494 }
495
496 /*if (!fFFB_Flag && (fFFB_holdoff_Counter>0 && fFFB_holdoff_Counter<=fFFB_holdoff) ){
497 fFFB_ErrorFlag = (kGlobalCut+kBModErrorFlag+kBModFFBErrorFlag+kEventCutMode3);
498 fFFB_holdoff_Counter--;
499 if(fFFB_holdoff_Counter==0)
500 fFFB_Flag=kTRUE;
501 }
502 else
503 fFFB_ErrorFlag=0;*/
504
505 //std::cout << "WordValue: " << fWord[fBmwObj_Index].GetValue() << std::endl;
506 //std::cout << "ErrorFlag: " << fFFB_ErrorFlag << std::endl;
507 return test_Mod;
508
509}
510
512{
513 for(size_t i=0;i<fModChannel.size();i++){
514 fModChannel[i]->IncrementErrorCounters();
515 }
516}
517
518void QwBeamMod::PrintErrorCounters() const{//inherited from the VQwSubsystemParity; this will display the error summary
519
520 std::cout<<"*********QwBeamMod Error Summary****************"<<std::endl;
522 for(size_t i=0;i<fModChannel.size();i++){
523 fModChannel[i]->PrintErrorCounters();
524 }
526}
527
529{
530 const QwBeamMod* input = dynamic_cast<const QwBeamMod*> (ev_error);
531
532 for(size_t i=0;i<fModChannel.size();i++){
533 fModChannel[i]->UpdateErrorFlag(input->fModChannel[i]->GetEventcutErrorFlag());
534 }
535 // fFFB_ErrorFlag = input->fFFB_ErrorFlag;
536
537}
538
539UInt_t QwBeamMod::GetEventcutErrorFlag(){//return the error flag
540 UInt_t ErrorFlag;
541 ErrorFlag=0;
542 for(size_t i=0;i<fModChannel.size();i++){
543 ErrorFlag |= fModChannel[i]->GetEventcutErrorFlag();
544 }
545 ErrorFlag |= fFFB_ErrorFlag;
546
547 return ErrorFlag;
548}
549
550
551
552
554 //std::cout << "Here in ProcessEvent() " << std::endl;
555 for (size_t i = 0; i < fModChannel.size(); i++) {
556 // First apply HW checks and update HW error flags.
557 fModChannel[i]->ApplyHWChecks();
558 fModChannel[i]->ProcessEvent();
559 //fModChannel[i]->PrintValue();
560 }
561}
562
564{
565 // Fill histograms here to bypass event cuts
567}
568
569Int_t QwBeamMod::ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t* buffer, UInt_t num_words)
570{
571 return 0;
572}
573
574//*****************************************************************
575
577 // std::cout << "Here in ClearEventData " << std::endl;
578 for(size_t i=0;i<fModChannel.size();i++)
580
581 for (size_t i=0;i<fWord.size();i++)
583 //std::cout << "Error Flag: " << fFFB_ErrorFlag << std::endl;
585}
586
587//*****************************************************************
589{
590 Bool_t ldebug=kFALSE;
591 if(ldebug)
592 {
593 std::cout<<"QwBeamMod::GetDetectorIndex\n";
594 std::cout<<fModChannelID.size()<<" already registered detector\n";
595 }
596
597 Int_t result=-1;
598 for(size_t i=0;i<fModChannelID.size();i++)
599 {
600
601 if(fModChannelID[i].fmodulename==name){
602 result=fModChannelID[i].fIndex;
603 }
604
605 if(ldebug)
606 {
607 std::cout<<"testing against ("<<fModChannelID[i].fTypeID
608 <<","<<fModChannelID[i].fmodulename<<")=>"<<result<<"\n";
609 }
610 }
611 //do it over the words, fWords
612 for(size_t i=0;i<fWord.size();i++) {
613
614 if(fWord[i].fWordName==name){
615 result=i;
616 }
617 }
618
619 return result;
620}
621//*****************************************************************
623{
624
625 if(Compare(value))
626 {
627
629 QwBeamMod* input = dynamic_cast<QwBeamMod*>(value);
630
631 for(size_t i=0;i<input->fModChannel.size();i++){
632 //input->fModChannel[i]->PrintValue();
633 this->fModChannel[i]->AssignValueFrom(input->fModChannel[i]);
634 //fModChannel[i]->PrintValue();
635 }
636 for(size_t i=0;i<input->fWord.size();i++)
637 this->fWord[i].fValue=input->fWord[i].fValue;
638 this->fFFB_ErrorFlag=input->fFFB_ErrorFlag;
639
640 }
641 return *this;
642}
643
645{
646 if(Compare(value)){
647 // std::cout << "Here in operator+= " << std::endl;
648 //QwBeamMod* input= (QwBeamMod*)value ;
649 QwBeamMod* input = dynamic_cast<QwBeamMod*>(value);
650
651 for(size_t i=0;i<input->fModChannel.size();i++)
652 *(this->fModChannel[i]) += *(input->fModChannel[i]);
653// for(size_t i=0;i<input->fWord.size();i++)
654// this->fWord[i]+=input->fWord[i];
655 this->fFFB_ErrorFlag |= input->fFFB_ErrorFlag;
656 // std::cout <<"Error Flag: " << this->fFFB_ErrorFlag <<std::endl;
657 }
658 return *this;
659}
660
662
663 if(Compare(value))
664 {
665 //QwBeamMod* input= (QwBeamMod*)value;
666 QwBeamMod* input = dynamic_cast<QwBeamMod*>(value);
667
668 for(size_t i=0;i<input->fModChannel.size();i++)
669 *(this->fModChannel[i]) -= *(input->fModChannel[i]);
670// for(size_t i=0;i<input->fWord.size();i++)
671// this->fWord[i]-=input->fWord[i];
672 this->fFFB_ErrorFlag |= input->fFFB_ErrorFlag;
673
674 }
675 return *this;
676}
677
679{
680 if(Compare(numer)&&Compare(denom))
681 {
682 //QwBeamMod* innumer= (QwBeamMod*)numer ;
683 QwBeamMod* innumer = dynamic_cast<QwBeamMod*>(numer);
684 //QwBeamMod* indenom= (QwBeamMod*)denom ;
685 QwBeamMod* indenom = dynamic_cast<QwBeamMod*>(denom);
686
687 for(size_t i=0;i<innumer->fModChannel.size();i++)
688 this->fModChannel[i]->Ratio(innumer->fModChannel[i],indenom->fModChannel[i]);
689 for(size_t i=0;i<innumer->fWord.size();i++)
690 this->fWord[i].fValue=innumer->fWord[i].fValue;
691
692 }
693}
694
695
696void QwBeamMod::Scale(Double_t factor)
697{
698 for(size_t i=0;i<fModChannel.size();i++)
699 fModChannel[i]->Scale(factor);
700}
701
703
704void QwBeamMod::AccumulateRunningSum(VQwSubsystem*, Int_t count, Int_t ErrorMask) { }
705
707{
708 Bool_t res = kTRUE;
709 if(typeid(*value)!=typeid(*this))
710 {
711 res=kFALSE;
712 // std::cout<<" types are not ok \n";
713 // std::cout<<" this is bypassed just for now but should be fixed eventually \n";
714 }
715 else
716 {
717 QwBeamMod* input = dynamic_cast<QwBeamMod*>(value);
718
719 if(input->fModChannel.size()!=fModChannel.size())
720 {
721 res=kFALSE;
722 // std::cout<<" not the same number of bcms \n";
723 }
724 }
725 return res;
726}
727
728
729//*****************************************************************
730
731void QwBeamMod::ConstructHistograms(TDirectory *folder, TString &prefix)
732{
733 // Go to the specified folder
734 if (folder != NULL) folder->cd();
735
736 // No histogram creation for asym, yield, diff, etc
737 if (prefix != "") return;
738
739 // Beam modulation channels
740 //for (size_t i = 0; i < fModChannel.size(); i++)
741 // fModChannel[i]->ConstructHistograms(folder,prefix);
742}
743
745{
746 // No data loaded
747 if (! HasDataLoaded()) return;
748
749 // No histograms defined
750 if (fHistograms.size() == 0) return;
751
752 // Beam modulation channels
753 //for (size_t i = 0; i < fModChannel.size(); i++)
754 // fModChannel[i]->FillHistograms();
755
756 // Do we have a ramp channel?
758
759 // Do we have a pattern word?
760 if (fPatternWordIndex < 0 || fPatternWordIndex > Int_t(fWord.size())) return;
761
762 // Determine the ramp/phase
763 Double_t ramp = fModChannel[fRampChannelIndex]->GetValue();
764 if (ramp < 0) return;
765
766 // Determine the pattern number -- the pattern number for single coil is
767 // between [0, 4] so we need to check for this.
768 Int_t pattern = -1;
769
770 if(fWord[fPatternWordIndex].fValue < 11 && fWord[fPatternWordIndex].fValue > 0)
771 pattern = fWord[fPatternWordIndex].fValue;
772 else
773 pattern = fWord[fPatternWordIndex].fValue - 11;
774
775 if (pattern < 0 || pattern > 4) return;
776
777 // Fill histograms for all BPMs and each of the modulation patterns
778 //
779 // Due to the the way the ADC averages the ramp signal we want to filter
780 // out events at the edged of the signal.
781 //
782 // Separated the ramp cut here because it is ridiculously long...
783 //
784
785 [[maybe_unused]]
786 Double_t ramp_block_41 = fModChannel[fRampChannelIndex]->GetValue(4) + fModChannel[fRampChannelIndex]->GetValue(1);
787 [[maybe_unused]]
788 Double_t ramp_block_32 = fModChannel[fRampChannelIndex]->GetValue(3) + fModChannel[fRampChannelIndex]->GetValue(2);
789 [[maybe_unused]]
790 Double_t ramp_block = ramp_block_41 - ramp_block_32;
791
792}
793
794
795void QwBeamMod::ConstructBranchAndVector(TTree *tree, TString & prefix, QwRootTreeBranchVector &values)
796{
797 TString basename;
798
799 for(size_t i = 0; i < fModChannel.size(); i++){
800 fModChannel[i]->ConstructBranchAndVector(tree, prefix, values);
801 }
802 // for (size_t i=0;i<fWord.size();i++)
803 // fWord[i].ConstructBranchAndVector(tree, prefix, values);
804 fTreeArrayIndex = values.size();
805 for (size_t i=0; i<fWord.size(); i++) {
806 // basename = fWord[i].fWordName;
807 basename = prefix(0, (prefix.First("|") >= 0)? prefix.First("|"): prefix.Length());
808 basename += fWord[i].fWordName;
809 values.push_back(basename.Data(), 'D');
810 tree->Branch(basename, &(values[fTreeArrayIndex + i]), values.LeafList(fTreeArrayIndex + i).c_str());
811 }
812}
813
815{
816 //std::cout << "inside FillTreeVector"<< std::endl;
817 //std::cout << "fTreeArrayIndex: " << fTreeArrayIndex << std::endl;
818 size_t index = fTreeArrayIndex;
819
820 for (size_t i = 0; i < fModChannel.size(); i++){
821 fModChannel[i]->FillTreeVector(values);
822 // fModChannel[i]->PrintValue();
823 }
824 for (size_t i = 0; i < fWord.size(); i++){
825 values.SetValue(index++, fWord[i].fValue);
826 }
827 //for (size_t i=0; i<values.size(); i++){
828 // std::cout << values[i] << " ";
829 // if ((i+1)%8==0) std::cout << std::endl;
830 // }
831
832}
833
834#ifdef HAS_RNTUPLE_SUPPORT
835void QwBeamMod::ConstructNTupleAndVector(std::unique_ptr<ROOT::RNTupleModel>& model, TString& prefix, std::vector<Double_t>& values, std::vector<std::shared_ptr<Double_t>>& fieldPtrs)
836{
837 TString basename;
838
839 for(size_t i = 0; i < fModChannel.size(); i++){
840 fModChannel[i]->ConstructNTupleAndVector(model, prefix, values, fieldPtrs);
841 }
842
843 fTreeArrayIndex = values.size();
844 for (size_t i=0; i<fWord.size(); i++) {
845 basename = prefix(0, (prefix.First("|") >= 0)? prefix.First("|"): prefix.Length());
846 basename += fWord[i].fWordName;
847 values.push_back(0.0);
848 fieldPtrs.push_back(model->MakeField<Double_t>(basename.Data()));
849 }
850}
851
852void QwBeamMod::FillNTupleVector(std::vector<Double_t>& values) const
853{
854 size_t index = fTreeArrayIndex;
855
856 for (size_t i = 0; i < fModChannel.size(); i++){
857 fModChannel[i]->FillNTupleVector(values);
858 }
859 for (size_t i = 0; i < fWord.size(); i++){
860 values[index++] = fWord[i].fValue;
861 }
862}
863#endif // HAS_RNTUPLE_SUPPORT
864
865
866//*****************************************************************
868{
869 QwMessage << "Name of the subsystem =" << fSystemName << QwLog::endl;
870
871 QwMessage << "there are " << fModChannel.size() << " mods" << QwLog::endl;
872
873 QwMessage << " Printing Running AVG and other channel info for fModChannel" << QwLog::endl;
874 for(size_t i=0;i<fModChannel.size();i++)
876 for(size_t i=0;i<fWord.size();i++)
877 fWord[i].Print();
878}
879
881{
882 for (size_t i=0;i<fModChannelID.size();i++)
883 {
884 std::cout<<"============================="<<std::endl;
885 std::cout<<" Detector ID="<<i<<std::endl;
886 fModChannelID[i].Print();
887 }
888}
889
890
892{
893 QwMessage <<std::endl<<"Detector name= "<<fmodulename<<QwLog::endl;
894 QwMessage <<"SubbankkIndex= "<<fSubbankIndex<<QwLog::endl;
895 QwMessage <<"word index in subbank= "<<fWordInSubbank<<QwLog::endl;
896 QwMessage <<"module type= "<<fmoduletype<<QwLog::endl;
897// QwMessage <<"detector type= "<<fdetectortype<<" that is index="<<fTypeID<<QwLog::endl;
898 QwMessage <<"Index of this detector in the vector of similar detector= "<<
900//QwMessage <<"Subelement index= "<< fSubelement<<QwLog::endl;
901
902
903 QwMessage << "---------------------------------------------------" << QwLog::endl;
905}
906
907
908#ifdef __USE_DATABASE__
909//*****************************************************************
910void QwBeamMod::FillDB_MPS(QwParityDB *db, TString datatype)
911{
912 Bool_t local_print_flag = false;
913
914 if(local_print_flag) {
915 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
916 QwMessage << " QwBeamMod::FillDB_MPS " << QwLog::endl;
917 QwMessage << " --------------------------------------------------------------- " << QwLog::endl;
918 }
919
920 std::vector<QwParitySchema::beam_optics_row> entrylist;
921
922 UInt_t analysis_id = db->GetAnalysisID();
923
924 for(size_t bpm = 0; bpm < fModChannel.size(); bpm++) {
925 for(size_t pattern = 0; pattern < 5; pattern++) {
926 // Explicitly zero the beam optics ID to ensure a non-sensical default
927 // is not picked up.
928 QwParitySchema::beam_optics table;
929 QwParitySchema::beam_optics_row row;
930 row[table.analysis_id] = analysis_id;
931 row[table.monitor_id] = 0; // placeholder
932 row[table.modulation_type_id] = pattern;
933 row[table.n] = 0; // placeholder
934 row[table.offset] = 0.0; // placeholder
935 row[table.amplitude] = 0.0; // placeholder
936 row[table.phase] = 0.0; // placeholder
937 row[table.o_error] = 0.0; // placeholder
938 row[table.a_error] = 0.0; // placeholder
939 row[table.p_error] = 0.0; // placeholder
940 row[table.gof_para] = 0.0; // placeholder
941
942 entrylist.push_back(row);
943 }
944 }
945
946 if(local_print_flag) {
947 QwMessage << QwColor(Qw::kGreen) << "Entrylist Size : "
948 << QwColor(Qw::kBoldRed) << entrylist.size()
949 << QwColor(Qw::kNormal) << QwLog::endl;
950 }
951
952 if( entrylist.size() ) {
953 auto c = db->GetScopedConnection();
954 for (const auto& entry : entrylist) {
955 c->QueryExecute(entry.insert_into());
956 }
957 }
958 else {
959 QwMessage << "QwBeamMod::FillDB_MPS :: Nothing to insert in database." << QwLog::endl;
960 }
961}
962
963void QwBeamMod::FillDB(QwParityDB *db, TString datatype)
964{
965 return;
966}
967
968void QwBeamMod::FillErrDB(QwParityDB *db, TString datatype)
969{
970 return;
971}
972#endif // __USE_DATABASE__
973
975{
976 return;
977};
Helper functions and utilities for ROOT histogram management.
Base and derived classes for scaler channel data handling.
class QwScaler_Channel< 0x00ffffff, 0 > QwSIS3801D24_Channel
class QwScaler_Channel< 0xffffffff, 0 > QwSIS3801D32_Channel
A logfile class, based on an identical class in the Hermes analyzer.
#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.
static const UInt_t kGlobalCut
Definition QwTypes.h:182
static const UInt_t kBModErrorFlag
Definition QwTypes.h:173
ULong64_t BankID_t
Definition QwTypes.h:21
static const UInt_t kEventCutMode3
Definition QwTypes.h:174
static const UInt_t kBModFFBErrorFlag
Definition QwTypes.h:172
UInt_t GetGlobalErrorFlag(TString evtype, Int_t evMode, Double_t stabilitycut)
Definition QwTypes.cc:132
UInt_t ROCID_t
Definition QwTypes.h:20
Beam modulation subsystem for parity analysis.
@ kBoldRed
Definition QwColor.h:79
@ kGreen
Definition QwColor.h:77
@ kNormal
Definition QwColor.h:81
std::vector< TH1_ptr > fHistograms
Histograms associated with this data element.
static std::ostream & endl(std::ostream &)
End of the line.
Definition QwLog.cc:297
Command-line and configuration file options processor.
Definition QwOptions.h:141
Configuration file parser with flexible tokenization and search capabilities.
T GetTypedNextToken()
Get next token into specific type.
void TrimWhitespace(TString::EStripType head_tail=TString::kBoth)
Bool_t HasVariablePair(const std::string &separatorchars, std::string &varname, std::string &varvalue)
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 ReturnValue(const std::string keyname, T &retvalue)
A helper class to manage a vector of branch entries for ROOT trees.
Definition QwRootFile.h:55
size_type size() const noexcept
Definition QwRootFile.h:83
std::string LeafList(size_type start_index=0) const
Definition QwRootFile.h:230
void push_back(const std::string &name, const char type='D')
Definition QwRootFile.h:197
void SetValue(size_type index, Double_t val)
Definition QwRootFile.h:110
static Int_t GetBufferOffset(Int_t scalerindex, Int_t wordindex, UInt_t header=1)
Concrete hardware channel for VQWK ADC modules (6x32-bit words)
static void PrintErrorCounterTail()
static Int_t GetBufferOffset(Int_t moduleindex, Int_t channelindex)
static void PrintErrorCounterHead()
Word-level data manipulation and bit operations.
Definition QwWord.h:31
Int_t fWordInSubbank
Definition QwWord.h:38
Int_t fSubbankIndex
Definition QwWord.h:37
TString fModuleType
Definition QwWord.h:39
TString fWordName
Definition QwWord.h:40
Abstract base for concrete hardware channels implementing dual-operator pattern.
void InitializeChannel(TString name)
Initialize the fields in this object.
BankID_t fCurrentBank_ID
Bank ID (and Marker word) that is currently being processed;.
TString fSystemName
Name of this subsystem.
Int_t GetSubbankIndex() const
virtual VQwSubsystem & operator=(VQwSubsystem *value)
Assignment Note: Must be called at the beginning of all subsystems routine call to operator=(VQwSubsy...
void RegisterRocBankMarker(QwParameterFile &mapstr)
TString GetName() const
std::map< TString, TString > fDetectorMaps
Map of file name to full path or content.
VQwSubsystem(const TString &name)
Constructor with name.
void SetDataLoaded(Bool_t flag)
ROCID_t fCurrentROC_ID
ROC ID that is currently being processed.
Bool_t HasDataLoaded() const
Mapping information for beam modulation channels.
Definition QwBeamMod.h:45
Int_t fWordInSubbank
Definition QwBeamMod.h:60
QwModChannelID(Int_t subbankid, Int_t wordssofar, TString name, TString modtype, QwBeamMod *obj)
TString fmodulename
Definition QwBeamMod.h:67
Int_t fSubbankIndex
Definition QwBeamMod.h:59
TString fmoduletype
Definition QwBeamMod.h:66
void ClearEventData() override
Definition QwBeamMod.cc:576
void LoadEventCuts_Line(QwParameterFile &mapstr, TString &varvalue, Int_t &eventcut_flag) override
Definition QwBeamMod.cc:282
static const Bool_t bDEBUG
Definition QwBeamMod.h:222
void LoadEventCuts_Init() override
Definition QwBeamMod.cc:278
VQwSubsystem & operator-=(VQwSubsystem *value) override
Definition QwBeamMod.cc:661
Bool_t CheckForBurpFail(const VQwSubsystem *subsys) override
Report the number of events failed due to HW and event cut failures.
Definition QwBeamMod.cc:322
void PrintErrorCounters() const override
Definition QwBeamMod.cc:518
Int_t fRampChannelIndex
Definition QwBeamMod.h:224
UInt_t fBmwObj_Index
Definition QwBeamMod.h:226
void PrintModChannelID()
Definition QwBeamMod.cc:880
VQwSubsystem & operator=(VQwSubsystem *value) override
Assignment Note: Must be called at the beginning of all subsystems routine call to operator=(VQwSubsy...
Definition QwBeamMod.cc:622
VQwSubsystem & operator+=(VQwSubsystem *value) override
Definition QwBeamMod.cc:644
Int_t fBMWObj_UL
Definition QwBeamMod.h:228
Int_t ProcessConfigurationBuffer(const ROCID_t roc_id, const BankID_t bank_id, UInt_t *buffer, UInt_t num_words) override
Definition QwBeamMod.cc:569
Int_t fBMWObj_LL
Definition QwBeamMod.h:227
UInt_t fFFB_holdoff
Definition QwBeamMod.h:218
void ProcessOptions(QwOptions &options) override
Process the command line options.
Definition QwBeamMod.cc:39
void AccumulateRunningSum(VQwSubsystem *, Int_t count=0, Int_t ErrorMask=0xFFFFFFF) override
Update the running sums for devices.
Definition QwBeamMod.cc:704
Int_t LoadChannelMap(TString mapfile) override
Mandatory map file definition.
Definition QwBeamMod.cc:43
void FillHistograms() override
Fill the histograms for this subsystem.
Definition QwBeamMod.cc:744
void ProcessEvent_2() override
Process the event data again, including data from other subsystems. Not all derived classes will requ...
Definition QwBeamMod.cc:563
Int_t LoadInputParameters(TString pedestalfile) override
Mandatory parameter file definition.
Definition QwBeamMod.cc:336
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
UInt_t fFFB_ErrorFlag
Definition QwBeamMod.h:220
void CalculateRunningAverage() override
Calculate the average for all good events.
Definition QwBeamMod.cc:702
void WritePromptSummary(QwPromptSummary *ps, TString type) override
Definition QwBeamMod.cc:974
void Ratio(VQwSubsystem *numer, VQwSubsystem *denom) override
Definition QwBeamMod.cc:678
void FillTreeVector(QwRootTreeBranchVector &values) const override
Fill the tree vector.
Definition QwBeamMod.cc:814
Int_t fPatternWordIndex
Definition QwBeamMod.h:225
std::vector< VQwHardwareChannel * > fModChannel
Definition QwBeamMod.h:208
UInt_t GetEventcutErrorFlag() override
Return the error flag to the top level routines related to stability checks and ErrorFlag updates.
Definition QwBeamMod.cc:539
Int_t GetDetectorIndex(TString name)
Definition QwBeamMod.cc:588
void Scale(Double_t factor) override
Definition QwBeamMod.cc:696
Bool_t Compare(VQwSubsystem *source)
Definition QwBeamMod.cc:706
std::vector< QwWord > fWord
Definition QwBeamMod.h:210
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...
Definition QwBeamMod.cc:381
void ConstructBranchAndVector(TTree *tree, TString &prefix, QwRootTreeBranchVector &values) override
Construct the branch and tree vector.
Definition QwBeamMod.cc:795
void ProcessEvent() override
Definition QwBeamMod.cc:553
void Print()
Definition QwBeamMod.cc:867
UInt_t fFFB_Index
Definition QwBeamMod.h:217
void IncrementErrorCounters() override
Increment the error counters.
Definition QwBeamMod.cc:511
std::vector< std::pair< Int_t, Int_t > > fWordsPerSubbank
Definition QwBeamMod.h:211
std::vector< QwModChannelID > fModChannelID
Definition QwBeamMod.h:209
Bool_t ApplySingleEventCuts() override
Apply the single event cuts.
Definition QwBeamMod.cc:447
Int_t fTreeArrayIndex
Definition QwBeamMod.h:206
void LoadEventCuts_Fin(Int_t &eventcut_flag) override
Definition QwBeamMod.cc:315
virtual void PrintValue() const
Print values of all channels.
virtual void FillDB_MPS(QwParityDB *, TString)
Fill the database with MPS-based variables Note that most subsystems don't need to do this.
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)