Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcHelicityScaler.cxx
Go to the documentation of this file.
1
19#include "THcHelicityScaler.h"
20//#include "GenScaler.h"
21#include "Scaler3801.h"
22#include "THaEvData.h"
23#include "THcGlobals.h"
24#include "THaGlobals.h"
25#include "THcParmList.h"
26#include "THcHelicity.h"
27#include "TMath.h"
28#include "TString.h"
29#include "TROOT.h"
30#include <cstring>
31#include <cstdio>
32#include <cstdlib>
33#include <iostream>
34#include <map>
35#include <bitset>
36#include <iterator>
37#include "THaVarList.h"
38#include "VarDef.h"
39#include "Helper.h"
40#include "Textvars.h" // for Podd::vsplit
41#include "THaString.h"
42
43using namespace std;
44using namespace Decoder;
46
47static const UInt_t ICOUNT = 1;
48static const UInt_t IRATE = 2;
49static const UInt_t ICURRENT = 3;
50static const UInt_t ICHARGE = 4;
51static const UInt_t ITIME = 5;
52static const UInt_t ICUT = 6;
53static const UInt_t MAXCHAN = 32;
54static const UInt_t defaultDT = 4;
55
56
57THcHelicityScaler::THcHelicityScaler(const char *name, const char* description)
58 : THaEvtTypeHandler(name,description),
59 fBankID(9801),
60 fUseFirstEvent(kTRUE), fDelayedType(-1),fHelicityCycleOffset(-3),
61 fBCM_Gain(0), fBCM_Offset(0), fBCM_SatOffset(0), fBCM_SatQuadratic(0), fBCM_delta_charge(0),
62 evcount(0), evcountR(0.0), ifound(0),
63 fNormIdx(-1), fNormSlot(-1),
64 dvars(0), dvarsFirst(0),
65 fScalerTree(0), fOnlyBanks(kFALSE),
66 fClockChan(-1), fLastClock(0)
67{
68
69 fRocSet.clear();
70 fModuleSet.clear();
71
72 //---------------------------------------------------------------------------
73
74 fROC=-1;
75 fNScalerChannels = 32;
76
77 AddEvtType(1);
78 AddEvtType(2);
79 AddEvtType(4);
80 AddEvtType(5);
81 AddEvtType(6);
82 AddEvtType(7);
83 AddEvtType(129);
84 SetDelayedType(129);
85
86}
87
89{
90
91 if (!TROOT::Initialized()) {
92 delete fScalerTree;
93 }
94 Podd::DeleteContainer(scalers);
95 Podd::DeleteContainer(scalerloc);
96 delete [] dvars;
97 delete [] dvarsFirst;
98 delete [] fBCM_SatOffset;
99 delete [] fBCM_SatQuadratic;
100 delete [] fBCM_delta_charge;
101 delete [] fBCM_Gain;
102 delete [] fBCM_Offset;
103
104 for( auto& fDelayedEvent: fDelayedEvents )
105 delete [] fDelayedEvent;
106 fDelayedEvents.clear();
107}
108
110{
111
112 // Process any delayed events in order received
113
114 for( auto rdata: fDelayedEvents ) {
115 evNumberR += 1;
116 AnalyzeBuffer(rdata);
117 }
118
119 for( auto& fDelayedEvent: fDelayedEvents )
120 delete [] fDelayedEvent;
121 fDelayedEvents.clear();
122
123 //Write the helicity variables to scaler tree
124 if (fScalerTree) { fScalerTree->Write(); }
125
126
127 //Compute Scaler Asymmetries
128 /*
129 for(Int_t i=0;i<fNScalerChannels;i++) {
130 if(fScalerSums[i]>0.5) {
131 fScalAsymmetry[i] = (fHScalers[0][i]-fHScalers[1][i])/fScalerSums[i];
132 fScalAsymmetryError[i] = 2*TMath::Sqrt(fHScalers[0][i]*fHScalers[1][i]
133 *fScalerSums[i])
134 /(fScalerSums[i]*fScalerSums[i]);
135 } else {
136 fScalAsymmetry[i] = 0.0;
137 fScalAsymmetryError[i] = 0.0;
138 }
139
140 }
141 */
142
143 //------Compute Time Asymmetries-------
144
145 //C.Y. 12/15/2020 Time Asymmetry / Error Calculations (with scaler_current cut)
146
147 if(fQuartetCount <= 1) {
148 fTimeAsymmetry = -1000.;
150 } else {
151 fTimeAsymmetry = fTimeAsymSum/fQuartetCount; //normalize asymmetry to total number of quartets
156 } else {
158 }
159 }
160
161
162
163 //------Compute Scaler Asymmetries-------
164
165 //C.Y. 12/15/2020 Scaler Asymmetry / Error Calculations (with scaler_current cut)
166
167 for(Int_t i=0; i<fNScalerChannels; i++) {
168
169 if(fQuartetCount <= 1) {
170 fScalAsymmetry[i] = -1000.;
171 fScalAsymmetryError[i] = 0.0;
172 } else {
173 fScalAsymmetry[i] = fScalAsymSum[i]/fQuartetCount; //normalize asymmetry to total number of quartets
178 } else {
179 fScalAsymmetryError[i] = 0.0;
180 }
181 }
182
183 }
184
185
186 //------Compute Charge Asymmetries-------
187
188 //C.Y. 12/14/2020 Charge Asymmetry / Error Calculations (with scaler_current cut)
189
190 //Set the helicity scaler module channels for each BCM
191 std::map<std::string, Int_t> bcmindex;
192 bcmindex["BCM1_Hel.scal"] = 0;
193 bcmindex["BCM2_Hel.scal"] = 2;
194 bcmindex["Unser_Hel.scal"] = 6;
195 bcmindex["BCM4A_Hel.scal"] = 10;
196 bcmindex["BCM4B_Hel.scal"] = 4;
197 bcmindex["BCM4C_Hel.scal"] = 12;
198 //bcmindex["1MHz_Hel.scal"] = 8;
199
200
201 //cout << endl << "---------------------- Beam Charge Asymmetries ---------------------- " << endl;
202 //cout << " BCM Total Charge Beam ON Beam ON Asymmetry" << endl;
203 //cout << " Name Charge Asymmetry Charge Asymmetry Error" << endl;
204
205 for(Int_t i=0;i<fNumBCMs;i++) {
206
207 if(fQuartetCount <= 1) {
208 fChargeAsymmetry[i] = -1000.;
209 fChargeAsymmetryError[i] = 0.0;
210 } else {
211 fChargeAsymmetry[i] = fChargeAsymSum[i]/fQuartetCount; //normalize charge asymmetry to total number of quartets (as the sum is for every quartet)
212
217 } else {
218 fChargeAsymmetryError[i] = 0.0;
219 }
220 }
221
222 //printf("%6s %12.2f %12.8f %12.2f %12.8f %12.8f\n",fBCM_Name[i].c_str(),fCharge[i],
223 // fChargeAsymmetry[i],fChargeSum[i],asy,asyerr);
224 }
225
226
227 //Compute +/- helicity Times (no BCM cut)
228 Double_t pclock = fHScalers[0][fClockChan];
229 Double_t mclock = fHScalers[1][fClockChan];
230
231 fTimePlus = pclock/fClockFreq;
232 fTimeMinus = mclock/fClockFreq;
233 //fTime = (pclock+mclock)/fClockFreq;
234 //if(pclock+mclock>0) {
235 //fTimeAsymmetry = (pclock-mclock)/(pclock+mclock);
236 //} else {
237 // fTimeAsymmetry = 0.0;
238 //}
239 //printf("TIME(s)%12.2f %12.8f %12.2f\n",fTime, fTimeAsymmetry, fTimeSum);
240
241
242 //Compute Helicity Trigger Asymmetries (no BCM cut)
245 } else {
246 fTriggerAsymmetry = 0.0;
247 }
248
249 return 0;
250}
251
252
254{
255
256 char prefix[2];
257 prefix[0]='g';
258 prefix[1]='\0';
259 fNumBCMs = 0;
260
261
262 DBRequest list[]={
263 {"NumBCMs",&fNumBCMs, kInt, 0, true},
264 {nullptr}
265 };
266 gHcParms->LoadParmValues((DBRequest*)&list, prefix);
267 if(fNumBCMs > 0) {
273
274 string bcm_namelist;
275 DBRequest list2[]={
276 {"BCM_Gain", fBCM_Gain, kDouble, (UInt_t) fNumBCMs},
277 {"BCM_Offset", fBCM_Offset, kDouble,(UInt_t) fNumBCMs},
278 {"BCM_SatQuadratic", fBCM_SatQuadratic, kDouble,(UInt_t) fNumBCMs,true},
279 {"BCM_SatOffset", fBCM_SatOffset, kDouble,(UInt_t) fNumBCMs,true},
280 {"BCM_Names", &bcm_namelist, kString},
281 {"BCM_Current_threshold", &fbcm_Current_Threshold, kDouble,0, true},
282 {"BCM_Current_threshold_index", &fbcm_Current_Threshold_Index, kInt,0,true},
283 {nullptr}
284 };
285
286
289 for(Int_t i=0;i<fNumBCMs;i++) {
290 fBCM_SatOffset[i]=0.;
291 fBCM_SatQuadratic[i]=0.;
292 }
293 gHcParms->LoadParmValues((DBRequest*)&list2, prefix);
294 vector<string> bcm_names = Podd::vsplit(bcm_namelist);
295 for(Int_t i=0;i<fNumBCMs;i++) {
296 fBCM_Name.push_back(bcm_names[i]+"_Hel.scal");
297 fBCM_delta_charge[i]=0.;
298 }
299 }
300
301
302 DBRequest list3[]={
303 {"helicity_CycleOffset",&fHelicityCycleOffset, kInt, 0, true},
304 {nullptr}
305 };
306 gHcParms->LoadParmValues((DBRequest*)&list3,"");
307
308 cout << "fHelicityCycleOffset= " <<fHelicityCycleOffset<<endl;
309
310
311
312
313 fTotalTime=0.;
315 fDeltaTime=-1.;
316
317
318
319 return kOK;
320}
321
322void THcHelicityScaler::SetDelayedType(int evtype) {
323
332 fDelayedType = evtype;
333}
334
336{
337
338 //C.Y. | THcScalerEvtHandler::Analyze() uses this flag (which is forced to be 1),
339 //but as to why, it is beyond me. For consistency, I have also used it here.
340 Int_t lfirst=1;
341
342 if(evdata->GetEvNum() > 0) {
343 evNumber=evdata->GetEvNum();
345 }
346
347 if ( !IsMyEvent(evdata->GetEvType()) ) return -1;
348
349 if (fDebugFile) {
350 *fDebugFile << endl << "---------------------------------- "<<endl<<endl;
351 *fDebugFile << "\nEnter THcHelicityScaler for fName = "<<fName<<endl;
352 EvDump(evdata);
353 }
354
355
356 if (lfirst && !fScalerTree) {
357
358 lfirst = 0;
359
360 //Assign a name to the helicity scaler tree
361 TString sname1 = "TSHel";
362 TString sname2 = sname1 + fName;
363 TString sname3 = fName + " Scaler Data";
364
365 if (fDebugFile) {
366 *fDebugFile << "\nAnalyze 1st time for fName = "<<fName<<endl;
367 *fDebugFile << sname2 << " " <<sname3<<endl;
368 }
369
370 //Create Scaler Tree
371 fScalerTree = new TTree(sname2.Data(),sname3.Data());
372 fScalerTree->SetAutoSave(200000000);
373
374 TString name, tinfo;
375
376 name = "evcount";
377 tinfo = name + "/D";
378 fScalerTree->Branch(name.Data(), &evcountR, tinfo.Data(), 4000);
379
380 name = "evNumber";
381 tinfo = name + "/D";
382 fScalerTree->Branch(name.Data(), &evNumberR, tinfo.Data(), 4000);
383
384 //C.Y. 12/02/2020 Added actual helicity to be stored in scaler tree
385 name = "actualHelicity";
386 tinfo = name + "/D";
387 fScalerTree->Branch(name.Data(), &actualHelicityR, tinfo.Data(), 4000);
388
389 //C.Y. 12/09/2020 Added quartetphase to be stored in scaler tree
390 name = "quartetPhase";
391 tinfo = name + "/D";
392 fScalerTree->Branch(name.Data(), &quartetphaseR, tinfo.Data(), 4000);
393
394 for (size_t i = 0; i < scalerloc.size(); i++) {
395 name = scalerloc[i]->name;
396 tinfo = name + "/D";
397 fScalerTree->Branch(name.Data(), &dvars[i], tinfo.Data(), 4000);
398 }
399
400 }
401
402 auto *rdata = (UInt_t*) evdata->GetRawDataBuffer();
403
404 if((Int_t)evdata->GetEvType() == fDelayedType) { // Save this event for processing later
405 UInt_t evlen = evdata->GetEvLength();
406
407 auto *datacopy = new UInt_t[evlen];
408 fDelayedEvents.push_back(datacopy);
409 memcpy(datacopy,rdata,evlen*sizeof(UInt_t));
410 return 1;
411 }
412
413 else { // A normal event
414 if (fDebugFile) *fDebugFile<<"\n\nTHcHelicityScaler :: Debugging event type "<<dec<<evdata->GetEvType()<< " event num = " << evdata->GetEvNum() << endl<<endl;
415 Int_t ret;
416 if((ret=AnalyzeBuffer(rdata)))
417 {
418 if (fDebugFile) *fDebugFile << "scaler tree ptr 1 "<<fScalerTree<<endl;
419 if (fDebugFile) *fDebugFile << "ret = "<< ret <<endl;
420 }
421
422 return ret;
423 }
424
425}
427{
428
429 fNTrigsInBuf = 0;
430
431 // Parse the data, load local data arrays.
432 auto *p = (UInt_t*) rdata;
433
434 UInt_t *plast = p+*p; // Index to last word in the bank
435 Int_t roc = -1;
436 UInt_t evlen = *p+1;
437
438 Int_t ifound=0;
439
440 while(p<plast) {
441 UInt_t banklen = *p;
442 p++; // point to header
443
444 if (fDebugFile) {
445 *fDebugFile << "Bank: " << hex << *p << dec << " len: " << *(p-1) << endl;
446 }
447 if((*p & 0xff00) == 0x1000) { // Bank Containing banks
448 if(evlen-*(p-1) > 1) { // Don't use overall event header
449 roc = (*p>>16) & 0xf;
450 if(fDebugFile) *fDebugFile << "ROC: " << roc << " " << evlen << " " << *(p-1) << hex << " " << *p << dec << endl;
451 // cout << "ROC: " << roc << " " << evlen << " " << *(p-1) << hex << " " << *p << dec << endl;
452 if(roc != fROC) { // Not a ROC with helicity scaler
453 p+=*(p-1)-1; // Skip to next ROC
454 }
455 }
456 p++; // Now pointing to a bank in the bank
457 } else if (((*p & 0xff00) == 0x100) && (*p != 0xC0000100)) {
458 // Bank containing integers. Look for scalers
459 // This is either ROC bank containing integers or
460 // a bank within a ROC containing data from modules of a single type
461 // Look for banks with the helicity scaler bank ID (9801)
462 // Assume that very first word is a scaler header
463 // At any point in the bank where the word is not a matching
464 // header, we stop.
465 UInt_t tag = (*p>>16) & 0xffff; // Bank ID (ROC #)
466 // UInt_t num = (*p) & 0xff;
467 UInt_t *pnext = p+banklen; // Next bank
468 p++; // First data word
469 // If the bank is not a helicity scaler bank
470 // or it is not one of the ROC containing helcity scaler data
471 // skip to the next bank
472 //cout << "BankID=" << tag << endl;
473
474 if (tag != fBankID) {
475 p = pnext; // Fall through to end of the above else if
476 // cout << " Skipping to next bank" << endl;
477
478 } else {
479 // This is a helicity scaler bank
480 if (roc == fROC) {
481 UInt_t nevents = (banklen-2)/fNScalerChannels;
482 //cout << "# of helicity events in bank:" << " " << nevents << endl;
483 if (nevents > 100) {
484 cout << "Error! Beam off for too long" << endl;
485 }
486
487 fNTrigsInBuf = 0;
488
489 // Save helcitiy and quad info for THcHelicity
490 for (UInt_t iev = 0; iev < nevents; iev++) { // find number of helicity events in each bank
492
493 //C.Y. 11/26/2020 This methods extracts the raw helicity information
494 AnalyzeHelicityScaler(p+index);
495 }
496 }
497 }
498
499 while(p < pnext) {
500 Int_t nskip = 0;
501 if(fDebugFile) {
502 *fDebugFile << "Scaler Header: " << hex << *p << dec;
503 }
504
505 if(nskip == 0) {
506 if(fDebugFile) {
507 *fDebugFile << endl;
508 }
509 break; // Didn't find a matching header
510 }
511 p = p + nskip;
512 }
513 p = pnext;
514 } else {
515 p = p+*(p-1); // Skip to next bank
516 }
517
518 }
519
520 if (fDebugFile) {
521 *fDebugFile << "Finished with decoding. "<<endl;
522 *fDebugFile << " Found flag = "<<ifound<<endl;
523 }
524
525 if (!ifound) return 0;
526
527 return 1;
528}
529
531{
532
533
534 Int_t hbits = (p[0]>>30) & 0x3; // quartet and helcity bits in scaler word
535 Bool_t isquartet = (hbits&2) != 0;
536 Int_t ispos = hbits&1;
537 Int_t actualhelicity = 0;
539 fNTrigsInBuf++;
540 fNTriggers++;
541
542 Int_t quartetphase = (fNTriggers-fFirstCycle)%4;
543
544
545 if(fFirstCycle >= -10) {
546 if(quartetphase == 0) {
547 Int_t predicted = RanBit30(fRingSeed_reported);
548 fRingSeed_reported = ((fRingSeed_reported<<1) | ispos) & 0x3FFFFFFF;
549 // Check if ringseed_predicted agrees with reported if(fNBits>=30)
550 if(fNBits >= 30 && predicted != fRingSeed_reported) {
551 cout << "THcHelicityScaler: Helicity Prediction Failed" << endl;
552 cout << "Reported " << bitset<32>(fRingSeed_reported) << endl;
553 cout << "Predicted " << bitset<32>(predicted) << endl;
554 }
555 fNBits++;
556 if(fNBits==30) {
557 cout << "THcHelicityScaler: A " << bitset<32>(fRingSeed_reported) <<
558 " found at cycle " << fNTriggers << endl;
559 }
560 //First H+ and qrt come in same window
561 if(!isquartet && fHelicityCycleOffset==0){
562 cout << "THcHelicityScaler: Quartet bit expected but not set (" <<
563 fNTriggers << ")" << endl;
564 fNBits = 0;
567 fFirstCycle = -100;
568 }
569 //Pre NPS/CODA 3 the H+ came 1 hel. window after qrt
570 }else if (quartetphase == 3 && fHelicityCycleOffset==-3){
571 if(!isquartet) {
572 cout << "THcHelicityScaler: Quartet bit expected but not set (" <<
573 fNTriggers << ")" << endl;
574 fNBits = 0;
577 fFirstCycle = -100;
578 }
579 }
580 } else { // First cycle not yet identified
581 if(isquartet) { // Helicity and quartet signal for next set of scalers
583 quartetphase = (fNTriggers-fFirstCycle)%4;
584 // Helicity at start of quartet is same as last of quartet, so we can start filling the seed
585 fRingSeed_reported = ((fRingSeed_reported<<1) | ispos) & 0x3FFFFFFF;
586 fNBits++;
587 if(fNBits==30) {
588 cout << "THcHelicityScaler: B " << bitset<32>(fRingSeed_reported) <<
589 " found at cycle " << fNTriggers << endl;
590 }
591 }
592 }
593
594 if(fNBits>=30) {
597
598#define DELAY9
599#ifdef DELAY9
600 if(quartetphase == 3) {
602 actualhelicity = (fRingSeed_actual&1)?+1:-1;
603 } else {
604 actualhelicity = (fRingSeed_actual&1)?+1:-1;
605 if(quartetphase == 0 || quartetphase == 1) {
606 actualhelicity = -actualhelicity;
607 }
608 }
609 quartetphase = (quartetphase+1)%4;
610#else
611 actualhelicity = (fRingSeed_actual&1)?+1:-1;
612 if(quartetphase == 1 || quartetphase == 2) {
613 actualhelicity = -actualhelicity;
614 }
615#endif
616 } else {
618 }
619
620 //C.Y. 12/09/2020 : Pass quartet phase to scaler tree varable
621 quartetphaseR = quartetphase;
622
623 //C.Y. 11/26/2020 The block of code below is used to extract the helicity information from each
624 //channel of the helicity scaler onto a variable to be stored in the scaler tree. For each channel,
625 //each helicity state (+, -, or undefined) is stored in a single varibale. Each helicity state
626 //is tagged to the corresponding scaler read via 'actualHelicityR' which is stored below..
627
628 //C.Y. Assign actual helicity value to a variable to be written to tree (the value may be (0, +hel (+1), or -hel(-1))
629 //Where actualhelicity=0 is NOT the MPS. It is zero when the actual helicity has not been determined.
630 actualHelicityR = actualhelicity;
631
632 //C.Y. 11/26/2020 Loop over all 32 scaler channels for a specific helicity scaler module (SIS 3801)
633 for(Int_t i=0;i<fNScalerChannels;i++) {
634
635 //C.Y. 11/26/2020 the count expression below gets the scaler raw helicity information (+, -, or MPS helicity states) for the ith channel
636 Int_t count = p[i]&0xFFFFFF; // Bottom 24 bits equivalent of scalers->Decode()
637 fScalerChan[i] = count; //pass the helicity raw information to each helicity scaler channel array element
638 }
639
640
641
642 //C.Y. 11/26/2020 The block of code below is used to get a cumulative sum of +/- helicity used
643 //for calculation of the cumulative beam charge asymmetry and other quantities in the ::End() method
644 if(actualhelicity!=0) {
645
646 //index=0 (+ helicity), index=1 (- helicity)
647 Int_t hindex = (actualhelicity>0)?0:1;
648
649 //C.Y. 11/24/2020 increment the counter for either '+' or '-' helicity states
650 (actualhelicity>0)?(fNTriggersPlus++):(fNTriggersMinus++);
651
652 //C.Y. 11/24/2020 Loop over all 32 scaler channels for a specific helicity scaler module (SIS 3801)
653 for(Int_t i=0;i<fNScalerChannels;i++) {
654
655 Int_t count = p[i]&0xFFFFFF; // Bottom 24 bits equivalent of scalers->Decode()
656
657 //Increment either the '+' (hindex=0) or '-' (hindex=1) helicity counts for each [i] scaler channel channel of a given module
658 fHScalers[hindex][i] += count;
659 fScalerSums[i] += count;
660 }
661 }
662
663
664 //Set the helicity scaler clock to define the time
665 fDeltaTime = fScalerChan[fClockChan]/fClockFreq; //total clock counts / clock_frequency (1MHz) for a specific scaler read interval
666 fTotalTime = fPrevTotalTime + fDeltaTime; //cumulative scaler time
667
668 if (fDeltaTime==0) {
669 cout << " ******************* Severe Warning ****************************" << endl;
670 cout << " In THcHelicityScaler have found fDeltaTime is zero !! " << endl;
671 cout << " ******************* Alert DAQ experts ***************************" << endl;
672 if (fDebugFile) *fDebugFile << " In THcHelicityScaler have found fDeltaTime is zero !! " << endl;
673 }
674
675 fPrevTotalTime=fTotalTime; //set the current total time to the previous time for the upcoming read
676
677 //C.Y. Nov 27, 2020 : (below) code to write the helicity raw data to a variable
678 //and map the variable to the scaler location
679
680 Double_t scal_current=0;
681 //Loop over each scaler variable from the map
682 for (size_t i = 0; i < scalerloc.size(); i++) {
683 size_t ivar = scalerloc[i]->ivar;
684 size_t idx = scalerloc[i]->index;
685 size_t ichan = scalerloc[i]->ichan;
686
687 //ANALYZE 1ST SCALER READ SEPARATE (There is no previous read before this one)
688 if (evcount==0) {
689
690 //Loop over the scaler variable (ivar), helicity slot (idx), and slot channel (ichan) --> idx=0 (only one helicity module per spectrometer arm)
691 if ((ivar < scalerloc.size()) &&
692 (idx < scalers.size()) &&
693 (ichan < MAXCHAN)) {
694
695 //If fUseFirstEvent=kTRUE (do not skip 1st read)
696 if(fUseFirstEvent){
697
698 //Write scaler counts
699 if (scalerloc[ivar]->ikind == ICOUNT){
700 dvars[ivar] = fScalerChan[ichan];
701 dvarsFirst[ivar] = 0.0;
702 }
703
704 //Write scaler time
705 if (scalerloc[ivar]->ikind == ITIME){
706 dvars[ivar] = fDeltaTime;
707 dvarsFirst[ivar] = 0.0;
708 }
709
710 //Write scaler rate
711 if (scalerloc[ivar]->ikind == IRATE) {
712 dvars[ivar] = (fScalerChan[ichan])/fDeltaTime;
713 dvarsFirst[ivar] = 0.0;
714 }
715
716 //Write either scaler current or scaler charge
717 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE){
718
719 //Get BCM index
720 Int_t bcm_ind=-1;
721 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
722 {
723 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
724 if (match!=string::npos)
725 {
726 bcm_ind=itemp;
727 }
728 }
729
730 //Calculate and write scaler current
731 if (scalerloc[ivar]->ikind == ICURRENT) {
732
733 //set default to zero
734 dvars[ivar]=0.;
735 //Check bcm index and calculate current in a temporary placeholder, "cur_temp"
736 if (bcm_ind != -1) {
737 Double_t cur_temp=((fScalerChan[ichan])/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
738 cur_temp = cur_temp + fBCM_SatQuadratic[bcm_ind]*TMath::Max(cur_temp-fBCM_SatOffset[i],0.0);
739 //Assign the calculated scaler current to dvars
740 dvars[ivar]=cur_temp;
741 }
742 //Check which bcm index to place a bcm current threshold cut later on. Assign the the beam current read to "scal_current" for later use.
743 if (bcm_ind == fbcm_Current_Threshold_Index) {scal_current= dvars[ivar];}
744 }
745
746 //Calculate andd write scaler charge
747 if (scalerloc[ivar]->ikind == ICHARGE) {
748 //Check bcm index and calculate current in a temporary placeholder, "cur_temp", to determine the charge later on.
749 if (bcm_ind != -1) {
750 Double_t cur_temp=((fScalerChan[ichan])/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
751 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.0);
752 //Calculate the charge for this scaler read
753 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
754 //Assign the charge to dvars
755 dvars[ivar] = fBCM_delta_charge[bcm_ind];
756 }
757 }
758 }
759 }
760
761 //If user has set fUseFirstEvent=kFALSE (then use dvarsFirst to store the information for the 1st event, and leave dvars at 0.0)
762 //(The same calculations as above, but assign dvars=0.0, so that it does not count when written to the scaler tree)
763 else { //ifnotusefirstevent
764
765 if (scalerloc[ivar]->ikind == ICOUNT) {
766 dvarsFirst[ivar] = fScalerChan[ichan];
767 dvars[ivar] = 0.0;
768 }
769
770 if (scalerloc[ivar]->ikind == ITIME){
771 dvarsFirst[ivar] = fDeltaTime;
772 dvars[ivar] = 0.0;
773 }
774
775 if (scalerloc[ivar]->ikind == IRATE) {
776 dvarsFirst[ivar] = (fScalerChan[ichan])/fDeltaTime;
777 dvars[ivar] = 0.0;
778 }
779
780 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE)
781 {
782 Int_t bcm_ind=-1;
783 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
784 {
785 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
786 if (match!=string::npos)
787 {
788 bcm_ind=itemp;
789 }
790 }
791
792 if (scalerloc[ivar]->ikind == ICURRENT) {
793
794 dvarsFirst[ivar] = 0.0;
795 dvars[ivar] = 0.0;
796
797 if (bcm_ind != -1) {
798 Double_t cur_temp=((fScalerChan[ichan])/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
799 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[i],0.0),2.0);
800 dvarsFirst[ivar] = cur_temp;
801 }
802
803 if (bcm_ind == fbcm_Current_Threshold_Index) scal_current= dvarsFirst[ivar];
804 }
805
806 if (scalerloc[ivar]->ikind == ICHARGE) {
807
808 if (bcm_ind != -1) {
809 Double_t cur_temp=((fScalerChan[ichan])/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
810 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.0);
811 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
812 dvarsFirst[ivar] = fBCM_delta_charge[bcm_ind];
813 }
814 }
815 }
816 }
817 }
818 else {
819 cout << "THcHelicityScaler:: ERROR:: incorrect index "<<ivar<<" "<<idx<<" "<<ichan<<endl;
820 }
821 }
822
823 //Calculate Scaler Quantities for ALL OTHER SCALER READS (OTHER THAN THE FIRST READ)
824 else{ // evcount != 0
825
826 if ((ivar < scalerloc.size()) &&
827 (idx < scalers.size()) &&
828 (ichan < MAXCHAN)) {
829
830 if (scalerloc[ivar]->ikind == ICOUNT) {
831 dvars[ivar] = fScalerChan[ichan];
832 }
833
834 if (scalerloc[ivar]->ikind == ITIME) {
835 dvars[ivar] = fDeltaTime;
836 }
837
838 if (scalerloc[ivar]->ikind == IRATE) {
839 dvars[ivar] = fScalerChan[ichan]/fDeltaTime;
840 }
841
842 if(scalerloc[ivar]->ikind == ICURRENT || scalerloc[ivar]->ikind == ICHARGE)
843 {
844 Int_t bcm_ind=-1;
845 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
846 {
847 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
848 if (match!=string::npos)
849 {
850 bcm_ind=itemp;
851 }
852 }
853 if (scalerloc[ivar]->ikind == ICURRENT) {
854 dvars[ivar]=0.;
855
856 if (bcm_ind != -1) {
857
858 if (fDeltaTime>0) {
859 Double_t cur_temp=(fScalerChan[ichan]/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
860 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.0);
861 dvars[ivar]=cur_temp;
862 }
863 }
864 if (bcm_ind == fbcm_Current_Threshold_Index) scal_current= dvars[ivar];
865 }
866
867 if (scalerloc[ivar]->ikind == ICHARGE) {
868
869 if (bcm_ind != -1) {
870 dvars[ivar] = 0.0;
871
872 if (fDeltaTime>0) {
873 Double_t cur_temp=(fScalerChan[ichan]/fDeltaTime-fBCM_Offset[bcm_ind])/fBCM_Gain[bcm_ind];
874 cur_temp=cur_temp+fBCM_SatQuadratic[bcm_ind]*TMath::Power(TMath::Max(cur_temp-fBCM_SatOffset[bcm_ind],0.0),2.0);
875 fBCM_delta_charge[bcm_ind]=fDeltaTime*cur_temp;
876 }
877 dvars[ivar] = fBCM_delta_charge[bcm_ind];
878 }
879 }
880 }
881
882 } else {
883 cout << "THcHelicityScaler:: ERROR:: incorrect index "<<ivar<<" "<<idx<<" "<<ichan<<endl;
884 }
885 }
886 }
887
888 //Analyze Scaler Reads ONLY FOR SCAL_CURRENTS > BCM_CURRENT THRESHOLD
889 //(These variables have the name structure: varibaleName_Cut)
890
891 for( auto& i: scalerloc ) {
892 size_t ivar = i->ivar;
893 size_t ichan = i->ichan;
894
895 if (scalerloc[ivar]->ikind == ICUT+ICOUNT){
896 if ( scal_current > fbcm_Current_Threshold) {
897 dvars[ivar] = fScalerChan[ichan];
898 }
899 }
900
901 if (scalerloc[ivar]->ikind == ICUT+ICHARGE){
902 Int_t bcm_ind=-1;
903 for(Int_t itemp =0; itemp<fNumBCMs;itemp++)
904 {
905 size_t match = string(scalerloc[ivar]->name.Data()).find(string(fBCM_Name[itemp]));
906 if (match!=string::npos)
907 {
908 bcm_ind=itemp;
909 }
910 }
911 if ( scal_current > fbcm_Current_Threshold && bcm_ind != -1) {
912 dvars[ivar] = fBCM_delta_charge[bcm_ind];
913 }
914 }
915
916 if (scalerloc[ivar]->ikind == ICUT+ITIME){
917 if ( scal_current > fbcm_Current_Threshold) {
918 dvars[ivar] = fDeltaTime;
919 }
920 }
921
922 }
923
924 //--------------------------------------
925
926 //C.Y. 12/13/2020 : Calculate the asymmetries / errors at the end of each quartet (S.A. Wood approach)
927
928 if(actualhelicity!=0) {
929
930 //Reset fHaveCycle to kFALSE at the start of each quartet (e.g. - + + -, + - - +
931 if(quartetphase==0) {
933 }
934
935 //Check if BCM scaler current is above set threshold
936 if(scal_current > fbcm_Current_Threshold && (quartetphase==0 || fHaveCycle[max(quartetphase-1,0)])) {
937 fHaveCycle[quartetphase] = kTRUE;
938
939 //Loop over each BCM to get the charge for each cycle of a quartet
940 for(Int_t i=0; i<fNumBCMs; i++) {
941 fChargeCycle[quartetphase][i] = fBCM_delta_charge[i];
942 }
943
944 //Loop over each Scaler Channel to the the counts for each cycle of a quartet
945 for(Int_t i=0; i<fNScalerChannels; i++) {
946 fScalCycle[quartetphase][i] = fScalerChan[i];
947 }
948
949 //Set the time for each cycle of the quartet
950 fTimeCycle[quartetphase] = fDeltaTime;
951
952 }
953
954 // Compute asymmetries at the end of each quartet
955 if(quartetphase == 3 && fHaveCycle[3]) {
956
957 //Loop over BCMs
958 for(Int_t i=0;i<fNumBCMs;i++) {
959
960 //compute charge asymmetry for each quartet at the end of said quartet
961 Double_t asy = actualhelicity*(fChargeCycle[0][i]+fChargeCycle[3][i]
962 - fChargeCycle[1][i]-fChargeCycle[2][i]) /
963 (fChargeCycle[0][i]+fChargeCycle[3][i]+fChargeCycle[1][i]+fChargeCycle[2][i]);
964
965 fChargeSum[i] += fChargeCycle[0][i]+fChargeCycle[1][i]
966 +fChargeCycle[2][i]+fChargeCycle[3][i];
967
968 //keep track of sums for proper error calculation
969 fChargeAsymSum[i] += asy;
970 fChargeAsymSum2[i] += asy*asy;
971 }
972
973 //-------
974
975 //Loop over Scaler Channels
976 for(Int_t i=0; i<fNScalerChannels; i++) {
977
978 //compute scaler asymmetry for each quartet at the end of said quartet
979 Double_t asy = actualhelicity*(fScalCycle[0][i]+fScalCycle[3][i]
980 - fScalCycle[1][i]-fScalCycle[2][i]) /
981 (fScalCycle[0][i]+fScalCycle[3][i]+fScalCycle[1][i]+fScalCycle[2][i]);
982
983 fScalSum[i] += fScalCycle[0][i]+fScalCycle[1][i]
984 +fScalCycle[2][i]+fScalCycle[3][i];
985
986 //keep track of sums for proper error calculation
987 fScalAsymSum[i] += asy;
988 fScalAsymSum2[i] += asy*asy;
989 }
990
991 //-------
992
993 //compute time asymmetry for each quartet at the end of said quartet
994 Double_t asy = actualhelicity*(fTimeCycle[0]+fTimeCycle[3]
995 - fTimeCycle[1]-fTimeCycle[2]) /
997
998 //keep track of total time
1000 +fTimeCycle[2]+fTimeCycle[3];
1001
1002 //keep track of sums for proper error calculation
1003 fTimeAsymSum += asy;
1004 fTimeAsymSum2 += asy*asy;
1005
1006 //------
1007
1008 //keep track of the total number of quartets
1009 fQuartetCount++;
1010
1011
1012 }
1013
1014
1015 }
1016
1017
1018 //--------------------------------------
1019
1020 //increment scaler reads
1021 evcount = evcount + 1;
1022 evcountR = evcount;
1023
1024 //clear Genscaler scalers
1025 for( auto& scaler: scalers ) scaler->Clear("");
1026
1027
1028 //C.Y. 12/02/2020 Fill Scaler Tree Here
1029 if (fScalerTree) {
1030 fScalerTree->Fill();
1031 }
1032
1033 return(0);
1034}
1035
1036//_____________________________________________________________________________
1038{
1039 UInt_t bit7 = (ranseed & 0x00000040) != 0;
1040 UInt_t bit28 = (ranseed & 0x08000000) != 0;
1041 UInt_t bit29 = (ranseed & 0x10000000) != 0;
1042 UInt_t bit30 = (ranseed & 0x20000000) != 0;
1043
1044 UInt_t newbit = (bit30 ^ bit29 ^ bit28 ^ bit7) & 0x1;
1045
1046 ranseed = ( (ranseed<<1) | newbit ) & 0x3FFFFFFF;
1047
1048 return ranseed;
1049
1050}
1051//_____________________________________________________________________________
1053{
1054
1055 ReadDatabase(date);
1056 const int LEN = 200;
1057 char cbuf[LEN];
1058
1059 fStatus = kOK;
1060 fNormIdx = -1;
1061
1062 for( auto& fDelayedEvent: fDelayedEvents )
1063 delete [] fDelayedEvent;
1064 fDelayedEvents.clear();
1065
1066 cout << "Howdy ! We are initializing THcHelicityScaler !! name = "
1067 << fName << endl;
1068
1069 if(eventtypes.empty()) {
1070 eventtypes.push_back(0); // Default Event Type
1071 }
1072
1073 if(fROC < 0) {
1074 fROC = 8; // Default to SHMS crate
1075 }
1076
1077
1078 //C.Y. 11/26/2020 : Read In and Parse the variables in the helicity scaler map file
1079 TString dfile;
1080 dfile = fName + "scaler.txt";
1081
1082 TString sname0 = "Scalevt";
1083 TString sname;
1084 sname = "hel"+fName+sname0; //This should be: 'helPScalevt' or 'helHScalevt'
1085
1086 //Open helicity scaler .dat map file
1087 FILE *fi = Podd::OpenDBFile(sname.Data(), date);
1088 if ( !fi ) {
1089 cout << "Cannot find db file for "<<fName<<" scaler event handler"<<endl;
1090 return kFileError;
1091 }
1092
1093 size_t minus1 = -1;
1094 size_t pos1;
1095 string scomment = "#";
1096 string svariable = "variable";
1097 string smap = "map";
1098 vector<string> dbline;
1099
1100 while( fgets(cbuf, LEN, fi) ) {
1101 std::string sin(cbuf);
1102 std::string sinput(sin.substr(0,sin.find_first_of('#')));
1103 if (fDebugFile) *fDebugFile << "string input "<<sinput<<endl;
1104 dbline = Podd::vsplit(sinput);
1105 if (!dbline.empty()) {
1106 pos1 = FindNoCase(dbline[0],scomment);
1107 if (pos1 != minus1) continue;
1108 pos1 = FindNoCase(dbline[0],svariable);
1109
1110 if (pos1 != minus1 && dbline.size()>4) {
1111 string sdesc;
1112 for (size_t j=5; j<dbline.size(); j++) sdesc = sdesc+" "+dbline[j];
1113 UInt_t islot = atoi(dbline[1].c_str());
1114 UInt_t ichan = atoi(dbline[2].c_str());
1115 UInt_t ikind = atoi(dbline[3].c_str());
1116 if (fDebugFile)
1117 *fDebugFile << "add var "<<dbline[1]<<" desc = "<<sdesc<<" islot= "<<islot<<" "<<ichan<<" "<<ikind<<endl;
1118 TString tsname(dbline[4].c_str());
1119 TString tsdesc(sdesc.c_str());
1120 AddVars(tsname,tsdesc,islot,ichan,ikind);
1121 // add extra scaler which is cut on the current
1122 if (ikind == ICOUNT ||ikind == ITIME ||ikind == ICHARGE ) {
1123 tsname=tsname+"Cut";
1124 AddVars(tsname,tsdesc,islot,ichan,ICUT+ikind);
1125 }
1126 }
1127
1128 pos1 = FindNoCase(dbline[0],smap);
1129 if (fDebugFile) *fDebugFile << "map ? "<<dbline[0]<<" "<<smap<<" "<<pos1<<" "<<dbline.size()<<endl;
1130 if (pos1 != minus1 && dbline.size()>6) {
1131 Int_t imodel, icrate, islot, inorm;
1132 UInt_t header, mask;
1133 char cdum[20];
1134 sscanf(sinput.c_str(),"%s %d %d %d %x %x %d \n",cdum,&imodel,&icrate,&islot, &header, &mask, &inorm);
1135 if ((fNormSlot >= 0) && (fNormSlot != inorm)) cout << "THcHelicityScaler::WARN: contradictory norm slot "<<fNormSlot<<" "<<inorm<<endl;
1136 fNormSlot = inorm; // slot number used for normalization. This variable is not used but is checked.
1137 Int_t clkchan = -1;
1138 Double_t clkfreq = 1;
1139 if (dbline.size()>8) {
1140 clkchan = atoi(dbline[7].c_str());
1141 clkfreq = 1.0*atoi(dbline[8].c_str());
1142 fClockChan=clkchan;
1143 fClockFreq=clkfreq;
1144 }
1145 if (fDebugFile) {
1146 *fDebugFile << "map line "<<dec<<imodel<<" "<<icrate<<" "<<islot<<endl;
1147 *fDebugFile <<" header 0x"<<hex<<header<<" 0x"<<mask<<dec<<" "<<inorm<<" "<<clkchan<<" "<<clkfreq<<endl;
1148 }
1149
1150 if (imodel==3801) {
1151 //switch(imodel) {
1152 //case 3801: //Hall C Helicity Scalers (SIS 3801)
1153 scalers.push_back(new Scaler3801(icrate, islot));
1154 if(!fOnlyBanks) fRocSet.insert(icrate);
1155 fModuleSet.insert(imodel);
1156 // break;
1157 }
1158
1159 if (!scalers.empty()) {
1160 UInt_t idx = scalers.size()-1;
1161 // Headers must be unique over whole event, not
1162 // just within a ROC
1163 scalers[idx]->SetHeader(header, mask);
1164 // The normalization slot has the clock in it, so we automatically recognize it.
1165 // fNormIdx is the index in scaler[] and
1166 // fNormSlot is the slot#, checked for consistency
1167 if (clkchan >= 0) {
1168 scalers[idx]->SetClock(defaultDT, clkchan, clkfreq);
1169 cout << "Setting scaler clock ... channel = "<<clkchan<<" ... freq = "<<clkfreq<<endl;
1170 if (fDebugFile) *fDebugFile <<"Setting scaler clock ... channel = "<<clkchan<<" ... freq = "<<clkfreq<<endl;
1171 fNormIdx = idx;
1172 if (islot != fNormSlot) cout << "THcHelicityScaler:: WARN: contradictory norm slot ! "<<islot<<endl;
1173
1174 }
1175 }
1176
1177 }
1178 }
1179
1180 } //end while loop
1181
1182 // can't compare UInt_t to Int_t (compiler warning), so do this
1183 nscalers=0;
1184 for (size_t i=0; i<scalers.size(); i++) nscalers++;
1185 // need to do LoadNormScaler after scalers created and if fNormIdx found
1186 if (fDebugFile) *fDebugFile <<"fNormIdx = "<<fNormIdx<<endl;
1187 if ((fNormIdx >= 0) && fNormIdx < nscalers) {
1188 for (Int_t i = 0; i < nscalers; i++) {
1189 if (i==fNormIdx) continue;
1190 scalers[i]->LoadNormScaler(scalers[fNormIdx]);
1191 }
1192 }
1193
1194 //Called after AddVars() has been called
1195 DefVars();
1196
1197 // Verify that the slots are not defined twice
1198 for (UInt_t i1=0; i1 < scalers.size()-1; i1++) {
1199 for (UInt_t i2=i1+1; i2 < scalers.size(); i2++) {
1200 if (scalers[i1]->GetSlot()==scalers[i2]->GetSlot())
1201 cout << "THcHelicityScaler:: WARN: same slot defined twice"<<endl;
1202 }
1203 }
1204 // Identify indices of scalers[] vector to variables.
1205 for (UInt_t i=0; i < scalers.size(); i++) {
1206 for( auto& j: scalerloc ) {
1207 if( j->islot == static_cast<UInt_t>(scalers[i]->GetSlot()) )
1208 j->index = i;
1209 }
1210 }
1211
1212 if(fDebugFile) *fDebugFile << "THcHelicityScaler:: Name of scaler bank "<<fName<<endl;
1213 for (size_t i=0; i<scalers.size(); i++) {
1214 if(fDebugFile) {
1215 *fDebugFile << "Scaler # "<<i<<endl;
1216 scalers[i]->SetDebugFile(fDebugFile);
1217 scalers[i]->DebugPrint(fDebugFile);
1218 }
1219 }
1220
1221
1222 //------Initialize Helicity Variables / Arrays-----
1223 fNTriggers = 0;
1224 fNTrigsInBuf = 0;
1225 fFirstCycle = -100;
1227 fRingSeed_actual = 0;
1228 fNBits = 0;
1230 fHScalers[0] = new Double_t[fNScalerChannels]; //+ helicity scaler counts
1231 fHScalers[1] = new Double_t[fNScalerChannels]; //- helicity scaler counts
1232 fScalerSums = new Double_t[fNScalerChannels]; //total helicity scaler counts (hel=0, excluded)
1233 fScalerChan = new Double_t[fNScalerChannels]; //C.Y. 11/26/2020 : added array to store helicity information per channel
1234
1235 for(Int_t i=0;i<fNScalerChannels;i++) {
1236 fHScalers[0][i] = 0.0;
1237 fHScalers[1][i] = 0.0;
1238 fScalerSums[i] = 0.0;
1239 fScalerChan[i] = 0.0;
1240 }
1241
1242
1243 fTimePlus = fTimeMinus = 0.0;
1244 fTriggerAsymmetry = 0.0;
1245
1246
1247
1248 //------C.Y.: 12/13/2020 Initialize Variables for quartet-by-quartet asymmetries (include BCM current threshold cut)---------
1249 //(and error calculation)
1250
1251 for(Int_t i=0; i<4; i++) {
1252
1253 fChargeCycle[i] = new Double_t[fNumBCMs];
1255
1256 fHaveCycle[i] = kFALSE;
1257
1258 for(Int_t j=0; j<fNumBCMs; j++) {
1259 fChargeCycle[i][j] = 0.0;
1260 }
1261
1262 for(Int_t j=0; j<fNScalerChannels; j++) {
1263 fScalCycle[i][j] = 0.0;
1264 }
1265
1266 fTimeCycle[i] = 0.0;
1267
1268 }
1269
1270 //Initialize quartet counter
1271 fQuartetCount = 0.0;
1272
1273 //Initialize variables for time asymmetry calculation
1274 fTimeSum = 0.0;
1275 fTimeAsymmetry = 0.0;
1276 fTimeAsymmetryError = 0.0;
1277 fTimeAsymSum = 0.0;
1278 fTimeAsymSum2 = 0.0;
1279
1280
1281 //Initialize variables for charge asymmetry calculation
1287
1288 for(Int_t i=0;i<fNumBCMs;i++) {
1289 fChargeSum[i] = 0.0;
1290 fChargeAsymmetry[i] = 0.0;
1291 fChargeAsymSum[i] = 0.0;
1292 fChargeAsymSum2[i] = 0.0;
1293 }
1294
1295 //Initialize variables for scaler asymmetry calculation
1301
1302 for(Int_t i=0;i<fNScalerChannels;i++) {
1303 fScalSum[i] = 0.0;
1304 fScalAsymmetry[i] = 0.0;
1305 fScalAsymSum[i] = 0.0;
1306 fScalAsymSum2[i] = 0.0;
1307 }
1308
1309
1310 //------------------------------------------------------------------------------------------
1311
1312
1313 //Call MakeParms() to define variables to be used in report file
1314 MakeParms();
1315
1316 return kOK;
1317}
1318
1320{
1321
1322 //Put Various helicity scaler results in gHcParms so they can be included in results.
1323
1324 //no bcm cut required
1325 gHcParms->Define(Form("g%s_hscaler_plus[%d]",fName.Data(),fNScalerChannels),
1326 "Plus Helcity Scalers",*(fHScalers[0]));
1327
1328 gHcParms->Define(Form("g%s_hscaler_minus[%d]",fName.Data(),fNScalerChannels),
1329 "Minus Helcity Scalers",*(fHScalers[1]));
1330
1331 //gHcParms->Define(Form("g%s_hscaler_sum[%d]",fName.Data(),fNScalerChannels),
1332 // "Helcity Scalers Sum",*fScalerSums);
1333
1334 gHcParms->Define(Form("g%s_hscaler_triggers",fName.Data()),
1335 "Total Helicity Scaler Triggers",fNTriggers);
1336
1337 gHcParms->Define(Form("g%s_hscaler_triggers_plus",fName.Data()),
1338 "Positive Helicity Scaler Triggers",fNTriggersPlus);
1339
1340 gHcParms->Define(Form("g%s_hscaler_triggers_minus",fName.Data()),
1341 "Negative Helicity Scaler Triggers",fNTriggersMinus);
1342
1343 gHcParms->Define(Form("g%s_hscaler_trigger_asy",fName.Data()),
1344 "Helicity Trigger Asymmetry",fTriggerAsymmetry);
1345
1346 //gHcParms->Define(Form("g%s_hscaler_time_plus",fName.Data()),
1347 // "Positive Helicity Time",fTimePlus);
1348
1349 //gHcParms->Define(Form("g%s_hscaler_time_minus",fName.Data()),
1350 // "Negative Helicity Time",fTimeMinus);
1351
1352 //bcm cut
1353 gHcParms->Define(Form("g%s_hscaler_sum[%d]",fName.Data(),fNScalerChannels),
1354 "Helcity Scalers Sum",*fScalSum);
1355
1356 gHcParms->Define(Form("g%s_hscaler_asy[%d]",fName.Data(),fNScalerChannels),
1357 "Helicity Scaler Asymmetry[%d]",*fScalAsymmetry);
1358
1359 gHcParms->Define(Form("g%s_hscaler_asyerr[%d]",fName.Data(),fNScalerChannels),
1360 "Helicity Scaler Asymmetry Error[%d]",*fScalAsymmetryError);
1361
1362 gHcParms->Define(Form("g%s_hscaler_charge[%d]",fName.Data(),fNumBCMs),
1363 "Helicity Gated Charge",*fChargeSum);
1364
1365 gHcParms->Define(Form("g%s_hscaler_charge_asy[%d]",fName.Data(),fNumBCMs),
1366 "Helicity Gated Charge Asymmetry",*fChargeAsymmetry);
1367
1368 gHcParms->Define(Form("g%s_hscaler_charge_asyerr[%d]",fName.Data(),fNumBCMs),
1369 "Helicity Gated Charge Asymmetry Error",*fChargeAsymmetryError);
1370
1371 gHcParms->Define(Form("g%s_hscaler_time",fName.Data()),
1372 "Helicity Gated Time (sec)",fTimeSum);
1373
1374 gHcParms->Define(Form("g%s_hscaler_time_asy",fName.Data()),
1375 "Helicity Gated Time Asymmetry",fTimeAsymmetry);
1376
1377 gHcParms->Define(Form("g%s_hscaler_time_asyerr",fName.Data()),
1378 "Helicity Gated Time Asymmetry Error",fTimeAsymmetryError);
1379
1380
1381
1382}
1383
1384
1385
1386//--------------------C.Y. Sep 20, 2020 : Added Utility Functions--------------------
1387
1388void THcHelicityScaler::AddVars( const TString& name, const TString& desc,
1389 UInt_t islot, UInt_t ichan, UInt_t ikind )
1390{
1391 if (fDebugFile) *fDebugFile << "C.Y. | Calling THcHelicityScaler::AddVars() "<<endl;
1392 // need to add fName here to make it a unique variable. (Left vs Right HRS, for example)
1393 TString name1 = fName + name;
1394 TString desc1 = fName + desc;
1395 // We don't yet know the correspondence between index of scalers[] and slots.
1396 // Will put that in later.
1397 scalerloc.push_back( new HCScalerLoc(name1, desc1, 0, islot, ichan, ikind,
1398 scalerloc.size()) );
1399}
1400
1402{
1403 if (fDebugFile) *fDebugFile << "C.Y. | Calling THcHelicityScaler::DefVars() "<<endl;
1404 // called after AddVars has finished being called.
1405 Nvars = scalerloc.size();
1406 if (Nvars == 0) return;
1407 dvars = new Double_t[Nvars]; // dvars is a member of this class
1408 dvarsFirst = new Double_t[Nvars]; // dvarsFirst is a member of this class
1409 memset(dvars, 0, Nvars*sizeof(Double_t));
1410 memset(dvarsFirst, 0, Nvars*sizeof(Double_t));
1411 if (gHaVars) {
1412 if(fDebugFile) *fDebugFile << "THcScalerEVtHandler:: Have gHaVars "<<gHaVars<<endl;
1413 } else {
1414 cout << "No gHaVars ?! Well, that's a problem !!"<<endl;
1415 return;
1416 }
1417 if(fDebugFile) *fDebugFile << "THcHelicityScaler:: scalerloc size "<<scalerloc.size()<<endl;
1418 const Int_t* count = nullptr;
1419 for (size_t i = 0; i < scalerloc.size(); i++) {
1420 gHaVars->DefineByType(scalerloc[i]->name.Data(), scalerloc[i]->description.Data(),
1421 &dvars[i], kDouble, count);
1422 }
1423}
1424
1425//---------------------------------------------------------------------------------
1426
1427
int Int_t
unsigned int UInt_t
bool Bool_t
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
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 ITIME
static const UInt_t ICHARGE
static const UInt_t ICUT
char * Form(const char *fmt,...)
virtual void Clear(Option_t *="")
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
THaVar * Define(const char *name, const Byte_t &var, const Int_t *count=nullptr)
virtual THaVar * DefineByType(const char *name, const char *desc, const void *loc, VarType type, const Int_t *count, const char *errloc="DefineByType")
virtual Int_t End(THaRunBase *r=0)
std::set< UInt_t > fRocSet
Int_t fHelicityHistory[200]
std::vector< HCScalerLoc * > scalerloc
Int_t AnalyzeBuffer(UInt_t *rdata)
std::set< UInt_t > fModuleSet
THcHelicityScaler(const char *, const char *)
void AddVars(const TString &name, const TString &desc, UInt_t iscal, UInt_t ichan, UInt_t ikind)
Double_t * fScalAsymmetryError
virtual Int_t ReadDatabase(const TDatime &date)
virtual void SetDelayedType(int evtype)
Int_t Analyze(THaEvData *evdata)
Double_t * fHScalers[2]
Double_t * fChargeAsymmetryError
Int_t AnalyzeHelicityScaler(const UInt_t *p)
Double_t * fScalCycle[4]
std::vector< Decoder::GenScaler * > scalers
std::vector< std::string > fBCM_Name
Double_t * fChargeCycle[4]
std::vector< UInt_t * > fDelayedEvents
static Int_t RanBit30(Int_t ranseed)
virtual ~THcHelicityScaler()
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
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)
double max(double x, double y)
string::size_type FindNoCase(string data, string chunk)
Double_t Power(Double_t x, Double_t y)
Double_t Sqrt(Double_t x)
Double_t Max(Double_t a, Double_t b)
STL namespace.