1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-------------------------------------------------------
19 // Implementation of the TPC tracker
21 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
22 //-------------------------------------------------------
23 #include <TObjArray.h>
30 #include "AliTPCtracker.h"
31 #include "AliTPCcluster.h"
32 #include "AliTPCParam.h"
33 #include "AliClusters.h"
35 ClassImp(AliTPCtracker)
37 //_____________________________________________________________________________
38 AliTPCtracker::AliTPCtracker(const AliTPCParam *par):
39 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
41 //---------------------------------------------------------------------
42 // The main TPC tracker constructor
43 //---------------------------------------------------------------------
44 fInnerSec=new AliTPCSector[fkNIS];
45 fOuterSec=new AliTPCSector[fkNOS];
48 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
49 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
51 fParam = (AliTPCParam*) par;
56 //_____________________________________________________________________________
57 AliTPCtracker::~AliTPCtracker() {
58 //------------------------------------------------------------------
59 // TPC tracker destructor
60 //------------------------------------------------------------------
69 //_____________________________________________________________________________
70 Double_t f1(Double_t x1,Double_t y1,
71 Double_t x2,Double_t y2,
72 Double_t x3,Double_t y3)
74 //-----------------------------------------------------------------
75 // Initial approximation of the track curvature
76 //-----------------------------------------------------------------
77 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
78 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
79 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
80 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
81 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
83 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
85 return -xr*yr/sqrt(xr*xr+yr*yr);
89 //_____________________________________________________________________________
90 Double_t f2(Double_t x1,Double_t y1,
91 Double_t x2,Double_t y2,
92 Double_t x3,Double_t y3)
94 //-----------------------------------------------------------------
95 // Initial approximation of the track curvature times center of curvature
96 //-----------------------------------------------------------------
97 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
98 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
99 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
100 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
101 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
103 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
105 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
108 //_____________________________________________________________________________
109 Double_t f3(Double_t x1,Double_t y1,
110 Double_t x2,Double_t y2,
111 Double_t z1,Double_t z2)
113 //-----------------------------------------------------------------
114 // Initial approximation of the tangent of the track dip angle
115 //-----------------------------------------------------------------
116 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
119 //_____________________________________________________________________________
120 Int_t AliTPCtracker::LoadClusters(TTree *cTree) {
121 //-----------------------------------------------------------------
122 // This function loads TPC clusters.
123 //-----------------------------------------------------------------
124 TBranch *branch=cTree->GetBranch("Segment");
126 Error("LoadClusters","Can't get the branch !");
130 AliClusters carray, *addr=&carray;
131 carray.SetClass("AliTPCcluster");
133 branch->SetAddress(&addr);
135 Int_t nentr=(Int_t)cTree->GetEntries();
137 for (Int_t i=0; i<nentr; i++) {
140 Int_t ncl=carray.GetArray()->GetEntriesFast();
142 Int_t nir=fInnerSec->GetNRows(), nor=fOuterSec->GetNRows();
143 Int_t id=carray.GetID();
144 if ((id<0) || (id>2*(fkNIS*nir + fkNOS*nor))) {
145 Fatal("LoadClusters","Wrong index !");
147 Int_t outindex = 2*fkNIS*nir;
150 Int_t row = id - sec*nir;
152 AliTPCRow &padrow=fInnerSec[sec][row];
154 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
155 padrow.InsertCluster(c,sec,row);
160 Int_t row = id - sec*nor;
162 AliTPCRow &padrow=fOuterSec[sec][row];
164 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
165 padrow.InsertCluster(c,sec+fkNIS,row);
168 carray.GetArray()->Clear();
174 //_____________________________________________________________________________
175 void AliTPCtracker::UnloadClusters() {
176 //-----------------------------------------------------------------
177 // This function unloads TPC clusters.
178 //-----------------------------------------------------------------
180 for (i=0; i<fkNIS; i++) {
181 Int_t nr=fInnerSec->GetNRows();
182 for (Int_t n=0; n<nr; n++) fInnerSec[i][n].ResetClusters();
184 for (i=0; i<fkNOS; i++) {
185 Int_t nr=fOuterSec->GetNRows();
186 for (Int_t n=0; n<nr; n++) fOuterSec[i][n].ResetClusters();
190 //_____________________________________________________________________________
191 Int_t AliTPCtracker::FollowProlongation(AliTPCseed& t, Int_t rf) {
192 //-----------------------------------------------------------------
193 // This function tries to find a track prolongation.
194 //-----------------------------------------------------------------
195 Double_t xt=t.GetX();
196 const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? kRowsToSkip :
197 Int_t(0.5*fSectors->GetNRows());
198 Int_t tryAgain=kSKIP;
200 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
201 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
202 if (alpha < 0. ) alpha += 2.*TMath::Pi();
203 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
205 Int_t nrows=fSectors->GetRowNumber(xt)-1;
206 for (Int_t nr=nrows; nr>=rf; nr--) {
207 Double_t x=fSectors->GetX(nr), ymax=fSectors->GetMaxY(nr);
208 if (!t.PropagateTo(x)) return 0;
212 Double_t maxchi2=kMaxCHI2;
213 const AliTPCRow &krow=fSectors[s][nr];
214 Double_t pt=t.GetConvConst()/(100/0.299792458/0.2)/t.Get1Pt();
215 Double_t sy2=AliTPCcluster::SigmaY2(t.GetX(),t.GetTgl(),pt);
216 Double_t sz2=AliTPCcluster::SigmaZ2(t.GetX(),t.GetTgl());
217 Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
220 if (t.GetNumberOfClusters()>4)
221 Warning("FindProlongation","Too broad road !");
226 for (Int_t i=krow.Find(y-road); i<krow; i++) {
227 AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
228 if (c->GetY() > y+road) break;
229 if (c->IsUsed()) continue;
230 if ((c->GetZ()-z)*(c->GetZ()-z) > 16.*(t.GetSigmaZ2()+sz2)) continue;
231 Double_t chi2=t.GetPredictedChi2(c);
232 if (chi2 > maxchi2) continue;
235 index=krow.GetIndex(i);
239 Float_t l=fSectors->GetPadPitchWidth();
240 Float_t corr=1.; if (nr>63) corr=0.67; // new (third) pad response !
241 t.SetSampledEdx(cl->GetQ()/l*corr,t.GetNumberOfClusters());
242 if (!t.Update(cl,maxchi2,index)) {
243 if (!tryAgain--) return 0;
244 } else tryAgain=kSKIP;
246 if (tryAgain==0) break;
249 if (!t.Rotate(fSectors->GetAlpha())) return 0;
250 } else if (y <-ymax) {
252 if (!t.Rotate(-fSectors->GetAlpha())) return 0;
260 //_____________________________________________________________________________
262 Int_t AliTPCtracker::FollowRefitInward(AliTPCseed *seed, AliTPCtrack *track) {
264 // This function propagates seed inward TPC using old clusters
267 // Sylwester Radomski, GSI
273 Int_t nRows = fSectors->GetNRows();
274 for (Int_t iRow = nRows-1; iRow >= 0; iRow--) {
276 Double_t x = fSectors->GetX(iRow);
277 if (!seed->PropagateTo(x)) return 0;
279 // try to find an assigned cluster in this row
281 AliTPCcluster* cluster = NULL;
284 for (Int_t iCluster = 0; iCluster < track->GetNumberOfClusters(); iCluster++) {
285 idx = track->GetClusterIndex(iCluster);
286 sec = (idx&0xff000000)>>24;
287 Int_t row = (idx&0x00ff0000)>>16;
288 if (((fSectors == fInnerSec) && (sec >= fkNIS)) ||
289 ((fSectors == fOuterSec) && (sec < fkNIS))) continue;
291 cluster = (AliTPCcluster*) GetCluster(idx);
296 // update the track seed with the found cluster
299 Double_t dAlpha = fParam->GetAngle(sec) - seed->GetAlpha();
300 if (TMath::Abs(dAlpha) > 0.0001) {
301 if (!seed->Rotate(dAlpha)) return 0;
302 if (!seed->PropagateTo(x)) return 0;
305 seed->Update(cluster, seed->GetPredictedChi2(cluster), idx);
312 //_____________________________________________________________________________
313 Int_t AliTPCtracker::FollowBackProlongation
314 (AliTPCseed& seed, const AliTPCtrack &track) {
315 //-----------------------------------------------------------------
316 // This function propagates tracks back through the TPC
317 //-----------------------------------------------------------------
318 Double_t alpha=seed.GetAlpha() - fSectors->GetAlphaShift();
319 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
320 if (alpha < 0. ) alpha += 2.*TMath::Pi();
321 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
323 Int_t idx=-1, sec=-1, row=-1;
324 Int_t nc=seed.GetNumber();
327 idx=track.GetClusterIndex(nc);
328 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
330 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
331 else { if (sec < fkNIS) row=-1; }
333 Int_t nr=fSectors->GetNRows();
334 for (Int_t i=0; i<nr; i++) {
335 Double_t x=fSectors->GetX(i), ymax=fSectors->GetMaxY(i);
336 Double_t y=seed.GetYat(x);
340 if (!seed.Rotate(fSectors->GetAlpha())) return 0;
341 } else if (y <-ymax) {
343 if (!seed.Rotate(-fSectors->GetAlpha())) return 0;
346 if (!seed.PropagateTo(x)) return 0;
350 Double_t maxchi2=kMaxCHI2;
351 Double_t pt=seed.GetConvConst()/(100/0.299792458/0.2)/seed.Get1Pt();
352 Double_t sy2=AliTPCcluster::SigmaY2(seed.GetX(),seed.GetTgl(),pt);
353 Double_t sz2=AliTPCcluster::SigmaZ2(seed.GetX(),seed.GetTgl());
354 Double_t road=4.*sqrt(seed.GetSigmaY2() + sy2), z=seed.GetZ();
356 Warning("FollowBackProlongation","Too broad road !");
360 Int_t accepted=seed.GetNumberOfClusters();
362 //try to accept already found cluster
363 AliTPCcluster *c=(AliTPCcluster*)GetCluster(idx);
365 if ((chi2=seed.GetPredictedChi2(c))<maxchi2 || accepted<27) {
366 index=idx; cl=c; maxchi2=chi2;
367 } //else cerr<<"AliTPCtracker::FollowBackProlongation: oulier !\n";
370 idx=track.GetClusterIndex(nc);
371 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
373 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
374 else { if (sec < fkNIS) row=-1; }
378 //try to fill the gap
379 const AliTPCRow &krow=fSectors[s][i];
382 for (Int_t i=krow.Find(y-road); i<krow; i++) {
383 AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
384 if (c->GetY() > y+road) break;
385 if (c->IsUsed()) continue;
386 if ((c->GetZ()-z)*(c->GetZ()-z)>16.*(seed.GetSigmaZ2()+sz2)) continue;
387 Double_t chi2=seed.GetPredictedChi2(c);
388 if (chi2 > maxchi2) continue;
391 index=krow.GetIndex(i);
397 Float_t l=fSectors->GetPadPitchWidth();
398 Float_t corr=1.; if (i>63) corr=0.67; // new (third) pad response !
399 seed.SetSampledEdx(cl->GetQ()/l*corr,seed.GetNumberOfClusters());
400 seed.Update(cl,maxchi2,index);
410 //_____________________________________________________________________________
411 void AliTPCtracker::MakeSeeds(Int_t i1, Int_t i2) {
412 //-----------------------------------------------------------------
413 // This function creates track seeds.
414 //-----------------------------------------------------------------
415 Double_t x[5], c[15];
417 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
418 Double_t cs=cos(alpha), sn=sin(alpha);
420 Double_t x1 =fSectors->GetX(i1);
421 Double_t xx2=fSectors->GetX(i2);
423 for (Int_t ns=0; ns<fN; ns++) {
424 Int_t nl=fSectors[(ns-1+fN)%fN][i2];
425 Int_t nm=fSectors[ns][i2];
426 Int_t nu=fSectors[(ns+1)%fN][i2];
427 const AliTPCRow& kr1=fSectors[ns][i1];
428 for (Int_t is=0; is < kr1; is++) {
429 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
430 for (Int_t js=0; js < nl+nm+nu; js++) {
431 const AliTPCcluster *kcl;
433 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
436 const AliTPCRow& kr2=fSectors[(ns-1+fN)%fN][i2];
438 y2=kcl->GetY(); z2=kcl->GetZ();
443 const AliTPCRow& kr2=fSectors[ns][i2];
445 x2=xx2; y2=kcl->GetY(); z2=kcl->GetZ();
447 const AliTPCRow& kr2=fSectors[(ns+1)%fN][i2];
449 y2=kcl->GetY(); z2=kcl->GetZ();
454 Double_t zz=z1 - (z1-z3)/(x1-x3)*(x1-x2);
455 if (TMath::Abs(zz-z2)>5.) continue;
457 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
459 Warning("MakeSeeds","Straight seed !");
464 x[4]=f1(x1,y1,x2,y2,x3,y3);
465 if (TMath::Abs(x[4]) >= 0.0066) continue;
466 x[2]=f2(x1,y1,x2,y2,x3,y3);
467 //if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
468 x[3]=f3(x1,y1,x2,y2,z1,z2);
469 if (TMath::Abs(x[3]) > 1.2) continue;
470 Double_t a=asin(x[2]);
471 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
472 if (TMath::Abs(zv-z3)>10.) continue;
474 Double_t sy1=kr1[is]->GetSigmaY2(), sz1=kr1[is]->GetSigmaZ2();
475 Double_t sy2=kcl->GetSigmaY2(), sz2=kcl->GetSigmaZ2();
476 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
477 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
479 Double_t f40=(f1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
480 Double_t f42=(f1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
481 Double_t f43=(f1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
482 Double_t f20=(f2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
483 Double_t f22=(f2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
484 Double_t f23=(f2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
485 Double_t f30=(f3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
486 Double_t f31=(f3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
487 Double_t f32=(f3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
488 Double_t f34=(f3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
492 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
493 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
494 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
495 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
496 c[13]=f30*sy1*f40+f32*sy2*f42;
497 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
499 UInt_t index=kr1.GetIndex(is);
500 AliTPCseed *track=new AliTPCseed(index, x, c, x1, ns*alpha+shift);
501 Float_t l=fSectors->GetPadPitchWidth();
502 track->SetSampledEdx(kr1[is]->GetQ()/l,0);
504 Int_t rc=FollowProlongation(*track, i2);
505 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
506 else fSeeds->AddLast(track);
512 //_____________________________________________________________________________
513 Int_t AliTPCtracker::ReadSeeds(const TFile *inp) {
514 //-----------------------------------------------------------------
515 // This function reades track seeds.
516 //-----------------------------------------------------------------
517 TDirectory *savedir=gDirectory;
519 TFile *in=(TFile*)inp;
521 Error("ReadSeeds","Input file has not been open !");
526 TTree *seedTree=(TTree*)in->Get("Seeds");
528 Error("ReadSeeds","Can't get a tree with track seeds !");
531 AliTPCtrack *seed=new AliTPCtrack;
532 seedTree->SetBranchAddress("tracks",&seed);
534 if (fSeeds==0) fSeeds=new TObjArray(15000);
536 Int_t n=(Int_t)seedTree->GetEntries();
537 for (Int_t i=0; i<n; i++) {
538 seedTree->GetEvent(i);
539 seed->ResetClusters();
540 fSeeds->AddLast(new AliTPCseed(*seed));
545 delete seedTree; //Thanks to Mariana Bondila
552 //_____________________________________________________________________________
553 Int_t AliTPCtracker::Clusters2Tracks(AliESD *event) {
554 //-----------------------------------------------------------------
555 // This is a track finder.
556 // The clusters must be already loaded !
557 //-----------------------------------------------------------------
560 Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows();
561 Int_t nrows=nlow+nup;
563 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
564 fSectors=fOuterSec; fN=fkNOS;
565 fSeeds=new TObjArray(15000);
566 MakeSeeds(nup-1, nup-1-gap);
567 MakeSeeds(nup-1-shift, nup-1-shift-gap);
571 Int_t nseed=fSeeds->GetEntriesFast();
572 for (Int_t i=0; i<nseed; i++) {
573 //tracking in the outer sectors
574 fSectors=fOuterSec; fN=fkNOS;
576 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
577 if (!FollowProlongation(t)) {
578 delete fSeeds->RemoveAt(i);
582 //tracking in the inner sectors
583 fSectors=fInnerSec; fN=fkNIS;
585 Double_t alpha=t.GetAlpha() - fInnerSec->GetAlphaShift();
586 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
587 if (alpha < 0. ) alpha += 2.*TMath::Pi();
588 Int_t ns=Int_t(alpha/fInnerSec->GetAlpha())%fkNIS;
590 alpha=ns*fInnerSec->GetAlpha()+fInnerSec->GetAlphaShift()-t.GetAlpha();
592 if (t.Rotate(alpha)) {
593 if (FollowProlongation(t)) {
594 if (t.GetNumberOfClusters() >= Int_t(0.4*nrows)) {
596 CookLabel(pt,0.1); //For comparison only
597 pt->PropagateTo(fParam->GetInnerRadiusLow());
599 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
601 event->AddTrack(&iotrack);
607 delete fSeeds->RemoveAt(i);
610 Info("Clusters2Tracks","Number of found tracks : %d",
611 event->GetNumberOfTracks());
613 fSeeds->Clear(); delete fSeeds; fSeeds=0;
618 //_____________________________________________________________________________
619 Int_t AliTPCtracker::RefitInward(AliESD* event) {
621 // The function propagates tracks throught TPC inward
622 // using already associated clusters.
623 // The clusters must be already loaded !
626 Int_t nTracks = event->GetNumberOfTracks();
629 for (Int_t i = 0; i < nTracks; i++) {
630 AliESDtrack* track = event->GetTrack(i);
631 ULong_t status = track->GetStatus();
633 if ( (status & AliESDtrack::kTPCrefit) != 0 ) continue;
634 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
636 if ( (status & AliESDtrack::kTRDout ) != 0 )
637 if ( (status & AliESDtrack::kTRDrefit ) == 0 ) continue;
639 AliTPCtrack* tpcTrack = new AliTPCtrack(*track);
640 AliTPCseed* seed=new AliTPCseed(*tpcTrack); seed->ResetClusters();
642 if ( (status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance();
644 fSectors = fOuterSec;
646 Int_t res = FollowRefitInward(seed, tpcTrack);
648 Int_t nc = seed->GetNumberOfClusters();
650 fSectors = fInnerSec;
652 res = FollowRefitInward(seed, tpcTrack);
653 UseClusters(seed, nc);
656 seed->PropagateTo(fParam->GetInnerRadiusLow());
657 seed->SetLabel(tpcTrack->GetLabel());
658 seed->SetdEdx(tpcTrack->GetdEdx());
659 track->UpdateTrackParams(seed, AliESDtrack::kTPCrefit);
667 Info("RefitInward","Number of refitted tracks : %d",nRefited);
672 Int_t AliTPCtracker::PropagateBack(AliESD *event) {
673 //-----------------------------------------------------------------
674 // This function propagates tracks back through the TPC.
675 // The clusters must be already loaded !
676 //-----------------------------------------------------------------
677 Int_t nentr=event->GetNumberOfTracks();
678 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
681 for (Int_t i=0; i<nentr; i++) {
682 AliESDtrack *esd=event->GetTrack(i);
683 ULong_t status=esd->GetStatus();
685 if ( (status & AliESDtrack::kTPCin ) == 0 ) continue;
686 if ( (status & AliESDtrack::kTPCout) != 0 ) continue;
687 if ( (status & AliESDtrack::kITSin) != 0 )
688 if ( (status & AliESDtrack::kITSout) == 0 ) continue;
690 const AliTPCtrack t(*esd);
691 AliTPCseed s(t); s.ResetClusters();
693 if ( (status & AliESDtrack::kITSout) == 0 ) s.ResetCovariance();
695 Int_t nc=t.GetNumberOfClusters();
696 s.SetNumber(nc); //set number of the cluster to start with
699 fSectors=fInnerSec; fN=fkNIS;
701 Double_t alpha=s.GetAlpha() - fSectors->GetAlphaShift();
702 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
703 if (alpha < 0. ) alpha += 2.*TMath::Pi();
704 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
705 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
708 if (!s.Rotate(alpha)) continue;
709 if (!FollowBackProlongation(s,t)) continue;
714 fSectors=fOuterSec; fN=fkNOS;
716 nc=s.GetNumberOfClusters();
718 alpha=s.GetAlpha() - fSectors->GetAlphaShift();
719 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
720 if (alpha < 0. ) alpha += 2.*TMath::Pi();
721 ns=Int_t(alpha/fSectors->GetAlpha())%fN;
723 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
726 if (!s.Rotate(alpha)) continue;
727 if (!FollowBackProlongation(s,t)) continue;
729 Int_t nrows=fOuterSec->GetNRows()+fInnerSec->GetNRows();
730 if (s.GetNumberOfClusters() < Int_t(0.4*nrows)) continue;
732 s.PropagateTo(fParam->GetOuterRadiusUp());
734 CookLabel(&s,0.1); //For comparison only
736 esd->UpdateTrackParams(&s,AliESDtrack::kTPCout);
739 Info("PropagateBack","Number of back propagated tracks: %d",ntrk);
744 //_________________________________________________________________________
745 AliCluster *AliTPCtracker::GetCluster(Int_t index) const {
746 //--------------------------------------------------------------------
747 // Return pointer to a given cluster
748 //--------------------------------------------------------------------
749 Int_t sec=(index&0xff000000)>>24;
750 Int_t row=(index&0x00ff0000)>>16;
751 Int_t ncl=(index&0x0000ffff)>>00;
753 const AliTPCcluster *cl=0;
756 cl=fInnerSec[sec][row].GetUnsortedCluster(ncl);
759 cl=fOuterSec[sec][row].GetUnsortedCluster(ncl);
762 return (AliCluster*)cl;
765 //__________________________________________________________________________
766 void AliTPCtracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const {
767 //--------------------------------------------------------------------
768 //This function "cooks" a track label. If label<0, this track is fake.
769 //--------------------------------------------------------------------
770 Int_t noc=t->GetNumberOfClusters();
771 Int_t *lb=new Int_t[noc];
772 Int_t *mx=new Int_t[noc];
773 AliCluster **clusters=new AliCluster*[noc];
776 for (i=0; i<noc; i++) {
778 Int_t index=t->GetClusterIndex(i);
779 clusters[i]=GetCluster(index);
783 for (i=0; i<noc; i++) {
784 AliCluster *c=clusters[i];
785 lab=TMath::Abs(c->GetLabel(0));
787 for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
793 for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
795 for (i=0; i<noc; i++) {
796 AliCluster *c=clusters[i];
797 if (TMath::Abs(c->GetLabel(1)) == lab ||
798 TMath::Abs(c->GetLabel(2)) == lab ) max++;
801 if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
804 Int_t tail=Int_t(0.10*noc);
806 for (i=1; i<=tail; i++) {
807 AliCluster *c=clusters[noc-i];
808 if (lab == TMath::Abs(c->GetLabel(0)) ||
809 lab == TMath::Abs(c->GetLabel(1)) ||
810 lab == TMath::Abs(c->GetLabel(2))) max++;
812 if (max < Int_t(0.5*tail)) lab=-lab;
822 //_________________________________________________________________________
823 void AliTPCtracker::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
824 //-----------------------------------------------------------------------
825 // Setup inner sector
826 //-----------------------------------------------------------------------
828 fAlpha=par->GetInnerAngle();
829 fAlphaShift=par->GetInnerAngleShift();
830 fPadPitchWidth=par->GetInnerPadPitchWidth();
831 f1PadPitchLength=par->GetInnerPadPitchLength();
832 f2PadPitchLength=f1PadPitchLength;
834 fN=par->GetNRowLow();
835 fRow=new AliTPCRow[fN];
836 for (Int_t i=0; i<fN; i++) fRow[i].SetX(par->GetPadRowRadiiLow(i));
838 fAlpha=par->GetOuterAngle();
839 fAlphaShift=par->GetOuterAngleShift();
840 fPadPitchWidth=par->GetOuterPadPitchWidth();
841 f1PadPitchLength=par->GetOuter1PadPitchLength();
842 f2PadPitchLength=par->GetOuter2PadPitchLength();
845 fRow=new AliTPCRow[fN];
846 for (Int_t i=0; i<fN; i++){
847 fRow[i].SetX(par->GetPadRowRadiiUp(i));
852 //_________________________________________________________________________
854 AliTPCRow::InsertCluster(const AliTPCcluster* c, Int_t sec, Int_t row) {
855 //-----------------------------------------------------------------------
856 // Insert a cluster into this pad row in accordence with its y-coordinate
857 //-----------------------------------------------------------------------
858 if (fN==kMaxClusterPerRow) {
859 ::Error("InsertCluster","Too many clusters !");
863 Int_t index=(((sec<<8)+row)<<16)+fN;
866 fSize=kMaxClusterPerRow/8;
867 fClusterArray=new AliTPCcluster[fSize];
869 fClusterArray[0]=*c; fClusters[fN++]=fClusterArray;
875 AliTPCcluster *buff=new AliTPCcluster[size];
876 memcpy(buff,fClusterArray,fSize*sizeof(AliTPCcluster));
877 for (Int_t i=0; i<fN; i++)
878 fClusters[i]=buff+(fClusters[i]-fClusterArray);
879 delete[] fClusterArray;
884 Int_t i=Find(c->GetY());
885 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCcluster*));
886 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
888 fClusters[i]=fClusterArray+fN; fClusterArray[fN++]=*c;
891 //___________________________________________________________________
892 Int_t AliTPCtracker::AliTPCRow::Find(Double_t y) const {
893 //-----------------------------------------------------------------------
894 // Return the index of the nearest cluster
895 //-----------------------------------------------------------------------
897 if (y <= fClusters[0]->GetY()) return 0;
898 if (y > fClusters[fN-1]->GetY()) return fN;
899 Int_t b=0, e=fN-1, m=(b+e)/2;
900 for (; b<e; m=(b+e)/2) {
901 if (y > fClusters[m]->GetY()) b=m+1;
907 //_____________________________________________________________________________
908 void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
909 //-----------------------------------------------------------------
910 // This funtion calculates dE/dX within the "low" and "up" cuts.
911 //-----------------------------------------------------------------
913 Int_t nc=GetNumberOfClusters();
915 Int_t swap;//stupid sorting
918 for (i=0; i<nc-1; i++) {
919 if (fdEdxSample[i]<=fdEdxSample[i+1]) continue;
920 Float_t tmp=fdEdxSample[i];
921 fdEdxSample[i]=fdEdxSample[i+1]; fdEdxSample[i+1]=tmp;
926 Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
928 for (i=nl; i<=nu; i++) dedx += fdEdxSample[i];