//$Id$ // Author: Anders Vestbo //*-- Copyright © ASV #include #include #include #include "AliL3Modeller.h" #include "AliL3MemHandler.h" #include "AliL3TrackArray.h" #include "AliL3ModelTrack.h" #include "AliL3DigitData.h" #include "AliL3Transform.h" //_____________________________________________________________ // AliL3Modeller // // Class for modeling TPC data. // // This performs the cluster finding, based on track parameters. // Basically it propagates the tracks to all padrows, and looks // for a corresponding cluster. For the moment only cog is calculated, // and no deconvolution is done. ClassImp(AliL3Modeller) AliL3Modeller::AliL3Modeller() { fMemHandler=0; fTracks=0; fTrackThreshold=0; fPadOverlap=0; fTimeOverlap=0; fRowData=0; } AliL3Modeller::~AliL3Modeller() { if(fMemHandler) delete fMemHandler; if(fTracks) delete fTracks; } void AliL3Modeller::Init(Int_t slice,Int_t patch,Char_t *trackdata,Char_t *path,Bool_t houghtracks) { fSlice = slice; fPatch = patch; fPadOverlap=6; fTimeOverlap=8; sprintf(fPath,"%s",path); AliL3Transform::Init(fPath); fTracks = new AliL3TrackArray("AliL3ModelTrack"); Char_t fname[100]; AliL3MemHandler *file = new AliL3MemHandler(); sprintf(fname,"%s/tracks_tr_%d_0.raw",trackdata,fSlice); //output tracks from the tracker (no merging) if(!file->SetBinaryInput(fname)) { cerr<<"AliL3Modeller::Init : Error opening trackfile: "<Binary2TrackArray(fTracks); file->CloseBinaryInput(); delete file; if(houghtracks) cout<<"AliL3Modeller is assuming local hough tracksegments!"<QSort(); for(Int_t i=0; iGetNTracks(); i++) { AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i); if(!track) continue; track->Init(fSlice,fPatch); //Only if the tracks has been merged across sector boundaries: //if(!houghtracks) //track->Rotate(fSlice,kTRUE); //!!!!!!!!!!!!!!!!!!! track->CalculateHelix(); } CalculateCrossingPoints(); if(!houghtracks) CheckForOverlaps(); fMemHandler = new AliL3MemHandler(); sprintf(fname,"%sdigits_%d_%d.raw",fPath,fSlice,fPatch); if(!fMemHandler->SetBinaryInput(fname)) { cerr<<"AliL3Modeller::Init : Error opening file "<CompBinary2Memory(ndigits); SetInputData(digits); } void AliL3Modeller::FindClusters() { if(!fTracks) { cerr<<"AliL3Modeller::Process : No tracks"<fDigitData; for(UInt_t j=0; jfNDigit; j++) { pad = digPt[j].fPad; time = digPt[j].fTime; charge = digPt[j].fCharge; row[ntimes*pad+time].fCharge = charge; row[ntimes*pad+time].fUsed = kFALSE; //cout<<"Row "<GetNTracks(); k++) { AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(k); if(!track) continue; //if(track->GetOverlap(i)>=0) continue;//Track is overlapping if(track->GetPadHit(i)<0 || track->GetTimeHit(i)<0 || track->GetOverlap(i)>=0) { track->SetCluster(i,0,0,0,0,0,0); //The track has left the patch. continue; } Int_t hitpad = (Int_t)rint(track->GetPadHit(i)); Int_t hittime = (Int_t)rint(track->GetTimeHit(i)); //cout<<"Checking track on row "<= AliL3Transform::GetNPads(i)) { //cout<<"Pad = "<= bounds) { cerr<<"AliL3Modeller::FindClusters : Index out of range : "< 0) { pad--; //In this case, we haven't found anything yet, } //so we will try to expand our search within the natural boundaries. else { pad=hitpad+1; padsign=1; } continue; } else if(padsign==1) { if(cluster.fCharge==0 && abs(pad-hitpad) <= fPadOverlap/2 && pad < AliL3Transform::GetNPads(i)-2) { pad++; //In this case, we haven't found anything yet, continue; //so we will try to expand our search within the natural boundaries. } else //We are out of range, or cluster if finished. { //cout<<"Outof range; charge "<UpdateRowPointer(rowPt); } delete [] row; cout<<"done processing"<GetNTracks(); i++) { AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i); if(!track) continue; if(track->GetNClusters() != AliL3Transform::GetNRows(fPatch)) cerr<GetNClusters()<<" nrows "<fCharge==0) { track->SetCluster(row,0,0,0,0,0,0); return; } Float_t fcharge = (Float_t)cluster->fCharge; Float_t fpad = ((Float_t)cluster->fPad/fcharge); Float_t ftime = ((Float_t)cluster->fTime/fcharge); Float_t sigmaY2,sigmaZ2; CalcClusterWidth(cluster,sigmaY2,sigmaZ2); track->SetCluster(row,fpad,ftime,fcharge,sigmaY2,sigmaZ2,npads); } void AliL3Modeller::FillZeros(AliL3DigitRowData *rowPt,Digit *row) { //Fill zero where data has been used. Int_t ntimes = AliL3Transform::GetNTimeBins()+1; AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData; for(UInt_t j=0; jfNDigit; j++) { Int_t pad = digPt[j].fPad; Int_t time = digPt[j].fTime; if(row[ntimes*pad+time].fUsed==kTRUE) digPt[j].fCharge = 0; } } void AliL3Modeller::WriteRemaining() { //Write remaining (nonzero) digits to file. AliL3DigitRowData *rowPt; rowPt = (AliL3DigitRowData*)fRowData; Int_t digitcount=0; Int_t ndigits[(AliL3Transform::GetNRows(fPatch))]; for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++) { AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData; ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]=0; for(UInt_t j=0; jfNDigit; j++) { if(digPt[j].fCharge==0) continue; digitcount++; ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]++; } //cout<<"Difference "<<(int)ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]<<" "<<(int)rowPt->fNDigit<UpdateRowPointer(rowPt); } Int_t size = digitcount*sizeof(AliL3DigitData) + AliL3Transform::GetNRows(fPatch)*sizeof(AliL3DigitRowData); Byte_t *data = new Byte_t[size]; memset(data,0,size); AliL3DigitRowData *tempPt = (AliL3DigitRowData*)data; rowPt = (AliL3DigitRowData*)fRowData; for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++) { Int_t localcount=0; tempPt->fRow = i; tempPt->fNDigit = ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]; AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData; for(UInt_t j=0; jfNDigit; j++) { if(digPt[j].fCharge==0) continue; if(localcount >= ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]) { cerr<<"AliL3Modeller::WriteRemaining : Digitarray out of range!!"<fDigitData[localcount].fCharge = digPt[j].fCharge; tempPt->fDigitData[localcount].fPad = digPt[j].fPad; tempPt->fDigitData[localcount].fTime = digPt[j].fTime; localcount++; } if(ndigits[(i-AliL3Transform::GetFirstRow(fPatch))] != localcount) { cerr<<"AliL3Modeller::WriteRemaining : Mismatch in digitcount"<UpdateRowPointer(rowPt); Byte_t *tmp = (Byte_t*)tempPt; Int_t size = sizeof(AliL3DigitRowData) + ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]*sizeof(AliL3DigitData); tmp += size; tempPt = (AliL3DigitRowData*)tmp; } Char_t fname[100]; AliL3MemHandler *mem = new AliL3MemHandler(); sprintf(fname,"%s/remains_%d_%d.raw",fPath,fSlice,fPatch); mem->SetBinaryOutput(fname); mem->Memory2CompBinary((UInt_t)AliL3Transform::GetNRows(fPatch),(AliL3DigitRowData*)data); mem->CloseBinaryOutput(); delete mem; delete [] data; } void AliL3Modeller::CalculateCrossingPoints() { cout<<"Calculating crossing points on "<GetNTracks()<<" tracks"<GetNTracks(); i++) { AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i); if(!track) continue; if(track->GetFirstPointX() > AliL3Transform::Row2X(AliL3Transform::GetLastRow(fPatch)) || track->GetPt()<0.1) fTracks->Remove(i); if(track->GetNHits() < fTrackThreshold) { fTracks->Remove(i); continue; } } Int_t sector,row; for(Int_t i=AliL3Transform::GetLastRow(fPatch); i>=AliL3Transform::GetFirstRow(fPatch); i--) { for(Int_t j=0; jGetNTracks(); j++) { AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(j); if(!track) continue; if(!track->GetCrossingPoint(i,hit)) { cerr<<"AliL3Modeller::CalculateCrossingPoints : Track "<GetFirstPointX()<< " nhits "<GetNHits()<GetTgl()<<" psi "<GetPsi()<<" charge "<GetCharge()<GetCenterX()<<" "<GetCenterY()<Remove(j); continue; } //cout<<" x "<GetPointX()<<" y "<GetPointY()<<" z "<GetPointZ()<AliL3Transform::GetNPads(i) || hit[2]<0 || hit[2]>AliL3Transform::GetNTimeBins()) {//Track is leaving the patch, so flag the track hits (<0) track->SetPadHit(i,-1); track->SetTimeHit(i,-1); continue; } track->SetPadHit(i,hit[1]); track->SetTimeHit(i,hit[2]); //if(hit[1]<0 || hit[2]>445) //if(hit[2]<0 || hit[2]>445) //cout<<"pad "<GetPt()<<" psi "<GetPsi()<<" tgl "<GetTgl()<<" firstpoint "<GetFirstPointX()<<" "<GetFirstPointY()<<" "<GetFirstPointZ()<Compress(); cout<<"And there are "<GetNTracks()<<" tracks remaining"<GetNTracks(); i++) { AliL3ModelTrack *track1 = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i); if(!track1) continue; for(Int_t j=i+1; jGetNTracks(); j++) { AliL3ModelTrack *track2 = (AliL3ModelTrack*)fTracks->GetCheckedTrack(j); if(!track2) continue; for(Int_t k=AliL3Transform::GetFirstRow(fPatch); k<=AliL3Transform::GetLastRow(fPatch); k++) { if(track1->GetPadHit(k)<0 || track1->GetTimeHit(k)<0 || track2->GetPadHit(k)<0 || track2->GetTimeHit(k)<0) continue; if(fabs(track1->GetPadHit(k)-track2->GetPadHit(k)) <= fPadOverlap && fabs(track1->GetTimeHit(k)-track2->GetTimeHit(k)) <= fTimeOverlap) { track2->SetOverlap(k,i); track1->SetOverlap(k,j); } } } } cout<<"done"<fCharge; Float_t pad = (Float_t)cl->fPad/charge; Float_t time = (Float_t)cl->fTime/charge; Float_t s2 = (Float_t)cl->fSigmaY2/charge - pad*pad; //Save the sigmas in pad and time: sigmaY2 = (s2);// + 1./12);//*padw*padw; /*Constants added by offline if(s2 != 0) { sigmaY2 = sigmaY2*0.108; if(fPatch<3) sigmaY2 = sigmaY2*2.07; } */ s2 = (Float_t)cl->fSigmaZ2/charge - time*time; timew = AliL3Transform::GetZWidth(); sigmaZ2 = (s2);// +1./12);//*timew*timew; /*Constants added by offline if(s2 != 0) { sigmaZ2 = sigmaZ2*0.169; if(fPatch < 3) sigmaZ2 = sigmaZ2*1.77; } */ }