#include "THaScintillator.h"
#include "THaEvData.h"
#include "THaDetMap.h"
#include "VarDef.h"
#include "VarType.h"
#include "THaTrack.h"
#include "TClonesArray.h"
#include "TMath.h"
#include "THaTrackProj.h"
#include <cstring>
#include <cstdio>
#include <cstdlib>
#include <iostream>
using namespace std;
THaScintillator::THaScintillator( const char* name, const char* description,
THaApparatus* apparatus ) :
THaNonTrackingDetector(name,description,apparatus)
{
fTWalkPar = 0;
fTrackProj = new TClonesArray( "THaTrackProj", 5 );
}
THaScintillator::THaScintillator( ) :
THaNonTrackingDetector()
{
fTWalkPar = NULL;
fTrackProj = NULL;
fRA_c = fRA_p = fRA = fLA_c = fLA_p = fLA = NULL;
fRT_c = fRT = fLT_c = fLT = NULL;
fRGain = fLGain = fRPed = fLPed = fROff = fLOff = NULL;
fTrigOff = fTime = fdTime = fYt = fYa = NULL;
fHitPad = NULL;
}
THaAnalysisObject::EStatus THaScintillator::Init( const TDatime& date )
{
if( THaNonTrackingDetector::Init( date ) )
return fStatus;
const DataDest tmp[NDEST] = {
{ &fRTNhit, &fRANhit, fRT, fRT_c, fRA, fRA_p, fRA_c, fROff, fRPed, fRGain },
{ &fLTNhit, &fLANhit, fLT, fLT_c, fLA, fLA_p, fLA_c, fLOff, fLPed, fLGain }
};
memcpy( fDataDest, tmp, NDEST*sizeof(DataDest) );
return fStatus = kOK;
}
Int_t THaScintillator::ReadDatabase( const TDatime& date )
{
static const char* const here = "ReadDatabase()";
const int LEN = 200;
char buf[LEN];
Int_t nelem;
FILE* fi = OpenFile( date );
if( !fi ) return kFileError;
while ( ReadComment( fi, buf, LEN ) ) {}
fscanf ( fi, "%d", &nelem );
fgets ( buf, LEN, fi );
if( fIsInit && nelem != fNelem ) {
Error( Here(here), "Cannot re-initalize with different number of paddles. "
"(was: %d, now: %d). Detector not re-initialized.", fNelem, nelem );
fclose(fi);
return kInitError;
}
fNelem = nelem;
while ( ReadComment( fi, buf, LEN ) ) {}
int i = 0;
fDetMap->Clear();
while (1) {
int pos;
Int_t first_chan, model;
Int_t crate, slot, first, last;
fgets ( buf, LEN, fi );
sscanf( buf, "%d%d%d%d%d%n", &crate, &slot, &first, &last, &first_chan, &pos );
if( crate < 0 ) break;
model=atoi(buf+pos);
if( fDetMap->AddModule( crate, slot, first, last, first_chan, model ) < 0 ) {
Error( Here(here), "Too many DetMap modules (maximum allowed - %d).",
THaDetMap::kDetMapSize);
fclose(fi);
return kInitError;
}
}
while ( ReadComment( fi, buf, LEN ) ) {}
Float_t x,y,z;
fscanf ( fi, "%f%f%f", &x, &y, &z );
fgets ( buf, LEN, fi );
fOrigin.SetXYZ( x, y, z );
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
fscanf ( fi, "%f%f%f", fSize, fSize+1, fSize+2 );
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
Float_t angle;
fscanf ( fi, "%f", &angle );
fgets ( buf, LEN, fi );
const Float_t degrad = TMath::Pi()/180.0;
tan_angle = TMath::Tan(angle*degrad);
sin_angle = TMath::Sin(angle*degrad);
cos_angle = TMath::Cos(angle*degrad);
DefineAxes(angle*degrad);
if( !fIsInit ) {
fLOff = new Double_t[ fNelem ];
fROff = new Double_t[ fNelem ];
fLPed = new Double_t[ fNelem ];
fRPed = new Double_t[ fNelem ];
fLGain = new Double_t[ fNelem ];
fRGain = new Double_t[ fNelem ];
fTrigOff = new Double_t[ fNelem ];
fLT = new Double_t[ fNelem ];
fLT_c = new Double_t[ fNelem ];
fRT = new Double_t[ fNelem ];
fRT_c = new Double_t[ fNelem ];
fLA = new Double_t[ fNelem ];
fLA_p = new Double_t[ fNelem ];
fLA_c = new Double_t[ fNelem ];
fRA = new Double_t[ fNelem ];
fRA_p = new Double_t[ fNelem ];
fRA_c = new Double_t[ fNelem ];
fNTWalkPar = 2*fNelem;
fTWalkPar = new Double_t[ fNTWalkPar ];
fHitPad = new Int_t[ fNelem ];
fTime = new Double_t[ fNelem ];
fdTime = new Double_t[ fNelem ];
fAmpl = new Double_t[ fNelem ];
fYt = new Double_t[ fNelem ];
fYa = new Double_t[ fNelem ];
fIsInit = true;
}
memset(fTrigOff,0,fNelem*sizeof(fTrigOff[0]));
fTdc2T = 0.1e-9;
fResolution = fTdc2T;
fCn = 1.7e+8;
fAttenuation = 0.7;
fAdcMIP = 1.e10;
for (int i=0; i<fNTWalkPar; i++) fTWalkPar[i]=0;
for (int i=0; i<fNelem; i++) fTrigOff[i]=0;
DBRequest list[] = {
{ "TDC_offsetsL", fLOff, kDouble, fNelem },
{ "TDC_offsetsR", fROff, kDouble, fNelem },
{ "ADC_pedsL", fLPed, kDouble, fNelem },
{ "ADC_pedsR", fRPed, kDouble, fNelem },
{ "ADC_coefL", fLGain, kDouble, fNelem },
{ "ADC_coefR", fRGain, kDouble, fNelem },
{ "TDC_res", &fTdc2T },
{ "TransSpd", &fCn },
{ "AdcMIP", &fAdcMIP },
{ "NTWalk", &fNTWalkPar, kInt },
{ "Timewalk", fTWalkPar, kDouble, 2*fNelem },
{ "ReTimeOff", fTrigOff, kDouble, fNelem },
{ "AvgRes", &fResolution },
{ "Atten", &fAttenuation },
{ 0 }
};
#if 0
if ( gHaDB && gHaDB->LoadValues(GetPrefix(),list,date) ) {
goto exit;
}
#endif
if( SeekDBdate( fi, date ) == 0 && fConfig.Length() > 0 &&
SeekDBconfig( fi, fConfig.Data() )) {}
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf(fi,"%lf",fLOff+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf(fi,"%lf",fROff+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf(fi,"%lf",fLPed+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf(fi,"%lf",fRPed+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf (fi,"%lf",fLGain+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
for (i=0;i<fNelem;i++)
fscanf (fi,"%lf",fRGain+i);
fgets ( buf, LEN, fi );
while ( ReadComment( fi, buf, LEN ) ) {}
if ( ! fgets(buf, LEN, fi) || strchr(buf,'[') ) goto exit;
sscanf(buf,"%lf",&fTdc2T);
fResolution = 3.*fTdc2T;
while ( ReadComment( fi, buf, LEN ) ) {}
if ( !fgets(buf, LEN, fi) || strchr(buf,'[') ) goto exit;
sscanf(buf,"%lf",&fCn);
while ( ReadComment( fi, buf, LEN ) ) {}
if ( !fgets ( buf, LEN, fi ) || strchr(buf,'[') ) goto exit;
sscanf(buf,"%lf",&fAttenuation);
while ( ReadComment( fi, buf, LEN ) ) {}
if ( !fgets(buf, LEN, fi) || strchr(buf,'[') ) goto exit;
sscanf(buf,"%lf",&fAdcMIP);
while ( ReadComment( fi, buf, LEN ) ) {}
{
int cnt=0;
while ( cnt<fNTWalkPar && fgets( buf, LEN, fi ) && ! strchr(buf,'[') ) {
char *ptr = buf;
int pos=0;
while ( cnt < fNTWalkPar && sscanf(ptr,"%lf%n",&fTWalkPar[cnt],&pos)>0 ) {
ptr += pos;
cnt++;
}
}
}
while ( ReadComment( fi, buf, LEN ) ) {}
{
int cnt=0;
while ( cnt<fNelem && fgets( buf, LEN, fi ) && ! strchr(buf,'[') ) {
char *ptr = buf;
int pos=0;
while ( cnt < fNelem && sscanf(ptr,"%lf%n",&fTrigOff[cnt],&pos)>0 ) {
ptr += pos;
cnt++;
}
}
}
exit:
fclose(fi);
if ( fDebug > 1 ) {
cout << '\n' << GetPrefix() << " calibration parameters: " << endl;;
for ( DBRequest *li = list; li->name; li++ ) {
cout << " " << li->name;
int maxc = li->nelem;
if (maxc==0)maxc=1;
for (int i=0; i<maxc; i++) {
if (li->type==kDouble) cout << " " << ((Double_t*)li->var)[i];
if (li->type==kInt) cout << " " << ((Int_t*)li->var)[i];
}
cout << endl;
}
}
return kOK;
}
Int_t THaScintillator::DefineVariables( EMode mode )
{
if( mode == kDefine && fIsSetup ) return kOK;
fIsSetup = ( mode == kDefine );
RVarDef vars[] = {
{ "nlthit", "Number of Left paddles TDC times", "fLTNhit" },
{ "nrthit", "Number of Right paddles TDC times", "fRTNhit" },
{ "nlahit", "Number of Left paddles ADCs amps", "fLANhit" },
{ "nrahit", "Number of Right paddles ADCs amps", "fRANhit" },
{ "lt", "TDC values left side", "fLT" },
{ "lt_c", "Corrected times left side", "fLT_c" },
{ "rt", "TDC values right side", "fRT" },
{ "rt_c", "Corrected times right side", "fRT_c" },
{ "la", "ADC values left side", "fLA" },
{ "la_p", "Corrected ADC values left side", "fLA_p" },
{ "la_c", "Corrected ADC values left side", "fLA_c" },
{ "ra", "ADC values right side", "fRA" },
{ "ra_p", "Corrected ADC values right side", "fRA_p" },
{ "ra_c", "Corrected ADC values right side", "fRA_c" },
{ "nthit", "Number of paddles with l&r TDCs", "fNhit" },
{ "t_pads", "Paddles with l&r coincidence TDCs", "fHitPad" },
{ "y_t", "y-position from timing (m)", "fYt" },
{ "y_adc", "y-position from amplitudes (m)", "fYa" },
{ "time", "Time of hit at plane (s)", "fTime" },
{ "dtime", "Est. uncertainty of time (s)", "fdTime" },
{ "dedx", "dEdX-like deposited in paddle", "fAmpl" },
{ "troff", "Trigger offset for paddles", "fTrigOff"},
{ "trn", "Number of tracks for hits", "GetNTracks()" },
{ "trx", "x-position of track in det plane", "fTrackProj.THaTrackProj.fX" },
{ "try", "y-position of track in det plane", "fTrackProj.THaTrackProj.fY" },
{ "trpath", "TRCS pathlen of track to det plane","fTrackProj.THaTrackProj.fPathl" },
{ "trdx", "track deviation in x-position (m)", "fTrackProj.THaTrackProj.fdX" },
{ "trpad", "paddle-hit associated with track", "fTrackProj.THaTrackProj.fChannel" },
{ 0 }
};
return DefineVarsFromList( vars, mode );
}
THaScintillator::~THaScintillator()
{
if( fIsSetup )
RemoveVariables();
if( fIsInit )
DeleteArrays();
if (fTrackProj) {
fTrackProj->Clear();
delete fTrackProj; fTrackProj = 0;
}
}
void THaScintillator::DeleteArrays()
{
delete [] fRA_c; fRA_c = NULL;
delete [] fRA_p; fRA_p = NULL;
delete [] fRA; fRA = NULL;
delete [] fLA_c; fLA_c = NULL;
delete [] fLA_p; fLA_p = NULL;
delete [] fLA; fLA = NULL;
delete [] fRT_c; fRT_c = NULL;
delete [] fRT; fRT = NULL;
delete [] fLT_c; fLT_c = NULL;
delete [] fLT; fLT = NULL;
delete [] fRGain; fRGain = NULL;
delete [] fLGain; fLGain = NULL;
delete [] fRPed; fRPed = NULL;
delete [] fLPed; fLPed = NULL;
delete [] fROff; fROff = NULL;
delete [] fLOff; fLOff = NULL;
delete [] fTWalkPar; fTWalkPar = NULL;
delete [] fTrigOff; fTrigOff = NULL;
delete [] fHitPad; fHitPad = NULL;
delete [] fTime; fTime = NULL;
delete [] fdTime; fdTime = NULL;
delete [] fYt; fYt = NULL;
delete [] fYa; fYa = NULL;
}
inline
void THaScintillator::ClearEvent()
{
const int lf = fNelem*sizeof(Double_t);
fLTNhit = 0;
memset( fLT, 0, lf );
memset( fLT_c, 0, lf );
fRTNhit = 0;
memset( fRT, 0, lf );
memset( fRT_c, 0, lf );
fLANhit = 0;
memset( fLA, 0, lf );
memset( fLA_p, 0, lf );
memset( fLA_c, 0, lf );
fRANhit = 0;
memset( fRA, 0, lf );
memset( fRA_p, 0, lf );
memset( fRA_c, 0, lf );
fNhit = 0;
memset( fHitPad, 0, fNelem*sizeof(fHitPad[0]) );
memset( fTime, 0, lf );
memset( fdTime, 0, lf );
memset( fYt, 0, lf );
memset( fYa, 0, lf );
fTrackProj->Clear();
}
Int_t THaScintillator::Decode( const THaEvData& evdata )
{
ClearEvent();
for( Int_t i = 0; i < fDetMap->GetSize(); i++ ) {
THaDetMap::Module* d = fDetMap->GetModule( i );
bool adc = ( d->model ? fDetMap->IsADC(d) : (i < fDetMap->GetSize()/2) );
for( Int_t j = 0; j < evdata.GetNumChan( d->crate, d->slot ); j++) {
Int_t chan = evdata.GetNextChan( d->crate, d->slot, j );
if( chan < d->lo || chan > d->hi ) continue;
Int_t data = evdata.GetData( d->crate, d->slot, chan, 0 );
Int_t k = d->first + chan - d->lo - 1;
#ifdef WITH_DEBUG
if( k<0 || k>NDEST*fNelem ) {
Warning( Here("Decode()"), "Illegal detector channel: %d", k );
continue;
}
#endif
DataDest* dest = fDataDest + k/fNelem;
k = k % fNelem;
if( adc ) {
dest->adc[k] = static_cast<Double_t>( data );
dest->adc_p[k] = data - dest->ped[k];
dest->adc_c[k] = dest->adc_p[k] * dest->gain[k];
(*dest->nahit)++;
} else {
dest->tdc[k] = static_cast<Double_t>( data );
dest->tdc_c[k] = (data - dest->offset[k])*fTdc2T;
(*dest->nthit)++;
}
}
}
if ( fDebug > 3 ) {
printf("\n\nEvent %d Trigger %d Scintillator %s\n:",
evdata.GetEvNum(), evdata.GetEvType(), GetPrefix() );
printf(" paddle Left(TDC ADC ADC_p) Right(TDC ADC ADC_p)\n");
for ( int i=0; i<fNelem; i++ ) {
printf(" %2d %5.0f %5.0f %5.0f %5.0f %5.0f %5.0f\n",
i+1,fLT[i],fLA[i],fLA_p[i],fRT[i],fRA[i],fRA_p[i]);
}
}
return fLTNhit+fRTNhit;
}
Int_t THaScintillator::ApplyCorrections( void )
{
Int_t nlt=0, nrt=0, nla=0, nra=0;
for (Int_t i=0; i<fNelem; i++) {
if (fLA[i] != 0.) {
fLA_p[i] = fLA[i] - fLPed[i];
fLA_c[i] = fLA_p[i]*fLGain[i];
nla++;
}
if (fRA[i] != 0.) {
fRA_p[i] = fRA[i] - fRPed[i];
fRA_c[i] = fRA_p[i]*fRGain[i];
nra++;
}
if (fLT[i] != 0.) {
fLT_c[i] = (fLT[i] - fLOff[i])*fTdc2T - TimeWalkCorrection(i,kLeft);
nlt++;
}
if (fRT[i] != 0.) {
fRT_c[i] = (fRT[i] - fROff[i])*fTdc2T - TimeWalkCorrection(i,kRight);
nrt++;
}
}
return !(fLTNhit==nlt && fLANhit==nla && fRTNhit==nrt && fRANhit==nra );
}
Double_t THaScintillator::TimeWalkCorrection(const Int_t& paddle,
const ESide side)
{
Double_t adc=0;
if (side == kLeft)
adc = fLA_p[paddle];
else
adc = fRA_p[paddle];
if (fNTWalkPar<=0 || !fTWalkPar) return 0.;
if ( adc <=0. ) return 0.;
Double_t ref = fAdcMIP;
Double_t tw(0), tw_ref(0.);
int npar = fNTWalkPar/(2*fNelem);
Double_t *par = &(fTWalkPar[npar*(side*fNelem+paddle)]);
tw = par[0]*pow(adc,-.5);
tw_ref = par[0]*pow(ref,-.5);
return tw-tw_ref;
}
Int_t THaScintillator::CoarseProcess( TClonesArray& )
{
static const Double_t sqrt2 = TMath::Sqrt(2.);
ApplyCorrections();
fNhit = 0;
for (int i=0; i<fNelem; i++) {
if (fLT[i]>0 && fRT[i]>0) {
fHitPad[fNhit++] = i;
fTime[i] = .5*(fLT_c[i]+fRT_c[i])-fSize[1]/fCn;
fdTime[i] = fResolution/sqrt2;
fYt[i] = .5*fCn*(fRT_c[i]-fLT_c[i]);
}
if (fLA_c[i]>0&&fRA_c[i]>0) {
fYa[i] = TMath::Log(fLA_c[i]/fRA_c[i])/(2.*fAttenuation);
fAmpl[i] = TMath::Sqrt(fLA_c[i]*fRA_c[i]*TMath::Exp(fAttenuation*2*fSize[1]))
/ fSize[2];
}
}
return 0;
}
Int_t THaScintillator::FineProcess( TClonesArray& tracks )
{
int n_track = tracks.GetLast()+1;
Double_t dpadx = (2.*fSize[0])/(fNelem);
Double_t padx0 = -dpadx*(fNelem-1)*.5;
for ( int i=0; i<n_track; i++ ) {
THaTrack* theTrack = static_cast<THaTrack*>( tracks[i] );
Double_t pathl=kBig, xc=kBig, yc=kBig, dx=kBig;
Int_t pad=-1;
if ( ! CalcTrackIntercept(theTrack, pathl, xc, yc) ) {
new ( (*fTrackProj)[i] )
THaTrackProj(xc,yc,pathl,dx,pad,this);
continue;
}
dx = kBig;
for ( Int_t j=0; j<fNhit; j++ ) {
Double_t dx2 = ( padx0 + fHitPad[j]*dpadx) - xc;
if (TMath::Abs(dx2) < TMath::Abs(dx) ) {
pad = fHitPad[j];
dx = dx2;
}
else if (pad>=0) break;
}
new ( (*fTrackProj)[i] )
THaTrackProj(xc,yc,pathl,dx,pad,this);
}
return 0;
}
ClassImp(THaScintillator)
Last change: Sat Nov 7 21:26:51 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.