Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcHelicity.cxx
Go to the documentation of this file.
1
2//
3// THcHelicity
4//
5// Helicity of the beam in delayed mode using quartets
6// +1 = plus, -1 = minus, 0 = unknown
7//
8// Also supports in-time mode with delay = 0
9//
11
12#include "THcHelicity.h"
13
14#include "THaApparatus.h"
15#include "THaEvData.h"
16#include "THcGlobals.h"
17#include "THcParmList.h"
18#include "THcHelicityScaler.h"
19#include "TH1F.h"
20#include "TMath.h"
21#include <iostream>
22#include <bitset>
23
24using namespace std;
25
26//_____________________________________________________________________________
27THcHelicity::THcHelicity( const char* name, const char* description,
28 THaApparatus* app ):
29 THaHelicityDet( name, description, app ),
30 fnQrt(-1), fHelDelay(8), fMAXBIT(30)
31{
32 // for( Int_t i = 0; i < NHIST; ++i )
33 // fHisto[i] = 0;
34 // memset(fHbits, 0, sizeof(fHbits));
37}
38
39//_____________________________________________________________________________
41 : fnQrt(-1), fHelDelay(8), fMAXBIT(30)
42{
43 // Default constructor for ROOT I/O
44
45 // for( Int_t i = 0; i < NHIST; ++i )
46 // fHisto[i] = 0;
47}
48
49//_____________________________________________________________________________
51{
53
54 // for( Int_t i = 0; i < NHIST; ++i ) {
55 // delete fHisto[i];
56 // }
57}
58
59//_____________________________________________________________________________
61 // Call `Setup` before everything else.
63
67 fLastMPSTime = 0;
69
70 // Call initializer for base class.
71 // This also calls `ReadDatabase` and `DefineVariables`.
72 EStatus status = THaDetector::Init(date);
73 if (status) {
74 fStatus = status;
75 return fStatus;
76 }
77
78 fEvNumCheck = 0;
80
87 fNBits = 0;
88 lastispos = 0;
92 fNQRTProblems = 0;
93 fPeriodCheck = 0.0;
95 fCycle = 0.0;
96 fRecommendedFreq = -1.0;
97
98 fStatus = kOK;
99 return fStatus;
100}
101
102//_____________________________________________________________________________
103void THcHelicity::Setup(const char* name, const char* description) {
104 // Prefix for parameters in `param` file.
105 string kwPrefix = string(GetApparatus()->GetName()) + "_" + name;
106 std::transform(kwPrefix.begin(), kwPrefix.end(), kwPrefix.begin(), ::tolower);
107 fKwPrefix = kwPrefix;
108}
109
110//_____________________________________________________________________________
112{
113
114 cout << "In THcHelicity::ReadDatabase" << endl;
115 // Read general HelicityDet database values (e.g. fSign)
116 // Int_t st = THaHelicityDet::ReadDatabase( date );
117 // if( st != kOK )
118 // return st;
119
120 // Read readout parameters (ROC addresses etc.)
122 date, fQWEAKDebug );
123 if( st != kOK )
124 return st;
125
126 fSign = 1; // Default helicity sign
127 fRingSeed_reported_initial = 0; // Initial see that should predict reported
128 // helicity of first quartet.
129 fFirstCycle = -1; // First Cycle that starts a quad (0 to 3)
130 fNLastQuartet = -1;
131 fNQuartet = 0;
132 // fFreq = 29.5596;
133 fFreq = 120.0007547169;
134 fHelDelay=8;
135
136 DBRequest list[]= {
137 // {"_hsign", &fSign, kInt, 0, 1},
138 {"helicity_delay", &fHelDelay, kInt, 0, 1},
139 {"helicity_freq", &fFreq, kDouble, 0, 1},
140 // {"helicity_seed", &fRingSeed_reported_initial, kInt, 0, 1},
141 // {"helicity_cycle", &fFirstCycle, kInt, 0, 1},
142 {0}
143 };
144
145 gHcParms->LoadParmValues(list, "");
146
147 fMAXBIT=30;
148
149 fTIPeriod = 250000000.0/fFreq;
150
151 // maximum of event in the pattern, for now we are working with quartets
152 // Int_t localpattern[4]={1,-1,-1,1};
153 // careful, the first value here should always +1
154 // for(int i=0;i<fQWEAKNPattern;i++)
155 // {
156 // fPatternSequence.push_back(localpattern[i]);
157 // }
158 HWPIN=kTRUE;
159
160 fQuartet[0]=fQuartet[1]=fQuartet[2]=fQuartet[3]=0;
161
163 // Set the seed for predicted reported and predicted actual
164 } else {
165 // Initialize mode to find quartets and then seed
166 }
167
168 cout << "Helicity decoder initialized with frequency of " << fFreq
169 << " Hz and reporting delay of " << fHelDelay << " cycles." << endl;
170
171 return kOK;
172}
173
174//_____________________________________________________________________________
175
180//_____________________________________________________________________________
182{
183 // Initialize global variables
184
185 cout << "Called THcHelicity::DefineVariables with mode == "
186 << mode << endl;
187
188 if( mode == kDefine && fIsSetup ) return kOK;
189 fIsSetup = ( mode == kDefine );
190
191 // Define standard variables from base class
193
194 const RVarDef var[] = {
195 { "nqrt", "position of cycle in quartet", "fnQrt" },
196 { "hel", "actual helicity for event", "fActualHelicity" },
197 { "helrep", "reported helicity for event", "fReportedHelicity" },
198 { "helpred", "predicted reported helicity for event", "fPredictedHelicity" },
199 { "mps", "In MPS blanking period", "fMPS"},
200 { "pcheck", "Period check", "fPeriodCheck"},
201 { "cycle", "Helicity Cycle", "fCycle"},
202 { "qrt", "Last cycle of quartet", "fQrt"},
203 { 0 }
204 };
205 cout << "Calling THcHelicity DefineVarsFromList" << endl;
206 return DefineVarsFromList( var, mode );
207}
208//_____________________________________________________________________________
209
211{
212
213 cout<<" ++++++ THcHelicity::Print ++++++\n";
214
215 cout<<" +++++++++++++++++++++++++++++++++++++\n";
216
217 return;
218}
219
220//_____________________________________________________________________________
222{
224
225 // fHisto[0] = new TH1F("hel.seed","hel.seed",32,-1.5,30.5);
226 // fHisto[1] = new TH1F("hel.error.code","hel.error.code",35,-1.5,33.5);
227
228 return 0;
229}
230
231//_____________________________________________________________________________
232//void THcHelicity::FillHisto()
233//{
234// fHisto[0]->Fill(fRing_NSeed);
235// fHisto[1]->Fill(fErrorCode);
236// return;
237//}
238//_____________________________________________________________________________
240{
241 // used as a control for the helciity computation
242 // 2^0: if the reported number of events in a pattern is larger than fQWEAKNPattern
243 // 2^1: if the offset between the ring reported value and TIR value is not fOffsetTIRvsRing
244 // 2^2: if the reported time in the ring is 0
245 // 2^3: if the predicted reported helicity doesn't match the reported helicity in the ring
246 // 2^4: if the helicity cannot be computed using the SetHelicity routine
247 // 2^5: if seed is being gathered
248
249 if(fErrorCode==0)
250 fErrorCode=(1<<error);
251 // only one reported error at the time
252
253 return;
254}
255
256//_____________________________________________________________________________
258{
259 // Clear event-by-event data
260
263 fEvtype = 0;
264
265 fQrt=0;
266 fErrorCode=0;
267
268 return;
269}
270//_____________________________________________________________________________
276
277//_____________________________________________________________________________
279{
280
281 // Decode Helicity data.
282 // Return 1 if helicity was assigned, 0 if not, <0 if error.
283 evnum = evdata.GetEvNum();
284
285 Int_t err = ReadData( evdata ); // from THcHelicityReader class
286 if( err ) {
287 Error( Here("THcHelicity::Decode"), "Error decoding helicity data." );
288 return err;
289 }
290
292 fMPS = fIsMPS?1:0;
293 fQrt = fIsQrt?1:0; // Last of quartet
294
295#if 0
299 if(nhelev >0) {
300 for(Int_t i=0;i<nhelev;i++) {
301 fScaleQuartet = (fHelicityHistory[i] & 2)!=0;
302 Int_t ispos = fHelicityHistory[i]&1;
303 if(fScaleQuartet) {
304 fScalerSeed = ((fScalerSeed<<1) | ispos) & 0x3FFFFFFF;
305 if(fNBits >= fMAXBIT+0) {
306 Int_t seedscan = fScalerSeed;
307 Int_t nbehind;
308 for(nbehind=0;nbehind<4;nbehind++) {
309 if(seedscan == fRingSeed_reported) {
310 if(nbehind>1) {
311 cout << "Scaler seed behind " << nbehind
312 << " quartets" << endl;
313 cout << "Ev seed " << bitset<32>(fRingSeed_reported) <<endl;
314 cout << "Scaler Seed " << bitset<32>(fScalerSeed) << endl;
315 }
316 break;
317 }
318 seedscan = RanBit30(seedscan);
319 }
320 if(nbehind>4) {
321 cout << "Scaler seed does not match" << endl;
322 cout << "Ev seed " << bitset<32>(fRingSeed_reported) <<endl;
323 cout << "Scaler Seed " << bitset<32>(fScalerSeed) << endl;
324 }
325 }
326 }
327 }
328 }
329 }
330#endif
331
332 if(fHelDelay == 0) { // If no delay actual=reported (but zero if in MPS)
334 return 0;
335 }
336
337 if(fDisabled) {
339 return 0;
340 }
341
342 fEvNumCheck++;
343 Int_t evnum = evdata.GetEvNum();
344
345 // For split replay support
346 if(fEvNumCheck == 1) fEvNumCheck = evnum;
347
348 if(fEvNumCheck!=evnum) {
349 cout << "THcHelicity: Missed " << evnum-fEvNumCheck << " events at event " << evnum << endl;
350 cout << " Disabling helicity decoding for rest of run." << endl;
351 cout << " Make sure \"RawDecode_master in cuts file accepts all physics events." <<endl;
354 return 0;
355 }
356
357 fActualHelicity = -10.0;
358 if(fFirstEvProcessed) { // Normal processing
359 // cout << evnum << " " << fNCycle << " " << fIsMPS << " " << fFoundMPS << " " << fTITime << " "
360 // << fLastMPSTime << " " << fNBits << endl;
361 Int_t missed = 0;
362 // Double_t elapsed_time = (fTITime - fFirstEvTime)/250000000.0;
363 if(fIsMPS) {
364 fPeriodCheck = fmod(fTITime/(250000000.0/fFreq)-fPeriodCheckOffset,1.0);
365
369 if(fFoundMPS) {
370 if(fRecommendedFreq < 0.0) {
371 if(TMath::Abs(fPeriodCheck-0.5) > 0.25) {
373 }
374 }
376 if(missed < 1) { // was <=1
381 } else {
383 }
384 // If there is a skip, pass it off to next non MPS event
385 // Need to also check here for missed MPS's
386 // cout << "Found MPS" << endl;
387 // check for Nint((time-last)/period) > 1
388 } else {
392 }
393 } else if (fFoundMPS) { //
394 if(fTITime - fLastMPSTime > fTIPeriod) { // We missed MPS periods
396 if(missed > 1) {
397 // cout << "Missed " << missed << " MPSes" << endl;
398 Int_t newNCycle = fNCycle + missed -1; // How many cycles really missed
399 Int_t quartets_missed = (newNCycle-fFirstCycle)/4 - (fNCycle-fFirstCycle)/4;
400 int quartetphase = (newNCycle-fFirstCycle)%4;
401 // cout << " " << fNCycle << " " << newNCycle << " " << fFirstCycle << " " << quartets_missed << " " << quartetphase << endl;
402 // cout << "Cycles " << fNCycle << " " << newNCycle << " " << fFirstCycle
403 // << " skipped " << quartets_missed << " quartets" << endl;
404 fNCycle = newNCycle;
405 // Need to reset fQuartet to reflect where we are based on the current
406 // reported helicity. So we don't fail quartet testing.
407 // But only do this if we are calibrated.
408 if(fNBits >= fMAXBIT+0) {
409
410 for(Int_t i=0;i<quartets_missed;i++) { // Advance the seeds.
411 // cout << "Advancing seed A " << fNBits << endl;
414 }
415
418 fActualHelicity = (quartetphase==0||quartetphase==3)?
420 fPredictedHelicity = (quartetphase==0||quartetphase==3)?
422
423 if (((fNCycle - fFirstCycle)%2)==1) {
425 fQuartet[1] = fQuartet[2] = -fQuartet[0];
426 } else {
428 fQuartet[2] = -fQuartet[1];
429 }
430 } else {
433 fQuartet[1] = 0;
434 }
435 }
436 fLastMPSTime += missed*fTIPeriod;
439 } else { // No missed periods. Get helicities from rings
440 if(fNBits>=fMAXBIT+0) {
441 int quartetphase = (fNCycle-fFirstCycle)%4;
444 fActualHelicity = (quartetphase==0||quartetphase==3)?
446 fPredictedHelicity = (quartetphase==0||quartetphase==3)?
448 } else {
449 fActualHelicity = 0;
450 }
451 }
452 if(fIsNewCycle) {
453 // cout << "Scaler Seed " << bitset<32>(fScalerSeed) << endl;
454 //Int_t predictedScalerSeed = RanBit30(fScalerSeed);
455 // cout << "Predct Seed " << bitset<32>(predictedScalerSeed) << endl;
456 //cout << "Event Seed " << bitset<32>(fRingSeed_reported) << endl;
457
460 fNCycle++;
461 // cout << "XX: " << fNCycle << " " << fReportedHelicity << " " << fIsQrt <<endl;
462
463 if(fIsQrt) { //
466 if(fHaveQRT) { // Should already have phase set
467 if((fNCycle-fFirstCycle)%4 != 0) {// Test if first in a quartet
469 if(fNQRTProblems > 10) {
470 cout << "QRT Problem resetting" << endl;
473 fNQRTProblems = 0;
474 } else {
475 cout << "Ignored " << fNQRTProblems << " problems" << endl;
476 }
477 } else {
478 fNQRTProblems = 0;
479 }
480 }
481 if(!fHaveQRT) {
482 fHaveQRT = kTRUE;
484
485 fFirstCycle = fNCycle; //
487 fNLastQuartet = fNQuartet - 1; // Make sure LoadHelicity uses
488 fNBits = 0;
489 fNQRTProblems = 0;
490 cout << "Phase found from QRT signal" << endl;
491 cout << "fFirstcycle = " << fFirstCycle << endl;
492 }
493 } else {
494 if(fHaveQRT) { // Using qrt signal.
497 if((fNCycle-fFirstCycle)%4 == 0) { // Shouldn't happen
499 if(fNQRTProblems > 10) {
500 cout << "Shouldn't happen, cycle=" << fNCycle << "/" << fFirstCycle << endl;
501 fHaveQRT = kFALSE; // False until a new QRT seen
502 fNBits = 0; // Reset
503 fNLastQuartet = fNQuartet; // Make sure LoadHelicity does not use
504 } else {
505 cout << "Ignored " << fNQRTProblems << " problems" << endl;
506 }
507 }
508 } else { // Presumable pre qrt signal data
509 if((fNCycle-fFirstCycle)%4 == 3) {// Test if last in a quartet
510 if((abs(fQuartet[0]+fQuartet[3]-fQuartet[1]-fQuartet[2])==4)) {
511 if(!fFoundQuartet) {
512 // fFirstCycle = fNCycle - 3;
513 cout << "Quartet potentially found, starting at cycle " << fFirstCycle
514 << endl;
516 fNLastQuartet = fNQuartet - 1; // Make sure LoadHelicity uses
518 }
519 } else {
520 if(fNCycle - fFirstCycle > 4) { // Not at start of run. Reset
521 cout << "Lost quartet sync at cycle " << fNCycle << endl;
522 cout << fQuartet[0] << " " << fQuartet[1] << " " << fQuartet[2] << " "
523 << fQuartet[3] << endl;
524 fFirstCycle += 4*((fNCycle-fFirstCycle)/4); // Update, but don't change phase
525 }
527 fNBits = 0;
528 cout << "Searching for first of a quartet at cycle " << " " << fFirstCycle
529 << endl;
530 cout << fQuartet[0] << " " << fQuartet[1] << " " << fQuartet[2] << " "
531 << fQuartet[3] << endl;
532 fFirstCycle++;
533 }
534 } else {
537 }
538 }
539 }
540 // Load the actual helicity. Calibrate if not calibrated.
542 // Here if we know we missed some earlier in the quartet, we need
543 // to make sure we get here and call LoadHelicity for the missing
544 // Cycles, reducing the missed count for each extra loadhelicity
545 // call that we make.
549 // cout << fTITime/250000000.0 << " " << fNCycle << " " << fReportedHelicity << endl;
550 // cout << fNCycle << ": " << fReportedHelicity << " "
551 // << fPredictedHelicity << " " << fActualHelicity << endl;
552 }
553 // Ignore until a MPS Is found
554
555 } else { // No MPS found yet
557 }
558 } else {
559 cout << "Initializing" << endl;
565 fLastMPSTime = fTITime; // Not necessarily during the MPS
566 fNCycle = 0;
571 fNBits = 0;
572 }
573
574 // Some sanity checks
575 if(fActualHelicity < -5) {
576 cout << "Actual Helicity never got defined" << endl;
577 }
578 if(fNBits < fMAXBIT+0) {
579 if(fActualHelicity == -1 || fActualHelicity == 1) {
580 cout << "Helicity of " << fActualHelicity << " reported prematurely at cycle " << fNCycle << endl;
581 }
582 }
584 return 0;
585}
586//_____________________________________________________________________________
588{
589 // End of run processing. Write histograms.
591
592 // for( Int_t i = 0; i < NHIST; ++i )
593 // fHisto[i]->Write();
594
595 if(fRecommendedFreq < 0.0) {
597 }
598 if(TMath::Abs(1-fRecommendedFreq/fFreq) >= 0.5e-6) {
599 cout << "------------- HELICITY DECODING ----------------------" << endl;
600 cout << "Actual helicity reversal frequency differs from \"helicity_freq\" value" << endl;
601 cout << "If there are helicity decoding errors beyond the start of the run, " << endl;
602 streamsize ss = cout.precision();
603 cout.precision(10);
604 cout << "try replacing helicity_freq value of " << fFreq << " with " << fRecommendedFreq << endl;
605 cout << "If that still gives helicity errors, try " << 0.9999999*fRecommendedFreq << endl;
606 cout.precision(ss);
607 cout << "------------------------------------------------------" << endl;
608 }
609
610 return 0;
611}
612
613//_____________________________________________________________________________
615{
616 // Set debug level of this detector as well as the THcHelicityReader
617 // helper class.
618
620 fQWEAKDebug = level;
621}
622
623//_____________________________________________________________________________
624void THcHelicity::LoadHelicity(Int_t reportedhelicity, Int_t cyclecount, Int_t missedcycles)
625{
626 // static const char* const here = "THcHelicity::LoadHelicity";
627 int quartetphase = (cyclecount-fFirstCycle)%4;
628 fnQrt = quartetphase;
629
630 // if(!fFixFirstCycle) {
631 if(fNQuartet - fNLastQuartet > 1) { // If we missed a quartet
632 if(fNBits< fMAXBIT+0) { // and we haven't gotten the seed, start over
633 cout << "fNBits = 0, missedcycles=" << missedcycles <<
634 " " << fNLastQuartet << " " << fNQuartet << endl;
635 fNBits = 0;
636 return;
637 }
638 }
639 // }
640 if(!fFoundQuartet) { // Wait until we have found quad phase before starting
641 return; // to calibrate
642 }
643 // LoadHelicity is called first event of each cycle.
644 // But only do it's thing if it is first cycle of quartet.
645 // Need to instead have it do it's thing on the first event of a quartet.
646 // But only for seed building. Current logic is OK once seed is found.
647 if(fNBits < fMAXBIT+0) {
648 if(fNQuartet != fNLastQuartet) {
649 // Sanity check fNQuartet == fNLastQuartet+1
650 if(fNBits == 0) {
651 cout << "Start calibrating at cycle " << cyclecount << endl;
653 }
654 // Do phase stuff right here
655 if((fReportedHelicity == kPlus && (quartetphase==0 || quartetphase == 3))
656 || (fReportedHelicity == kMinus && (quartetphase==1 || quartetphase == 2))) {
657 fRingSeed_reported = ((fRingSeed_reported<<1) | 1) & 0x3FFFFFFF;
658 // cout << " " << fNQuartet << " 1" << endl;
659 } else {
660 fRingSeed_reported = (fRingSeed_reported<<1) & 0x3FFFFFFF;
661 // cout << " " << fNQuartet << " 0" << endl;
662 }
663 fNBits++;
665 fNBits = 0;
667 } else if (fNBits==fMAXBIT+0) {
668 cout << " Seed Found " << bitset<32>(fRingSeed_reported) << " at cycle " << cyclecount << " with first cycle " << fFirstCycle << endl;
670 cout << "Scaler Seed " << bitset<32>(fglHelicityScaler->GetReportedSeed()) << endl;
671 }
673 cout << "Seed at cycle " << fFirstCycle << " should be " << hex << backseed << dec << endl;
674 // Create the "actual seed"
676 for(Int_t i=0;i<fHelDelay/4; i++) {
678 }
681 } else if (fglHelicityScaler && fNBits>2) { // Try the scalers
684 fRingSeed_reported = RanBit30(scalerseed);
685 cout << " -- Getting seed from scalers -- " << endl;
686 cout << " Seed Found " << bitset<32>(fRingSeed_reported) << " at cycle " << cyclecount << " with first cycle " << fFirstCycle << endl;
687 cout << "Scaler Seed " << bitset<32>(scalerseed) << endl;
688 // Create the "actual seed"
690 for(Int_t i=0;i<fHelDelay/4; i++) {
692 }
695 fNBits = fMAXBIT+0;
696 }
697 }
699 } // Need to change this to build seed even when not at start of quartet
700 } else {
701 if(quartetphase == 0) {
702 // If quartetphase !=, the seeds will alread have been advanced
703 // except that we won't have made the initial
704 // cout << "Advancing seed B " << fNBits << endl;
707
710 // if(fTITime/250000000.0 > 380.0) cout << fTITime/250000000.0 << " " << fNCycle << " " << hex <<
711 // fRingSeed_reported << " " << fRingSeed_actual << dec << endl;
713 cout << "Helicity prediction failed " << fReportedHelicity << " "
714 << fPredictedHelicity << " " << fActualHelicity << endl;
715 cout << bitset<32>(fRingSeed_reported) << " " << bitset<32>(fRingSeed_actual) << endl;
716 fNBits = 0; // Need to reaquire seed
719 }
722 }
723 fActualHelicity = (quartetphase==0||quartetphase==3)?
725 fPredictedHelicity = (quartetphase==0||quartetphase==3)?
727 }
728 return;
729}
730//_____________________________________________________________________________
732{
733
734 UInt_t bit7 = (ranseed & 0x00000040) != 0;
735 UInt_t bit28 = (ranseed & 0x08000000) != 0;
736 UInt_t bit29 = (ranseed & 0x10000000) != 0;
737 UInt_t bit30 = (ranseed & 0x20000000) != 0;
738
739 UInt_t newbit = (bit30 ^ bit29 ^ bit28 ^ bit7) & 0x1;
740
741
742 if(ranseed<=0) {
743 if(fQWEAKDebug>1)
744 std::cerr<<"ranseed must be greater than zero!"<<"\n";
745
746 newbit = 0;
747 }
748 ranseed = ( (ranseed<<1) | newbit ) & 0x3FFFFFFF;
749 //here ranseed is changed
750 if(fQWEAKDebug>1)
751 {
752 cout<< "THcHelicity::RanBit30, newbit="<<newbit<<"\n";
753 }
754
755 return ranseed;
756
757}
758//_____________________________________________________________________________
760/* Back track the seed by 30 samples */
761{
762#if 1
763 Int_t seed = currentseed;
764 for(Int_t i=0;i<30;i++) {
765 UInt_t bit1 = (seed & 0x00000001) != 0;
766 UInt_t bit8 = (seed & 0x00000080) != 0;
767 UInt_t bit29 = (seed & 0x10000000) != 0;
768 UInt_t bit30 = (seed & 0x20000000) != 0;
769
770 UInt_t newbit30 = (bit30 ^ bit29 ^ bit8 ^ bit1) & 0x1;
771 seed = (seed >> 1) | (newbit30<<29);
772 }
773#else
774 Int_t bits = currentseed;
775 Int_t seed=0;
776 for(Int_t i=0;i<30;i++) {
777 Int_t val;
778 // XOR at virtual position 0 and 29
779 if(i==0) {
780 val = ((bits & (1<<(i)))!=0) ^ ((bits & (1<<(i+29)))!=0);
781 } else {
782 val = ((bits & (1<<(i)))!=0) ^ ((seed & (1<<(i-1)))!=0);
783 }
784 if(i<=1) {
785 val = ((bits & (1<<(1-i)))!=0) ^ val;
786 } else {
787 val = ((seed & (1<<(i-2)))!=0) ^ val;
788 }
789 if(i<=22) {
790 val = ((bits & (1<<(i-22)))!=0) ^ val;
791 } else {
792 val = ((seed & (1<<(i-23)))!=0) ^ val;
793 }
794 seed |= (val<<i);
795 }
796#endif
797 return seed;
798}
799
801
int Int_t
unsigned int UInt_t
#define f(i)
const Bool_t kFALSE
const Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char mode
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
virtual void SetDebug(Int_t level)
static Int_t DefineVarsFromList(const void *list, EType type, EMode mode, const char *def_prefix, const TObject *obj, const char *prefix, const char *here, const char *comment_subst="")
virtual const char * Here(const char *) const
const char * GetPrefix() const
virtual void MakePrefix()
THaApparatus * GetApparatus() const
UInt_t GetEvNum() const
virtual void Clear(Option_t *opt="")
virtual const char * GetDBFileName() const
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual void Clear(Option_t *opt="")
virtual Int_t ReadData(const THaEvData &evdata)
Int_t ReadDatabase(const char *dbfilename, const char *prefix, const TDatime &date, int debug_flag=0)
virtual Int_t GetReportedSeed()
virtual Bool_t IsSeedGood()
virtual Int_t * GetHelicityHistoryP()
virtual Int_t GetNcycles()
virtual Int_t GetNevents()
Long64_t fLastEvTime
Definition THcHelicity.h:64
THcHelicityScaler * fglHelicityScaler
Int_t fLastHelpCycle
virtual Int_t DefineVariables(EMode mode=kDefine)
Int_t fQuartet[4]
Definition THcHelicity.h:78
Int_t fScalerSeed
void SetErrorCode(Int_t error)
Bool_t fIsNewCycle
Definition THcHelicity.h:74
Int_t fNCycle
Definition THcHelicity.h:75
Bool_t fFoundQuartet
Definition THcHelicity.h:73
virtual Int_t Decode(const THaEvData &evdata)
Int_t fRingSeed_reported_initial
Definition THcHelicity.h:49
std::string fKwPrefix
Definition THcHelicity.h:41
Int_t fLastReportedHelicity
Definition THcHelicity.h:62
Int_t fLastActualHelicity
Int_t fRingSeed_actual
Definition THcHelicity.h:85
Int_t fReportedHelicity
Definition THcHelicity.h:66
Int_t lastispos
virtual void SetDebug(Int_t level)
Bool_t fFixFirstCycle
Definition THcHelicity.h:51
Int_t fQuadPattern[8]
virtual Int_t ReadDatabase(const TDatime &date)
Int_t fActualHelicity
Definition THcHelicity.h:69
Int_t fScaleQuartet
Int_t fQuartetStartHelicity
Definition THcHelicity.h:70
virtual void MakePrefix()
Int_t fMAXBIT
Definition THcHelicity.h:91
void LoadHelicity(Int_t reportedhelicity, Int_t cyclecount, Int_t missedcycles)
Double_t fPeriodCheckOffset
Definition THcHelicity.h:58
Bool_t fFirstEvProcessed
Definition THcHelicity.h:61
Int_t fFirstCycle
Definition THcHelicity.h:50
Double_t fCycle
Definition THcHelicity.h:59
Int_t fPredictedHelicity
Definition THcHelicity.h:68
Bool_t HWPIN
Definition THcHelicity.h:95
Int_t RanBit30(Int_t ranseed)
Int_t fHelperHistory
virtual void SetHelicityScaler(THcHelicityScaler *f)
Int_t fNLastQuartet
Definition THcHelicity.h:77
Int_t fEvNumCheck
Bool_t fHaveQRT
Definition THcHelicity.h:81
Int_t fHelperQuartetHistory
Long64_t fFirstEvTime
Definition THcHelicity.h:63
virtual ~THcHelicity()
Int_t fHelDelay
Definition THcHelicity.h:89
Int_t fNQuartet
Definition THcHelicity.h:76
Double_t fPeriodCheck
Definition THcHelicity.h:57
Long64_t fLastMPSTime
Definition THcHelicity.h:65
Bool_t fDisabled
Double_t fErrorCode
Double_t fRecommendedFreq
Definition THcHelicity.h:53
Int_t * fHelicityHistory
Double_t fTIPeriod
Definition THcHelicity.h:55
Int_t fRingSeed_reported
Definition THcHelicity.h:84
Int_t fNQRTProblems
Definition THcHelicity.h:82
Bool_t fFoundMPS
Definition THcHelicity.h:72
Double_t fFreq
Definition THcHelicity.h:52
Int_t fQuartetStartPredictedHelicity
Definition THcHelicity.h:71
virtual void Clear(Option_t *opt="")
Int_t GetSeed30(Int_t currentseed)
void Setup(const char *name, const char *description)
void PrintEvent(Int_t evtnum)
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
const char * GetName() const override
const char * GetTitle() const override
virtual void Error(const char *method, const char *msgfmt,...) const
RVec< PromoteTypes< T0, T1 > > fmod(const RVec< T0 > &v, const T1 &y)
RVec< PromoteType< T > > abs(const RVec< T > &v)
RVec< PromoteType< T > > floor(const RVec< T > &v)
Int_t Nint(T x)
Double_t Abs(Double_t d)
STL namespace.