Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcScalerEvtHandler.cxx
Go to the documentation of this file.
1
36#include "THcScalerEvtHandler.h"
37//#include "GenScaler.h"
38#include "Scaler3800.h"
39#include "Scaler3801.h"
40#include "Scaler1151.h"
41#include "Scaler560.h"
42#include "Scaler9001.h"
43#include "Scaler9250.h"
44#include "THaEvData.h"
45#include "THcParmList.h"
46#include "THcGlobals.h"
47#include "THaGlobals.h"
48#include "TMath.h"
49#include "TString.h"
50#include "TROOT.h"
51#include <cstring>
52#include <cstdio>
53#include <cstdlib>
54#include <iostream>
55#include <map>
56#include <iterator>
57#include "THaVarList.h"
58#include "VarDef.h"
59#include "Helper.h"
60#include "Textvars.h" // for Podd::vsplit
61#include "THaString.h"
62
63using namespace std;
64using namespace Decoder;
66
67static const UInt_t ICOUNT = 1;
68static const UInt_t IRATE = 2;
69static const UInt_t ICURRENT = 3;
70static const UInt_t ICHARGE = 4;
71static const UInt_t ITIME = 5;
72static const UInt_t ICUT = 6;
73static const UInt_t MAXCHAN = 32;
74static const UInt_t defaultDT = 4;
75
76THcScalerEvtHandler::THcScalerEvtHandler(const char *name, const char* description)
77 : THaEvtTypeHandler(name,description),
78 fBCM_Gain(0), fBCM_Offset(0), fBCM_SatOffset(0), fBCM_SatQuadratic(0), fBCM_delta_charge(0),
79 evcount(0), evcountR(0.0), ifound(0), fNormIdx(-1),
80 fNormSlot(-1),
81 dvars(0),dvars_prev_read(0), dvarsFirst(0), fScalerTree(0), fUseFirstEvent(kTRUE),
82 fOnlySyncEvents(kFALSE), fOnlyBanks(kFALSE), fDelayedType(-1),
83 fClockChan(-1), fLastClock(0), fClockOverflows(0)
84{
85 fRocSet.clear();
86 fModuleSet.clear();
87 scal_prev_read.clear();
88 scal_present_read.clear();
89 scal_overflows.clear();
90}
91
93{
94 // The tree object is owned by ROOT since it gets associated wth the output
95 // file, so DO NOT delete it here.
96 if (!TROOT::Initialized()) {
97 delete fScalerTree;
98 }
99 Podd::DeleteContainer(scalers);
100 Podd::DeleteContainer(scalerloc);
101 delete [] dvars_prev_read;
102 delete [] dvars;
103 delete [] dvarsFirst;
104 delete [] fBCM_Gain;
105 delete [] fBCM_Offset;
106 delete [] fBCM_SatOffset;
107 delete [] fBCM_SatQuadratic;
108 delete [] fBCM_delta_charge;
109
110 for( auto& fDelayedEvent: fDelayedEvents )
111 delete[] fDelayedEvent;
112 fDelayedEvents.clear();
113}
114
116{
117 // Process any delayed events in order received
118
119 cout << "THcScalerEvtHandler::End Analyzing " << fDelayedEvents.size() << " delayed scaler events" << endl;
120 for( auto rdata: fDelayedEvents ) {
121 AnalyzeBuffer(rdata, kFALSE);
122 }
123 if (fDebugFile) *fDebugFile << "scaler tree ptr "<<fScalerTree<<endl;
124 evNumber += 1;
127
128 for( auto& fDelayedEvent: fDelayedEvents )
129 delete[] fDelayedEvent;
130 fDelayedEvents.clear();
131
133 return 0;
134}
135
136
138{
139 char prefix[2];
140 prefix[0]='g';
141 prefix[1]='\0';
142 fNumBCMs = 0;
143 DBRequest list[]={
144 {"NumBCMs",&fNumBCMs, kInt, 0, true},
145 {nullptr}
146 };
147 gHcParms->LoadParmValues((DBRequest*)&list, prefix);
148 //cout << " NUmber of BCMs = " << fNumBCMs << endl;
149 //
150 if(fNumBCMs > 0) {
156 string bcm_namelist;
157 DBRequest list2[]={
158 {"BCM_Gain", fBCM_Gain, kDouble, (UInt_t) fNumBCMs},
159 {"BCM_Offset", fBCM_Offset, kDouble,(UInt_t) fNumBCMs},
160 {"BCM_SatQuadratic", fBCM_SatQuadratic, kDouble,(UInt_t) fNumBCMs,true},
161 {"BCM_SatOffset", fBCM_SatOffset, kDouble,(UInt_t) fNumBCMs,true},
162 {"BCM_Names", &bcm_namelist, kString},
163 {"BCM_Current_threshold", &fbcm_Current_Threshold, kDouble,0, true},
164 {"BCM_Current_threshold_index", &fbcm_Current_Threshold_Index, kInt,0,true},
165 {nullptr}
166 };
169 for(Int_t i=0;i<fNumBCMs;i++) {
170 fBCM_SatOffset[i]=0.;
171 fBCM_SatQuadratic[i]=0.;
172 }
173 gHcParms->LoadParmValues((DBRequest*)&list2, prefix);
174 vector<string> bcm_names = Podd::vsplit(bcm_namelist);
175 for(Int_t i=0;i<fNumBCMs;i++) {
176 fBCM_Name.push_back(bcm_names[i]+".scal");
177 fBCM_delta_charge[i]=0.;
178 }
179 }
180 fTotalTime=0.;
182 fDeltaTime=-1.;
183 //
184 //
185 return kOK;
186}
196 fDelayedType = evtype;
197}
198
200{
201 Int_t lfirst=1;
202
203 if(evdata->GetEvNum() > 0) {
204 evNumber=evdata->GetEvNum();
206 }
207 if ( !IsMyEvent(evdata->GetEvType()) ) return -1;
208
209 if (fDebugFile) {
210 *fDebugFile << endl << "---------------------------------- "<<endl<<endl;
211 *fDebugFile << "\nEnter THcScalerEvtHandler for fName = "<<fName<<endl;
212 EvDump(evdata);
213 }
214
215 if (lfirst && !fScalerTree) {
216
217
218 lfirst = 0; // Can't do this in Init for some reason
219
220 TString sname1 = "TS";
221 TString sname2 = sname1 + fName;
222 TString sname3 = fName + " Scaler Data";
223
224 if (fDebugFile) {
225 *fDebugFile << "\nAnalyze 1st time for fName = "<<fName<<endl;
226 *fDebugFile << sname2 << " " <<sname3<<endl;
227 }
228
229 fScalerTree = new TTree(sname2.Data(),sname3.Data());
230 fScalerTree->SetAutoSave(200000000);
231
232 TString name, tinfo;
233
234 name = "evcount";
235 tinfo = name + "/D";
236 fScalerTree->Branch(name.Data(), &evcountR, tinfo.Data(), 4000);
237
238 name = "evNumber";
239 tinfo = name + "/D";
240 fScalerTree->Branch(name.Data(), &evNumberR, tinfo.Data(), 4000);
241
242 for (size_t i = 0; i < scalerloc.size(); i++) {
243 name = scalerloc[i]->name;
244 tinfo = name + "/D";
245 fScalerTree->Branch(name.Data(), &dvars[i], tinfo.Data(), 4000);
246 }
247
248 } // if (lfirst && !fScalerTree)
249
250 auto *rdata = (UInt_t*) evdata->GetRawDataBuffer();
251
252 if( (Int_t)evdata->GetEvType() == fDelayedType) { // Save this event for processing later
253 UInt_t evlen = evdata->GetEvLength();
254
255 auto *datacopy = new UInt_t[evlen];
256 fDelayedEvents.push_back(datacopy);
257 memcpy(datacopy,rdata,evlen*sizeof(UInt_t));
258 return 1;
259 } else { // A normal event
260 if (fDebugFile) *fDebugFile<<"\n\nTHcScalerEvtHandler :: Debugging event type "<<dec<<evdata->GetEvType()<< " event num = " << evdata->GetEvNum() << endl<<endl;
261 Int_t ret;
262 if((ret=AnalyzeBuffer(rdata,fOnlySyncEvents))) {
263 if (fDebugFile) *fDebugFile << "scaler tree ptr "<<fScalerTree<<endl;
265 }
266 return ret;
267
268 }
269
270}
272{
273
274 // Parse the data, load local data arrays.
275 auto *p = (UInt_t*) rdata;
276
277 UInt_t *plast = p+*p; // Index to last word in the bank
278
279 ifound=0;
280 while(p<plast) {
281 p++; // point to header
282 if (fDebugFile) {
283 *fDebugFile << "Bank: " << hex << *p << dec << " len: " << *(p-1) << endl;
284 }
285 if((*p & 0xff00) == 0x1000) { // Bank Containing banks
286 p++; // Now pointing to a bank in the bank
287 } else if (((*p & 0xff00) == 0x100) && (*p != 0xC0000100)) {
288 // Bank containing integers. Look for scalers
289 // This is either ROC bank containing integers or
290 // a bank within a ROC containing data from modules of a single type
291 // Look for scaler data
292 // Assume that very first word is a scaler header
293 // At any point in the bank where the word is not a matching
294 // header, we stop.
295 UInt_t tag = (*p>>16) & 0xffff;
296 UInt_t num = (*p) & 0xff;
297 UInt_t *pnext = p+*(p-1); // Next bank
298 p++; // First data word
299
300 // Skip over banks that can't contain scalers
301 // If SetOnlyBanks(kTRUE) called, fRocSet will be empty
302 // so only bank tags matching module types will be considered.
303 if(fModuleSet.find(tag)!=fModuleSet.end()) {
304 if(onlysync && num==0) {
305 ifound = 0;
306 return 0;
307 }
308 } else if (fRocSet.find(tag)==fRocSet.end()) {
309 p = pnext; // Fall through to end of the above else if
310 }
311
312 // Look for normalization scaler module first.
313 if(fNormIdx >= 0) {
314 UInt_t *psave = p;
315 while(p < pnext) {
316 if(scalers[fNormIdx]->IsSlot(*p)) {
317 scalers[fNormIdx]->Decode(p);
318 ifound = 1;
319 break;
320 }
321 p += scalers[fNormIdx]->GetNumChan() + 1;
322 }
323 p = psave;
324 }
325 while(p < pnext) {
326 Int_t nskip = 0;
327 if(fDebugFile) {
328 *fDebugFile << "Scaler Header: " << hex << *p << dec;
329 }
330 for(size_t j=0; j<scalers.size(); j++) {
331 if(scalers[j]->IsSlot(*p)) {
332 nskip = scalers[j]->GetNumChan() + 1;
333 if((Int_t) j != fNormIdx) {
334 if(fDebugFile) {
335 *fDebugFile << " found (" << j << ") skip " << nskip << endl;
336 }
337 scalers[j]->Decode(p);
338 ifound = 1;
339 }
340 break;
341 }
342 }
343 if(nskip == 0) {
344 if(fDebugFile) {
345 *fDebugFile << endl;
346 }
347 break; // Didn't find a matching header
348 }
349 p = p + nskip;
350 }
351 p = pnext;
352 } else {
353 p = p+*(p-1); // Skip to next bank
354 }
355 }
356
357 if (fDebugFile) {
358 *fDebugFile << "Finished with decoding. "<<endl;
359 *fDebugFile << " Found flag = "<<ifound<<endl;
360 }
361
362 // HMS has headers which are different from SOS, but both are
363 // event type 0 and come here. If you found no headers, return.
364
365 if (!ifound) return 0;
366
367 // The correspondance between dvars and the scaler and the channel
368 // will be driven by a scaler.map file -- later
369 Double_t scal_current=0;
370 UInt_t thisClock = scalers[fNormIdx]->GetData(fClockChan);
371 if(thisClock < fLastClock) { // Count clock scaler wrap arounds
373 }
375 fLastClock = thisClock;
377 if (fDeltaTime==0) {
378 cout << " ******************* Severe Warning ****************************" << endl;
379 cout << " In THcScalerEvtHandler have found fDeltaTime is zero !! " << endl;
380 cout << " ******************* Alert DAQ experts ****************************" << endl;
381 }
383 Int_t nscal=0;
384 for (size_t i = 0; i < scalerloc.size(); i++) {
385 size_t ivar = scalerloc[i]->ivar;
386 size_t idx = scalerloc[i]->index;
387 size_t ichan = scalerloc[i]->ichan;
388 if (evcount==0) {
389 if (fDebugFile) *fDebugFile << "Debug dvarsFirst "<<i<<" "<<ivar<<" "<<idx<<" "<<ichan<<endl;
390 if ((ivar < scalerloc.size()) &&
391 (idx < scalers.size()) &&
392 (ichan < MAXCHAN)){
393 if(fUseFirstEvent){
394 if (scalerloc[ivar]->ikind == ICOUNT){
395 UInt_t scaldata = scalers[idx]->GetData(ichan);
396 dvars[ivar] = scaldata;
397 scal_present_read.push_back(scaldata);
398 scal_prev_read.push_back(0);
399 scal_overflows.push_back(0);
400 dvarsFirst[ivar] = 0.0;
401 }
402 if (scalerloc[ivar]->ikind == ITIME){
403 dvars[ivar] =fTotalTime;
404 dvarsFirst[ivar] = 0;
405 }
406 if (scalerloc[ivar]->ikind == IRATE) {
407 dvars[ivar] = (scalers[idx]->GetData(ichan))/fDeltaTime;
408 dvarsFirst[ivar] = dvars[ivar];
409 //printf("%s %f\n",scalerloc[ivar]->name.Data(),scalers[idx]->GetRate(ichan)); //checks
410 }
411 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE){
412 Int_t bcm_ind=-1;
413 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
414 {
415 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
416 if (match!=string::npos)
417 {
418 bcm_ind=itemp;
419 }
420 }
421 if (scalerloc[ivar]->ikind == ICURRENT) {
422 dvars[ivar]=0.;
423 if (bcm_ind != -1) {
424 dvars[ivar]=((scalers[idx]->GetData(ichan))/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
425 dvars[ivar]=dvars[ivar]+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(dvars[ivar]-fBCM_SatOffset[bcm_ind],0.0),2.0);
426
427 }
428 if (bcm_ind == fbcm_Current_Threshold_Index) scal_current= dvars[ivar];
429 }
430 if (scalerloc[ivar]->ikind == ICHARGE) {
431 if (bcm_ind != -1) {
432 Double_t cur_temp=((scalers[idx]->GetData(ichan))/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
433 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.0);
434 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
435 dvars[ivar]+=fBCM_delta_charge[bcm_ind];
436 }
437 }
438 // printf("1st event %i index %i fBCMname %s scalerloc %s offset %f gain %f computed %f\n",evcount, bcm_ind, fBCM_Name[bcm_ind],scalerloc[ivar]->name.Data(),fBCM_Offset[bcm_ind],fBCM_Gain[bcm_ind],dvars[ivar]);
439 }
440
441 if (fDebugFile) *fDebugFile << " dvarsFirst "<<scalerloc[ivar]->ikind<<" "<<dvarsFirst[ivar]<<endl;
442
443 } else { //ifnotusefirstevent
444 if (scalerloc[ivar]->ikind == ICOUNT) {
445 dvarsFirst[ivar] = scalers[idx]->GetData(ichan) ;
446 scal_present_read.push_back(dvarsFirst[ivar]);
447 scal_prev_read.push_back(0);
448 }
449 if (scalerloc[ivar]->ikind == ITIME){
450 dvarsFirst[ivar] = fTotalTime;
451 }
452 if (scalerloc[ivar]->ikind == IRATE) {
453 dvarsFirst[ivar] = (scalers[idx]->GetData(ichan))/fDeltaTime;
454 //printf("%s %f\n",scalerloc[ivar]->name.Data(),scalers[idx]->GetRate(ichan)); //checks
455 }
456 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE)
457 {
458 Int_t bcm_ind=-1;
459 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
460 {
461 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
462 if (match!=string::npos)
463 {
464 bcm_ind=itemp;
465 }
466 }
467 if (scalerloc[ivar]->ikind == ICURRENT) {
468 dvarsFirst[ivar]=0.0;
469 if (bcm_ind != -1) {
470 dvarsFirst[ivar]=((scalers[idx]->GetData(ichan))/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
471 dvarsFirst[ivar]=dvarsFirst[ivar]+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(dvars[ivar]-fBCM_SatOffset[bcm_ind],0.0),2.);
472 }
473 if (bcm_ind == fbcm_Current_Threshold_Index) scal_current= dvarsFirst[ivar];
474 }
475 if (scalerloc[ivar]->ikind == ICHARGE) {
476 if (bcm_ind != -1) {
477 Double_t cur_temp=((scalers[idx]->GetData(ichan))/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
478 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.);
479 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
480 dvarsFirst[ivar]+=fBCM_delta_charge[bcm_ind];
481 }
482 }
483 }
484 if (fDebugFile) *fDebugFile << " dvarsFirst "<<scalerloc[ivar]->ikind<<" "<<dvarsFirst[ivar]<<endl;
485 }
486 }
487 else {
488 cout << "THcScalerEvtHandler:: ERROR:: incorrect index "<<ivar<<" "<<idx<<" "<<ichan<<endl;
489 }
490 }else{ // evcount != 0
491 if (fDebugFile) *fDebugFile << "Debug dvars "<<i<<" "<<ivar<<" "<<idx<<" "<<ichan<<endl;
492 if ((ivar < scalerloc.size()) &&
493 (idx < scalers.size()) &&
494 (ichan < MAXCHAN)) {
495 if (scalerloc[ivar]->ikind == ICOUNT) {
496 UInt_t scaldata = scalers[idx]->GetData(ichan);
497 if(scaldata < scal_prev_read[nscal]) {
498 scal_overflows[nscal]++;
499 }
500 dvars[ivar] = scaldata + (1+((Double_t)kMaxUInt))*scal_overflows[nscal]
501 -dvarsFirst[ivar];
502 scal_present_read[nscal]=scaldata;
503 nscal++;
504 }
505 if (scalerloc[ivar]->ikind == ITIME) {
506 dvars[ivar] = fTotalTime;
507 }
508 if (scalerloc[ivar]->ikind == IRATE) {
509 UInt_t scaldata = scalers[idx]->GetData(ichan);
510 UInt_t diff;
511 if(scaldata < scal_prev_read[nscal-1]) {
512 diff = (kMaxUInt-(scal_prev_read[nscal-1] - 1)) + scaldata;
513 } else {
514 diff = scaldata - scal_prev_read[nscal-1];
515 }
516 dvars[ivar] = diff/fDeltaTime;
517 // printf("%s %f\n",scalerloc[ivar]->name.Data(),scalers[idx]->GetRate(ichan));//checks
518 }
519 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE)
520 {
521 Int_t bcm_ind=-1;
522 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
523 {
524 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
525 if (match!=string::npos)
526 {
527 bcm_ind=itemp;
528 }
529 }
530 if (scalerloc[ivar]->ikind == ICURRENT) {
531 dvars[ivar]=0;
532 if (bcm_ind != -1) {
533 UInt_t scaldata = scalers[idx]->GetData(ichan);
534 UInt_t diff;
535 if(scaldata < scal_prev_read[nscal-1]) {
536 diff = (kMaxUInt-(scal_prev_read[nscal-1] - 1)) + scaldata;
537 } else {
538 diff = scaldata - scal_prev_read[nscal-1];
539 }
540 dvars[ivar]=0.;
541 if (fDeltaTime>0) {
542 Double_t cur_temp=(diff/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
543 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.);
544
545 dvars[ivar]=cur_temp;
546 }
547 }
548 if (bcm_ind == fbcm_Current_Threshold_Index) scal_current= dvars[ivar];
549 }
550 if (scalerloc[ivar]->ikind == ICHARGE) {
551 if (bcm_ind != -1) {
552 UInt_t scaldata = scalers[idx]->GetData(ichan);
553 UInt_t diff;
554 if(scaldata < scal_prev_read[nscal-1]) {
555 diff = (kMaxUInt-(scal_prev_read[nscal-1] - 1)) + scaldata;
556 } else {
557 diff = scaldata - scal_prev_read[nscal-1];
558 }
559 fBCM_delta_charge[bcm_ind]=0;
560 if (fDeltaTime>0) {
561 Double_t cur_temp=(diff/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
562 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.);
563 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
564 }
565 dvars[ivar]+=fBCM_delta_charge[bcm_ind];
566 }
567 }
568 // printf("event %i index %i fBCMname %s scalerloc %s offset %f gain %f computed %f\n",evcount, bcm_ind, fBCM_Name[bcm_ind],scalerloc[ivar]->name.Data(),fBCM_Offset[bcm_ind],fBCM_Gain[bcm_ind],dvars[ivar]);
569 }
570 if (fDebugFile) *fDebugFile << " dvars "<<scalerloc[ivar]->ikind<<" "<<dvars[ivar]<<endl;
571 } else {
572 cout << "THcScalerEvtHandler:: ERROR:: incorrect index "<<ivar<<" "<<idx<<" "<<ichan<<endl;
573 }
574 }
575
576 }
577 //
578 for( auto& i: scalerloc ) {
579 size_t ivar = i->ivar;
580 size_t idx = i->index;
581 size_t ichan = i->ichan;
582 if (scalerloc[ivar]->ikind == ICUT+ICOUNT){
583 UInt_t scaldata = scalers[idx]->GetData(ichan);
584 if ( scal_current > fbcm_Current_Threshold) {
585 UInt_t diff;
586 if(scaldata < dvars_prev_read[ivar]) {
587 diff = (kMaxUInt-(dvars_prev_read[ivar] - 1)) + scaldata;
588 } else {
589 diff = scaldata - dvars_prev_read[ivar];
590 }
591 dvars[ivar] += diff;
592 }
593 dvars_prev_read[ivar] = scaldata;
594 }
595 if (scalerloc[ivar]->ikind == ICUT+ICHARGE){
596 Int_t bcm_ind=-1;
597 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
598 {
599 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
600 if (match!=string::npos)
601 {
602 bcm_ind=itemp;
603 }
604 }
605 if ( scal_current > fbcm_Current_Threshold && bcm_ind != -1) {
606 dvars[ivar] += fBCM_delta_charge[bcm_ind];
607 }
608 }
609 if (scalerloc[ivar]->ikind == ICUT+ITIME){
610 if ( scal_current > fbcm_Current_Threshold) {
611 dvars[ivar] += fDeltaTime;
612 }
613 }
614 }
615 //
616 evcount = evcount + 1;
618 //
619 for (size_t j=0; j<scal_prev_read.size(); j++) scal_prev_read[j]=scal_present_read[j];
620 //
621 for( auto& scaler: scalers ) scaler->Clear("");
622
623 return 1;
624}
625
626
628{
629 //
630 ReadDatabase(date);
631 const int LEN = 200;
632 char cbuf[LEN];
633
634 fStatus = kOK;
635 fNormIdx = -1;
636
637 for( auto& fDelayedEvent: fDelayedEvents )
638 delete [] fDelayedEvent;
639 fDelayedEvents.clear();
640
641 cout << "Howdy ! We are initializing THcScalerEvtHandler !! name = "
642 << fName << endl;
643
644 if(eventtypes.empty()) {
645 eventtypes.push_back(0); // Default Event Type
646 }
647
648 TString dfile;
649 dfile = fName + "scaler.txt";
650
651// Parse the map file which defines what scalers exist and the global variables.
652
653 TString sname0 = "Scalevt";
654 TString sname;
655 sname = fName+sname0;
656
657 FILE *fi = Podd::OpenDBFile(sname.Data(), date);
658 if ( !fi ) {
659 cout << "Cannot find db file for "<<fName<<" scaler event handler"<<endl;
660 return kFileError;
661 }
662
663 size_t minus1 = -1;
664 size_t pos1;
665 string scomment = "#";
666 string svariable = "variable";
667 string smap = "map";
668 vector<string> dbline;
669
670 while( fgets(cbuf, LEN, fi) ) {
671 std::string sin(cbuf);
672 std::string sinput(sin.substr(0,sin.find_first_of('#')));
673 if (fDebugFile) *fDebugFile << "string input "<<sinput<<endl;
674 dbline = Podd::vsplit(sinput);
675 if (!dbline.empty()) {
676 pos1 = FindNoCase(dbline[0],scomment);
677 if (pos1 != minus1) continue;
678 pos1 = FindNoCase(dbline[0],svariable);
679 if (pos1 != minus1 && dbline.size()>4) {
680 string sdesc;
681 for (size_t j=5; j<dbline.size(); j++) sdesc = sdesc+" "+dbline[j];
682 UInt_t islot = atoi(dbline[1].c_str());
683 UInt_t ichan = atoi(dbline[2].c_str());
684 UInt_t ikind = atoi(dbline[3].c_str());
685 if (fDebugFile)
686 *fDebugFile << "add var "<<dbline[1]<<" desc = "<<sdesc<<" islot= "<<islot<<" "<<ichan<<" "<<ikind<<endl;
687 TString tsname(dbline[4].c_str());
688 TString tsdesc(sdesc.c_str());
689 AddVars(tsname,tsdesc,islot,ichan,ikind);
690 // add extra scaler which is cut on the current
691 if (ikind == ICOUNT ||ikind == ITIME ||ikind == ICHARGE ) {
692 tsname=tsname+"Cut";
693 AddVars(tsname,tsdesc,islot,ichan,ICUT+ikind);
694 }
695 }
696 pos1 = FindNoCase(dbline[0],smap);
697 if (fDebugFile) *fDebugFile << "map ? "<<dbline[0]<<" "<<smap<<" "<<pos1<<" "<<dbline.size()<<endl;
698 if (pos1 != minus1 && dbline.size()>6) {
699 Int_t imodel, icrate, islot, inorm;
700 UInt_t header, mask;
701 char cdum[20];
702 sscanf(sinput.c_str(),"%s %d %d %d %x %x %d \n",cdum,&imodel,&icrate,&islot, &header, &mask, &inorm);
703 if ((fNormSlot >= 0) && (fNormSlot != inorm)) cout << "THcScalerEvtHandler::WARN: contradictory norm slot "<<fNormSlot<<" "<<inorm<<endl;
704 fNormSlot = inorm; // slot number used for normalization. This variable is not used but is checked.
705 Int_t clkchan = -1;
706 Double_t clkfreq = 1;
707 if (dbline.size()>8) {
708 clkchan = atoi(dbline[7].c_str());
709 clkfreq = 1.0*atoi(dbline[8].c_str());
710 fClockChan=clkchan;
711 fClockFreq=clkfreq;
712 }
713 if (fDebugFile) {
714 *fDebugFile << "map line "<<dec<<imodel<<" "<<icrate<<" "<<islot<<endl;
715 *fDebugFile <<" header 0x"<<hex<<header<<" 0x"<<mask<<dec<<" "<<inorm<<" "<<clkchan<<" "<<clkfreq<<endl;
716 }
717 switch (imodel) {
718 case 560:
719 scalers.push_back(new Scaler560(icrate, islot));
720 if(!fOnlyBanks) fRocSet.insert(icrate);
721 fModuleSet.insert(imodel);
722 break;
723 case 1151:
724 scalers.push_back(new Scaler1151(icrate, islot));
725 if(!fOnlyBanks) fRocSet.insert(icrate);
726 fModuleSet.insert(imodel);
727 break;
728 case 3800:
729 scalers.push_back(new Scaler3800(icrate, islot));
730 if(!fOnlyBanks) fRocSet.insert(icrate);
731 fModuleSet.insert(imodel);
732 break;
733 case 3801:
734 scalers.push_back(new Scaler3801(icrate, islot));
735 if(!fOnlyBanks) fRocSet.insert(icrate);
736 fModuleSet.insert(imodel);
737 break;
738 case 9001: // TI Scalers
739 scalers.push_back(new Scaler9001(icrate, islot));
740 if(!fOnlyBanks) fRocSet.insert(icrate);
741 fModuleSet.insert(imodel);
742 break;
743 case 9250: // FADC250 Scalers
744 scalers.push_back(new Scaler9250(icrate, islot));
745 if(!fOnlyBanks) fRocSet.insert(icrate);
746 fModuleSet.insert(imodel);
747 break;
748 default:
749 break;
750 }
751 if (!scalers.empty()) {
752 UInt_t idx = scalers.size()-1;
753 // Headers must be unique over whole event, not
754 // just within a ROC
755 scalers[idx]->SetHeader(header, mask);
756// The normalization slot has the clock in it, so we automatically recognize it.
757// fNormIdx is the index in scaler[] and
758// fNormSlot is the slot#, checked for consistency
759 if (clkchan >= 0) {
760 scalers[idx]->SetClock(defaultDT, clkchan, clkfreq);
761 cout << "Setting scaler clock ... channel = "<<clkchan<<" ... freq = "<<clkfreq<<endl;
762 if (fDebugFile) *fDebugFile <<"Setting scaler clock ... channel = "<<clkchan<<" ... freq = "<<clkfreq<<endl;
763 fNormIdx = idx;
764 if (islot != fNormSlot) cout << "THcScalerEvtHandler:: WARN: contradictory norm slot ! "<<islot<<endl;
765
766 }
767 }
768 }
769 }
770 }
771 // can't compare UInt_t to Int_t (compiler warning), so do this
772 nscalers=0;
773 for (size_t i=0; i<scalers.size(); i++) nscalers++;
774 // need to do LoadNormScaler after scalers created and if fNormIdx found
775 if (fDebugFile) *fDebugFile <<"fNormIdx = "<<fNormIdx<<endl;
776 if ((fNormIdx >= 0) && fNormIdx < nscalers) {
777 for (Int_t i = 0; i < nscalers; i++) {
778 if (i==fNormIdx) continue;
779 scalers[i]->LoadNormScaler(scalers[fNormIdx]);
780 }
781 }
782
783#ifdef HARDCODED
784 // This code is superseded by the parsing of a map file above. It's another way ...
785 if (fName == "Left") {
786 AddVars("TSbcmu1", "BCM x1 counts", 1, 4, ICOUNT);
787 AddVars("TSbcmu1r","BCM x1 rate", 1, 4, IRATE);
788 AddVars("TSbcmu3", "BCM u3 counts", 1, 5, ICOUNT);
789 AddVars("TSbcmu3r", "BCM u3 rate", 1, 5, IRATE);
790 } else {
791 AddVars("TSbcmu1", "BCM x1 counts", 0, 4, ICOUNT);
792 AddVars("TSbcmu1r","BCM x1 rate", 0, 4, IRATE);
793 AddVars("TSbcmu3", "BCM u3 counts", 0, 5, ICOUNT);
794 AddVars("TSbcmu3r", "BCM u3 rate", 0, 5, IRATE);
795 }
796#endif
797
798
799 DefVars();
800
801#ifdef HARDCODED
802 // This code is superseded by the parsing of a map file above. It's another way ...
803 if (fName == "Left") {
804 scalers.push_back(new Scaler1151(1,0));
805 scalers.push_back(new Scaler3800(1,1));
806 scalers.push_back(new Scaler3800(1,2));
807 scalers.push_back(new Scaler3800(1,3));
808 scalers[0]->SetHeader(0xabc00000, 0xffff0000);
809 scalers[1]->SetHeader(0xabc10000, 0xffff0000);
810 scalers[2]->SetHeader(0xabc20000, 0xffff0000);
811 scalers[3]->SetHeader(0xabc30000, 0xffff0000);
812 scalers[0]->LoadNormScaler(scalers[1]);
813 scalers[1]->SetClock(4, 7, 1024);
814 scalers[2]->LoadNormScaler(scalers[1]);
815 scalers[3]->LoadNormScaler(scalers[1]);
816 } else {
817 scalers.push_back(new Scaler3800(2,0));
818 scalers.push_back(new Scaler3800(2,0));
819 scalers.push_back(new Scaler1151(2,1));
820 scalers.push_back(new Scaler1151(2,2));
821 scalers[0]->SetHeader(0xceb00000, 0xffff0000);
822 scalers[1]->SetHeader(0xceb10000, 0xffff0000);
823 scalers[2]->SetHeader(0xceb20000, 0xffff0000);
824 scalers[3]->SetHeader(0xceb30000, 0xffff0000);
825 scalers[0]->SetClock(4, 7, 1024);
826 scalers[1]->LoadNormScaler(scalers[0]);
827 scalers[2]->LoadNormScaler(scalers[0]);
828 scalers[3]->LoadNormScaler(scalers[0]);
829 }
830#endif
831
832 // Verify that the slots are not defined twice
833 for (UInt_t i1=0; i1 < scalers.size()-1; i1++) {
834 for (UInt_t i2=i1+1; i2 < scalers.size(); i2++) {
835 if (scalers[i1]->GetSlot()==scalers[i2]->GetSlot())
836 cout << "THcScalerEvtHandler:: WARN: same slot defined twice"<<endl;
837 }
838 }
839 // Identify indices of scalers[] vector to variables.
840 for (UInt_t i=0; i < scalers.size(); i++) {
841 for( auto& j: scalerloc ) {
842 if (j->islot==static_cast<UInt_t>(scalers[i]->GetSlot()))
843 j->index = i;
844 }
845 }
846
847 if(fDebugFile) *fDebugFile << "THcScalerEvtHandler:: Name of scaler bank "<<fName<<endl;
848 for (size_t i=0; i<scalers.size(); i++) {
849 if(fDebugFile) {
850 *fDebugFile << "Scaler # "<<i<<endl;
851 scalers[i]->SetDebugFile(fDebugFile);
852 scalers[i]->DebugPrint(fDebugFile);
853 }
854 }
855
856
857 //
858 return kOK;
859}
860
861void THcScalerEvtHandler::AddVars(const TString& name, const TString& desc,
862 UInt_t islot, UInt_t ichan, UInt_t ikind )
863{
864 // need to add fName here to make it a unique variable. (Left vs Right HRS, for example)
865 TString name1 = fName + name;
866 TString desc1 = fName + desc;
867 // We don't yet know the correspondence between index of scalers[] and slots.
868 // Will put that in later.
869 scalerloc.push_back( new HCScalerLoc(name1, desc1, 0, islot, ichan, ikind,
870 scalerloc.size()) );
871}
872
874{
875 // called after AddVars has finished being called.
876 Nvars = scalerloc.size();
877 if (Nvars == 0) return;
878 dvars_prev_read = new UInt_t[Nvars]; // dvars_prev_read is a member of this class
879 dvars = new Double_t[Nvars]; // dvars is a member of this class
880 dvarsFirst = new Double_t[Nvars]; // dvarsFirst is a member of this class
881 memset(dvars, 0, Nvars*sizeof(Double_t));
882 memset(dvars_prev_read, 0, Nvars*sizeof(UInt_t));
883 memset(dvarsFirst, 0, Nvars*sizeof(Double_t));
884 if (gHaVars) {
885 if(fDebugFile) *fDebugFile << "THcScalerEVtHandler:: Have gHaVars "<<gHaVars<<endl;
886 } else {
887 cout << "No gHaVars ?! Well, that's a problem !!"<<endl;
888 return;
889 }
890 if(fDebugFile) *fDebugFile << "THcScalerEvtHandler:: scalerloc size "<<scalerloc.size()<<endl;
891 const Int_t* count = nullptr;
892 for (size_t i = 0; i < scalerloc.size(); i++) {
893 gHaVars->DefineByType(scalerloc[i]->name.Data(), scalerloc[i]->description.Data(),
894 &dvars[i], kDouble, count);
895 //gHaVars->DefineByType(scalerloc[i]->name.Data(), scalerloc[i]->description.Data(),
896 // &dvarsFirst[i], kDouble, count);
897 }
898}
899
int Int_t
unsigned int UInt_t
bool Bool_t
const UInt_t kMaxUInt
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t mask
R__EXTERN class THaVarList * gHaVars
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
static const UInt_t ICURRENT
static const UInt_t IRATE
static const UInt_t defaultDT
static const UInt_t ICOUNT
static const UInt_t ITIME
static const UInt_t MAXCHAN
static const UInt_t ICHARGE
static const UInt_t ICUT
UInt_t GetEvType() const
const UInt_t * GetRawDataBuffer() const
UInt_t GetEvNum() const
UInt_t GetEvLength() const
std::vector< UInt_t > eventtypes
virtual Bool_t IsMyEvent(UInt_t type) const
virtual void EvDump(THaEvData *evdata) const
std::ofstream * fDebugFile
virtual THaVar * DefineByType(const char *name, const char *desc, const void *loc, VarType type, const Int_t *count, const char *errloc="DefineByType")
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Event handler for Hall C scalers.
Int_t Analyze(THaEvData *evdata)
std::set< UInt_t > fModuleSet
std::vector< HCScalerLoc * > scalerloc
std::vector< UInt_t > scal_prev_read
std::vector< Decoder::GenScaler * > scalers
std::vector< UInt_t * > fDelayedEvents
virtual Int_t ReadDatabase(const TDatime &date)
std::vector< UInt_t > scal_present_read
std::vector< std::string > fBCM_Name
virtual void SetDelayedType(int evtype)
std::set< UInt_t > fRocSet
void AddVars(const TString &name, const TString &desc, UInt_t iscal, UInt_t ichan, UInt_t ikind)
virtual Int_t End(THaRunBase *r=nullptr)
Int_t AnalyzeBuffer(UInt_t *rdata, Bool_t onlysync)
std::vector< UInt_t > scal_overflows
THcScalerEvtHandler(const char *, const char *)
TString fName
static Bool_t Initialized()
const char * Data() const
virtual Int_t Fill()
virtual void SetAutoSave(Long64_t autos=-300000000)
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const override
virtual Int_t Branch(const char *folder, Int_t bufsize=32000, Int_t splitlevel=99)
RVec< PromoteType< T > > sin(const RVec< T > &v)
string::size_type FindNoCase(string data, string chunk)
Double_t Power(Double_t x, Double_t y)
Double_t Max(Double_t a, Double_t b)
STL namespace.