#include "THaG0Helicity.h"
#include "THaEvData.h"
#include "TH1F.h"
#include "TMath.h"
#include <iostream>
using namespace std;
static const Double_t kDefaultTdavg = 14050.;
static const Double_t kDefaultTtol = 40.;
static const Int_t kDefaultDelay = 8;
static const Int_t kDefaultMissQ = 30;
THaG0Helicity::THaG0Helicity( const char* name, const char* description,
THaApparatus* app ) :
THaHelicityDet( name, description, app ),
fG0delay(kDefaultDelay), fEvtype(0), fTdavg(kDefaultTdavg), fTdiff(0.),
fT0(0.), fT9(0.), fT0T9(kFALSE), fQuad_calibrated(kFALSE),
fValidHel(kFALSE), fRecovery_flag(kTRUE),
fTlastquad(0), fTtol(kDefaultTtol), fQuad(0),
fFirstquad(0), fLastTimestamp(0.0), fTimeLastQ1(0.0),
fT9count(0), fPredicted_reading(0), fQ1_reading(0),
fPresent_helicity(kUnknown), fSaved_helicity(kUnknown),
fQ1_present_helicity(kUnknown), fMaxMissed(kDefaultMissQ),
fNqrt(0), fHisto(NULL), fNB(0), fIseed(0), fIseed_earlier(0),
fInquad(0), fTET9Index(0), fTELastEvtQrt(-1), fTELastEvtTime(-1.),
fTELastTime(-1.), fTEPresentReadingQ1(-1), fTEStartup(3), fTETime(-1.),
fTEType9(kFALSE), fManuallySet(0)
{
memset(fHbits, 0, sizeof(fHbits));
}
THaG0Helicity::THaG0Helicity() : fHisto(NULL)
{
}
THaG0Helicity::~THaG0Helicity()
{
DefineVariables( kDelete );
delete fHisto;
}
Int_t THaG0Helicity::DefineVariables( EMode mode )
{
if( mode == kDefine && fIsSetup ) return kOK;
THaHelicityDet::DefineVariables( mode );
const RVarDef var[] = {
{ "qrt", "qrt from ROC", "fQrt" },
{ "nqrt", "number of qrts seen", "fNqrt" },
{ "quad", "quad (1, 2, or 4)", "fQuad" },
{ "tdiff", "time since quad start", "fTdiff" },
{ "gate", "Helicity Gate from ROC", "fGate" },
{ "pread", "Present helicity reading", "fPresentReading" },
{ "timestamp", "Timestamp from ROC", "fTimestamp" },
{ "validtime", "validtime flag", "fValidTime" },
{ "validHel", "validHel flag", "fValidHel" },
{ 0 }
};
return DefineVarsFromList( var, mode );
}
Int_t THaG0Helicity::ReadDatabase( const TDatime& date )
{
Int_t st = THaHelicityDet::ReadDatabase( date );
if( st != kOK )
return st;
st = THaG0HelicityReader::ReadDatabase( GetDBFileName(), GetPrefix(),
date, fG0Debug );
if( st != kOK )
return st;
FILE* file = OpenFile( date );
if( !file ) return kFileError;
Int_t delay = kDefaultDelay;
Double_t tdavg = kDefaultTdavg;
Double_t ttol = kDefaultTtol;
Int_t missqrt = kDefaultMissQ;
DBRequest req[] = {
{ "delay", &delay, kInt, 0, 1, -2 },
{ "tdavg", &tdavg, kDouble, 0, 1, -2 },
{ "ttol", &ttol, kDouble, 0, 1, -2 },
{ "missqrt", &missqrt, kInt, 0, 1, -2 },
{ 0 }
};
st = LoadDB( file, date, req, GetPrefix() );
fclose(file);
if( st )
return kInitError;
if( !TESTBIT(fManuallySet,0) )
fG0delay = delay;
if( !TESTBIT(fManuallySet,1) )
fTdavg = tdavg;
if( !TESTBIT(fManuallySet,2) )
fTtol = ttol;
if( !TESTBIT(fManuallySet,3) )
fMaxMissed = missqrt;
return kOK;
}
Int_t THaG0Helicity::Begin( THaRunBase* )
{
if (fDebug == -1) {
if( !fHisto )
fHisto = new TH1F("hahel1","Time diff for QRT",1000,11000,22000);
else
fHisto->Clear();
}
return 0;
}
void THaG0Helicity::Clear( Option_t* opt )
{
THaHelicityDet::Clear(opt);
fLastTimestamp = fTimestamp;
THaG0HelicityReader::Clear(opt);
fEvtype = fQuad = 0;
fValidHel = kFALSE;
fPresent_helicity = kUnknown;
}
Int_t THaG0Helicity::Decode( const THaEvData& evdata )
{
Int_t err = ReadData( evdata );
if( err ) {
Error( Here("Decode"), "Error decoding helicity data." );
return err;
}
fEvtype = evdata.GetEvType();
if (fDebug >= 3) {
cout << dec << "--> Data for spectrometer " << GetName() << endl;
cout << " evtype " << fEvtype;
cout << " helicity reading " << fPresentReading;
cout << " qrt " << fQrt;
cout << " gate " << fGate;
cout << " time stamp " << fTimestamp << endl;
}
if( fG0delay > 0 ) {
TimingEvent();
QuadCalib();
LoadHelicity();
}
else if( fGate ) {
fPresent_helicity = ( fPresentReading ) ? kPlus : kMinus;
}
if( fSign >= 0 )
fHelicity = fPresent_helicity;
else
fHelicity = ( fPresent_helicity == kPlus ) ? kMinus : kPlus;
if (fDebug >= 3)
cout << "G0 helicity "<< GetName() << " = " << fHelicity <<endl;
return HelicityValid() ? 1 : 0;
}
Int_t THaG0Helicity::End( THaRunBase* )
{
if( fHisto )
fHisto->Write();
return 0;
}
void THaG0Helicity::SetDebug( Int_t level )
{
THaHelicityDet::SetDebug( level );
fG0Debug = level;
}
void THaG0Helicity::SetG0Delay( Int_t value )
{
fG0delay = value;
fManuallySet |= BIT(0);
}
void THaG0Helicity::SetTdavg( Double_t value )
{
fTdavg = value;
fManuallySet |= BIT(1);
}
void THaG0Helicity::SetTtol( Double_t value )
{
fTtol = value;
fManuallySet |= BIT(2);
}
void THaG0Helicity::SetMaxMsQrt( Int_t value )
{
fMaxMissed = value;
fManuallySet |= BIT(3);
}
void THaG0Helicity::TimingEvent()
{
static const char* const here = "TimingEvent";
if (fEvtype == 9) {
fTEType9 = kTRUE;
fTETime = fTimestamp;
}
else if (fQrt && !fTELastEvtQrt
&& fTimestamp - fTELastTime > 8 * fTdavg) {
fTEType9 = kFALSE;
if (fTELastEvtTime>0.
&& (fTimestamp - fTELastEvtTime) < 0.1 * fTdavg) {
fTETime = (fTimestamp + fTELastEvtTime) * 0.5;
} else {
fTETime = fTimestamp;
}
}
else {
fTELastEvtQrt = fQrt;
fTELastEvtTime = fTimestamp;
return;
}
if (fTELastTime > 0)
{
Double_t t9diff = 0.25 * fTdavg;
Double_t tdiff = fTETime - fTELastTime;
Int_t nt9miss = (Int_t) (tdiff / t9diff - 0.5);
fTET9Index += nt9miss;
if (fTET9Index > 3) {
fTEPresentReadingQ1 = -1;
fTET9Index = fTET9Index % 4;
}
if (fTEType9 &&
TMath::Abs(tdiff - (nt9miss+1) * t9diff) > 10*(nt9miss + 1) )
Warning( Here(here), "Weird time difference between timing events: "
"%f at %f.", tdiff, fTETime);
}
++fTET9Index;
if (fQrt) {
if (fTET9Index != 4 && !fTEStartup)
Warning("THaG0Helicity", "QRT wrong spacing: %d at %f.",
fTET9Index, fTETime);
fTET9Index = 0;
fTEPresentReadingQ1 = fPresentReading;
fTEStartup = 0;
}
else {
if (fTEStartup)
--fTEStartup;
if (fTEPresentReadingQ1 > -1)
if ((fPresentReading == fTEPresentReadingQ1 &&
(fTET9Index == 1 || fTET9Index == 2))
|| (fPresentReading != fTEPresentReadingQ1 &&
fTET9Index == 3))
Warning( Here(here), "Wrong helicity pattern: "
"%d in window 1, %d in window %d at timestamp %f.",
fTEPresentReadingQ1, fPresentReading, (fTET9Index+1), fTETime);
}
fTELastTime = fTETime;
fTELastEvtQrt = fQrt;
fTELastEvtTime = fTimestamp;
return;
}
void THaG0Helicity::QuadCalib()
{
static const char* const here = "QuadCalib";
int delayed_qrt = 1;
if (fEvtype == 9) {
fT9count += 1;
if ( delayed_qrt && fQrt ) {
} else {
fT9 = fTimestamp;
}
}
if (fEvtype != 9 && !fGate) return;
fTdiff = fTimestamp - fT0;
if (fDebug >= 3) {
cout << here << " " << fTimestamp<<" "<<fT0;
cout << " "<<fTdiff
<< " " << fEvtype
<<" "<<fQrt<<endl;
}
if (fFirstquad == 0 &&
fTdiff > (1.25*fTdavg + fTtol)) {
Int_t nqmiss = (Int_t)(fTdiff/fTdavg);
if (fDebug >= 1)
Info( Here(here), "Recovering large DT, nqmiss = %d "
"at timestamp %f, tdiff %f", nqmiss, fTimestamp, fTdiff );
if (fQuad_calibrated && nqmiss < fMaxMissed) {
for (Int_t i = 0; i < nqmiss; i++) {
if (fQrt && fT9 > 0 && fTimestamp - fT9 < 8*fTdavg) {
fT0 = TMath::Floor((fTimestamp-fT9)/(.25*fTdavg))
*(.25*fTdavg) + fT9;
fT0T9 = kTRUE;
}
else {
fT0 += fTdavg;
fT0T9 = kFALSE;
}
QuadHelicity(1);
fNqrt++;
fQ1_reading = (fPredicted_reading == -1) ? 0 : 1;
fQ1_present_helicity = fPresent_helicity;
if (fDebug>=1) {
Info(Here(here)," %5d M M %1d %2d %10.0f %10.0f %10.0f Missing",
fNqrt,fQ1_reading,fQ1_present_helicity,fTimestamp,fT0,fTdiff);
}
}
fTdiff = fTimestamp - fT0;
} else {
Warning( Here(here), "Cannot recover: Too many skipped QRT (Low rate ?) "
"or uninitialized" );
fNqrt += nqmiss;
fT0 += nqmiss*fTdavg;
fRecovery_flag = kTRUE;
fQuad_calibrated = kFALSE;
fFirstquad = 1;
}
}
if (fQrt && fFirstquad == 1) {
fT0 = fTimestamp;
fT0T9 = kFALSE;
fFirstquad = 0;
}
if (fEvtype == 9 && fQrt) {
if (fTimeLastQ1 > 0) {
Double_t q1tdiff = fmod(fTimestamp - fTimeLastQ1,fTdavg);
if (q1tdiff > .1*fTdavg)
Warning( Here(here), "QRT==1 timing error -- Last time = %f "
"This time = %f\nHELICITY SIGNALS MAY BE CORRUPT",
fTimeLastQ1, fTimestamp );
}
fTimeLastQ1 = fTimestamp;
}
if (fQrt)
fT9count = 0;
if ( ( fTdiff > 0.5*fTdavg ) && fQrt ) {
if (fDebug >= 3)
Info(Here(here),"found qrt ");
if (fT9 > 0 && fTimestamp - fT9 < 8*fTdavg) {
fT0 = TMath::Floor((fTimestamp-fT9)/(.25*fTdavg))
*(.25*fTdavg) + fT9;
fT0T9 = kTRUE;
}
else {
fT0 = (fLastTimestamp == 0
|| fTimestamp - fLastTimestamp > 0.1*fTdavg)
? fTimestamp : 0.5 * (fTimestamp + fLastTimestamp);
fT0T9 = kFALSE;
}
fQ1_reading = fPresentReading;
QuadHelicity();
fNqrt++;
fQ1_present_helicity = fPresent_helicity;
if (fQuad_calibrated && fDebug >= 3)
cout << "------------ quad calibrated --------------------------"<<endl;
if (fQuad_calibrated && !CompHel() ) {
Warning( Here(here), "QuadCalib G0 prediction failed at timestamp %f.\n"
"HELICITY ASSIGNMENT IN PREVIOUS %d WINDOWS MAY BE INCORRECT",
fTimestamp, fG0delay);
if (fDebug >= 3) {
cout << "Qrt " << fQrt << " Tdiff "<<fTdiff;
cout<<" gate "<<fGate<<" evtype "<<fEvtype<<endl;
}
fRecovery_flag = kTRUE;
fQuad_calibrated = kFALSE;
fFirstquad = 1;
fPresent_helicity = kUnknown;
}
if (fDebug>=3) {
Info(Here(here), "%5d %1d %1d %1d %2d %10.0f %10.0f %10.0f",fNqrt,
fEvtype,fQrt,fQ1_reading,fQ1_present_helicity,fTimestamp,fT0,fTdiff);
}
if (fDebug==-1) {
cout << fTdiff<<endl;
if (fHisto) {
fHisto->Fill(fTdiff);
}
}
}
fTdiff = fTimestamp - fT0;
return;
}
void THaG0Helicity::LoadHelicity()
{
static const char* const here = "LoadHelicity";
if ( fQuad_calibrated )
fValidHel = kTRUE;
if ( !fQuad_calibrated || !fGate ) {
fPresent_helicity = kUnknown;
return;
}
if (fDebug >= 2)
cout << "Loading helicity ##########" << endl;
fInquad = 0;
if (fQrt) {
fInquad = 1;
if (fPresentReading != fQ1_reading) {
Warning( Here(here), "Invalid bit pattern !! "
"fPresentReading != fQ1_reading while QRT == 1 at timestamp "
"%f.\nHELICITY ASSIGNMENT MAY BE INCORRECT", fTimestamp);
}
}
if (!fQrt && fPresentReading != fQ1_reading) {
fInquad = 2;
}
if (!fQrt && fPresentReading == fQ1_reading) {
fInquad = 4;
}
if (fInquad == 0) {
fPresent_helicity = kUnknown;
Error( Here(here), "Quad assignment impossible !! "
" qrt = %d present read = %d Q1 read = %d at timestamp %f.\n"
"HELICITY SET TO UNKNOWN",
fQrt, fPresentReading, fQ1_reading, fTimestamp );
}
else if (fInquad == 1 || fInquad >= 4) {
fPresent_helicity = fQ1_present_helicity;
} else {
if (fQ1_present_helicity == kPlus)
fPresent_helicity = kMinus;
if (fQ1_present_helicity == kMinus)
fPresent_helicity = kPlus;
}
if (
(!fT0T9 &&
((fInquad == 1 && fTdiff > 0.5 * fTdavg) ||
(fInquad == 4 && fTdiff < 0.5 * fTdavg)))
||
(fT0T9 &&
((fInquad == 1 && fTdiff > 0.3 * fTdavg) ||
(fInquad == 2 && (fTdiff < 0.2 * fTdavg ||
fTdiff > 0.8 * fTdavg)) ||
(fInquad == 4 && fTdiff < 0.7 * fTdavg))))
{
Warning( Here(here), "Quad assignment inconsistent with timing !! "
"quad = %d tdiff = %f at timestamp %f. "
"HELICITY ASSIGNMENT MAY BE INCORRECT",
fInquad, fTdiff, fTimestamp );
}
if (fDebug >= 2) {
cout << dec << "Quad "<<fInquad;
cout << " Time diff "<<fTdiff;
cout << " Qrt "<<fQrt;
cout << " Helicity "<<fPresent_helicity<<endl;
}
fQuad = fInquad;
return;
}
void THaG0Helicity::QuadHelicity( Int_t cond )
{
static const char* const here = "QuadHelicity";
if (fRecovery_flag) {
fNB = 0;
fSaved_helicity = kUnknown;
}
fRecovery_flag = kFALSE;
fPresent_helicity = fSaved_helicity;
if (cond == 0 && fNB>0
&& fTlastquad > 0
&& fT0 - fTlastquad < 0.3 * fTdavg) {
if (fDebug >= 2) Warning(Here(here),"SKIP this quad, QRT's too close");
return;
}
fTlastquad = fT0;
if (fNB < kNbits) {
fHbits[fNB] = fPresentReading;
fNB++;
fQuad_calibrated = kFALSE;
} else if (fNB == kNbits) {
fIseed_earlier = GetSeed();
fIseed = fIseed_earlier;
for( int i = 0; i < kNbits+1; i++) {
fPredicted_reading = RanBit(0);
RanBit(1);
}
if( fPredicted_reading > 0 )
fPresent_helicity = kPlus;
else if( fPredicted_reading < 0 )
fPresent_helicity = kMinus;
else
fPresent_helicity = kUnknown;
for( int i = 0; i < fG0delay/4; i++)
fPresent_helicity = RanBit(1);
fNB++;
fSaved_helicity = fPresent_helicity;
fQuad_calibrated = kTRUE;
} else {
fPredicted_reading = RanBit(0);
fSaved_helicity = fPresent_helicity = RanBit(1);
if (fDebug >= 3) {
cout << "quad helicities " << dec;
cout << fPredicted_reading;
cout << " ?=? " << fPresentReading;
cout << " pred-> " << fPresent_helicity
<< " time " << fTimestamp
<< " <<<<<<<<<<<<<================================="
<<endl;
if (CompHel()) cout << "HELICITY OK "<<endl;
}
fQuad_calibrated = kTRUE;
}
return;
}
THaHelicityDet::EHelicity THaG0Helicity::RanBit( int which )
{
const int MASK = BIT(0)+BIT(2)+BIT(3)+BIT(23);
if (which == 0) {
if( TESTBIT(fIseed_earlier,23) ) {
fIseed_earlier = ((fIseed_earlier^MASK)<<1) | BIT(0);
return kPlus;
} else {
fIseed_earlier <<= 1;
return kMinus;
}
} else {
if( TESTBIT(fIseed,23) ) {
fIseed = ((fIseed^MASK)<<1) | BIT(0);
return kPlus;
} else {
fIseed <<= 1;
return kMinus;
}
}
}
UInt_t THaG0Helicity::GetSeed()
{
int seedbits[kNbits];
UInt_t ranseed = 0;
for (int i = 0; i < 20; i++)
seedbits[23-i] = fHbits[i];
seedbits[3] = fHbits[20]^seedbits[23];
seedbits[2] = fHbits[21]^seedbits[22]^seedbits[23];
seedbits[1] = fHbits[22]^seedbits[21]^seedbits[22];
seedbits[0] = fHbits[23]^seedbits[20]^seedbits[21]^seedbits[23];
for (int i=kNbits-1; i >= 0; i--)
ranseed = ranseed<<1|(seedbits[i]&1);
ranseed &= 0xFFFFFF;
return ranseed;
}
Bool_t THaG0Helicity::CompHel()
{
if (fPresentReading == 0 &&
fPredicted_reading == kMinus) return kTRUE;
if (fPresentReading == 1 &&
fPredicted_reading == kPlus) return kTRUE;
return kFALSE;
}
ClassImp(THaG0Helicity)
Last change: Sat Nov 7 21:26:47 2009
Last generated: 2009-11-07 21:26
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.