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 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
36 extern Double_t SigmaZ2(Double_t, Double_t);
38 ClassImp(AliTPCtracker)
40 //_____________________________________________________________________________
41 AliTPCtracker::AliTPCtracker(const AliTPCParam *par):
42 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
44 //---------------------------------------------------------------------
45 // The main TPC tracker constructor
46 //---------------------------------------------------------------------
47 fInnerSec=new AliTPCSector[fkNIS];
48 fOuterSec=new AliTPCSector[fkNOS];
51 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
52 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
54 fParam = (AliTPCParam*) par;
66 //_____________________________________________________________________________
67 AliTPCtracker::~AliTPCtracker() {
68 //------------------------------------------------------------------
69 // TPC tracker destructor
70 //------------------------------------------------------------------
85 //_____________________________________________________________________________
87 void AliTPCtracker::SetBarrelTree(const char *mode) {
89 // Creates a tree for BarrelTracks
90 // mode = "back" or "refit"
95 if (!IsStoringBarrel()) return;
97 TDirectory *sav = gDirectory;
98 if (!fBarrelFile) fBarrelFile = new TFile("AliBarrelTracks.root", "UPDATE");
101 sprintf(buff, "BarrelTPC_%d_%s", GetEventNumber(), mode);
104 fBarrelTree = new TTree(buff, "Barrel TPC tracks");
106 if (!fBarrelArray) fBarrelArray = new TClonesArray("AliBarrelTrack", 4);
107 for(Int_t i=0; i<4; i++) new((*fBarrelArray)[i]) AliBarrelTrack();
109 fBarrelTree->Branch("tracks", &fBarrelArray);
113 //_____________________________________________________________________________
115 void AliTPCtracker::StoreBarrelTrack(AliTPCtrack *ps, Int_t refPlane, Int_t isIn) {
117 // Stores Track at a given reference plane
120 // isIn: 1 - backward, 2 - refit
123 if (!IsStoringBarrel()) return;
124 if (refPlane < 0 || refPlane > 4) return;
125 if (isIn > 2) return;
127 static Int_t nClusters;
129 static Double_t chi2;
132 Int_t newClusters, newWrong;
135 if ( (refPlane == 1 && isIn == kTrackBack) ||
136 (refPlane == 4 && isIn == kTrackRefit) ) {
138 fBarrelArray->Clear();
139 nClusters = nWrong = 0;
146 if (refPlane == 1) refX = fParam->GetInnerRadiusLow();
147 if (refPlane == 2) refX = fParam->GetInnerRadiusUp();
148 if (refPlane == 3) refX = fParam->GetOuterRadiusLow();
149 if (refPlane == 4) refX = fParam->GetOuterRadiusUp();
151 ps->PropagateTo(refX);
153 fBarrelTrack = (AliBarrelTrack*)(*fBarrelArray)[index++];
154 ps->GetBarrelTrack(fBarrelTrack);
156 newClusters = ps->GetNumberOfClusters() - nClusters;
157 newWrong = ps->GetNWrong() - nWrong;
158 newChi2 = ps->GetChi2() - chi2;
160 nClusters = ps->GetNumberOfClusters();
161 nWrong = ps->GetNWrong();
162 chi2 = ps->GetChi2();
164 fBarrelTrack->SetNClusters(newClusters, newChi2);
165 fBarrelTrack->SetNWrongClusters(newWrong);
166 fBarrelTrack->SetRefPlane(refPlane, isIn);
170 //_____________________________________________________________________________
171 Double_t f1(Double_t x1,Double_t y1,
172 Double_t x2,Double_t y2,
173 Double_t x3,Double_t y3)
175 //-----------------------------------------------------------------
176 // Initial approximation of the track curvature
177 //-----------------------------------------------------------------
178 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
179 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
180 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
181 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
182 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
184 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
186 return -xr*yr/sqrt(xr*xr+yr*yr);
190 //_____________________________________________________________________________
191 Double_t f2(Double_t x1,Double_t y1,
192 Double_t x2,Double_t y2,
193 Double_t x3,Double_t y3)
195 //-----------------------------------------------------------------
196 // Initial approximation of the track curvature times center of curvature
197 //-----------------------------------------------------------------
198 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
199 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
200 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
201 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
202 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
204 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
206 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
209 //_____________________________________________________________________________
210 Double_t f3(Double_t x1,Double_t y1,
211 Double_t x2,Double_t y2,
212 Double_t z1,Double_t z2)
214 //-----------------------------------------------------------------
215 // Initial approximation of the tangent of the track dip angle
216 //-----------------------------------------------------------------
217 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
220 //_____________________________________________________________________________
221 Int_t AliTPCtracker::LoadClusters(TTree *cTree) {
222 //-----------------------------------------------------------------
223 // This function loads TPC clusters.
224 //-----------------------------------------------------------------
225 TBranch *branch=cTree->GetBranch("Segment");
227 Error("LoadClusters","Can't get the branch !");
231 AliClusters carray, *addr=&carray;
232 carray.SetClass("AliTPCcluster");
234 branch->SetAddress(&addr);
236 Int_t nentr=(Int_t)cTree->GetEntries();
238 for (Int_t i=0; i<nentr; i++) {
241 Int_t ncl=carray.GetArray()->GetEntriesFast();
243 Int_t nir=fInnerSec->GetNRows(), nor=fOuterSec->GetNRows();
244 Int_t id=carray.GetID();
245 if ((id<0) || (id>2*(fkNIS*nir + fkNOS*nor))) {
246 Fatal("LoadClusters","Wrong index !");
248 Int_t outindex = 2*fkNIS*nir;
251 Int_t row = id - sec*nir;
253 AliTPCRow &padrow=fInnerSec[sec][row];
255 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
256 padrow.InsertCluster(c,sec,row);
261 Int_t row = id - sec*nor;
263 AliTPCRow &padrow=fOuterSec[sec][row];
265 AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
266 padrow.InsertCluster(c,sec+fkNIS,row);
269 carray.GetArray()->Clear();
275 //_____________________________________________________________________________
276 void AliTPCtracker::UnloadClusters() {
277 //-----------------------------------------------------------------
278 // This function unloads TPC clusters.
279 //-----------------------------------------------------------------
281 for (i=0; i<fkNIS; i++) {
282 Int_t nr=fInnerSec->GetNRows();
283 for (Int_t n=0; n<nr; n++) fInnerSec[i][n].ResetClusters();
285 for (i=0; i<fkNOS; i++) {
286 Int_t nr=fOuterSec->GetNRows();
287 for (Int_t n=0; n<nr; n++) fOuterSec[i][n].ResetClusters();
291 //_____________________________________________________________________________
292 Int_t AliTPCtracker::FollowProlongation(AliTPCseed& t, Int_t rf) {
293 //-----------------------------------------------------------------
294 // This function tries to find a track prolongation.
295 //-----------------------------------------------------------------
296 Double_t xt=t.GetX();
297 const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? kRowsToSkip :
298 Int_t(0.5*fSectors->GetNRows());
299 Int_t tryAgain=kSKIP;
301 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
302 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
303 if (alpha < 0. ) alpha += 2.*TMath::Pi();
304 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
306 Int_t nrows=fSectors->GetRowNumber(xt)-1;
307 for (Int_t nr=nrows; nr>=rf; nr--) {
308 Double_t x=fSectors->GetX(nr), ymax=fSectors->GetMaxY(nr);
309 if (!t.PropagateTo(x)) return 0;
313 Double_t maxchi2=kMaxCHI2;
314 const AliTPCRow &krow=fSectors[s][nr];
315 Double_t pt=t.GetConvConst()/(100/0.299792458/0.2)/t.Get1Pt();
316 Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),pt);
317 Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl());
318 Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
321 if (t.GetNumberOfClusters()>4)
322 Warning("FindProlongation","Too broad road !");
327 for (Int_t i=krow.Find(y-road); i<krow; i++) {
328 AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
329 if (c->GetY() > y+road) break;
330 if (c->IsUsed()) continue;
331 if ((c->GetZ()-z)*(c->GetZ()-z) > 16.*(t.GetSigmaZ2()+sz2)) continue;
332 Double_t chi2=t.GetPredictedChi2(c);
333 if (chi2 > maxchi2) continue;
336 index=krow.GetIndex(i);
340 Float_t l=fSectors->GetPadPitchWidth();
341 Float_t corr=1.; if (nr>63) corr=0.67; // new (third) pad response !
342 t.SetSampledEdx(cl->GetQ()/l*corr,t.GetNumberOfClusters());
343 if (!t.Update(cl,maxchi2,index)) {
344 if (!tryAgain--) return 0;
345 } else tryAgain=kSKIP;
347 if (tryAgain==0) break;
350 if (!t.Rotate(fSectors->GetAlpha())) return 0;
351 } else if (y <-ymax) {
353 if (!t.Rotate(-fSectors->GetAlpha())) return 0;
361 //_____________________________________________________________________________
363 Int_t AliTPCtracker::FollowRefitInward(AliTPCseed *seed, AliTPCtrack *track) {
365 // This function propagates seed inward TPC using old clusters
368 // Sylwester Radomski, GSI
374 Int_t nRows = fSectors->GetNRows();
375 for (Int_t iRow = nRows-1; iRow >= 0; iRow--) {
377 Double_t x = fSectors->GetX(iRow);
378 if (!seed->PropagateTo(x)) return 0;
380 // try to find an assigned cluster in this row
382 AliTPCcluster* cluster = NULL;
385 for (Int_t iCluster = 0; iCluster < track->GetNumberOfClusters(); iCluster++) {
386 idx = track->GetClusterIndex(iCluster);
387 sec = (idx&0xff000000)>>24;
388 Int_t row = (idx&0x00ff0000)>>16;
389 if (((fSectors == fInnerSec) && (sec >= fkNIS)) ||
390 ((fSectors == fOuterSec) && (sec < fkNIS))) continue;
392 cluster = (AliTPCcluster*) GetCluster(idx);
397 // update the track seed with the found cluster
400 Double_t dAlpha = fParam->GetAngle(sec) - seed->GetAlpha();
401 if (TMath::Abs(dAlpha) > 0.0001) {
402 if (!seed->Rotate(dAlpha)) return 0;
403 if (!seed->PropagateTo(x)) return 0;
406 seed->Update(cluster, seed->GetPredictedChi2(cluster), idx);
413 //_____________________________________________________________________________
414 Int_t AliTPCtracker::FollowBackProlongation
415 (AliTPCseed& seed, const AliTPCtrack &track) {
416 //-----------------------------------------------------------------
417 // This function propagates tracks back through the TPC
418 //-----------------------------------------------------------------
419 Double_t alpha=seed.GetAlpha() - fSectors->GetAlphaShift();
420 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
421 if (alpha < 0. ) alpha += 2.*TMath::Pi();
422 Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
424 Int_t idx=-1, sec=-1, row=-1;
425 Int_t nc=seed.GetNumber();
428 idx=track.GetClusterIndex(nc);
429 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
431 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
432 else { if (sec < fkNIS) row=-1; }
434 Int_t nr=fSectors->GetNRows();
435 for (Int_t i=0; i<nr; i++) {
436 Double_t x=fSectors->GetX(i), ymax=fSectors->GetMaxY(i);
437 Double_t y=seed.GetYat(x);
441 if (!seed.Rotate(fSectors->GetAlpha())) return 0;
442 } else if (y <-ymax) {
444 if (!seed.Rotate(-fSectors->GetAlpha())) return 0;
447 if (!seed.PropagateTo(x)) return 0;
451 Double_t maxchi2=kMaxCHI2;
452 Double_t pt=seed.GetConvConst()/(100/0.299792458/0.2)/seed.Get1Pt();
453 Double_t sy2=SigmaY2(seed.GetX(),seed.GetTgl(),pt);
454 Double_t sz2=SigmaZ2(seed.GetX(),seed.GetTgl());
455 Double_t road=4.*sqrt(seed.GetSigmaY2() + sy2), z=seed.GetZ();
457 Warning("FollowBackProlongation","Too broad road !");
461 Int_t accepted=seed.GetNumberOfClusters();
463 //try to accept already found cluster
464 AliTPCcluster *c=(AliTPCcluster*)GetCluster(idx);
466 if ((chi2=seed.GetPredictedChi2(c))<maxchi2 || accepted<27) {
467 index=idx; cl=c; maxchi2=chi2;
468 } //else cerr<<"AliTPCtracker::FollowBackProlongation: oulier !\n";
471 idx=track.GetClusterIndex(nc);
472 sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
474 if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
475 else { if (sec < fkNIS) row=-1; }
479 //try to fill the gap
480 const AliTPCRow &krow=fSectors[s][i];
483 for (Int_t i=krow.Find(y-road); i<krow; i++) {
484 AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
485 if (c->GetY() > y+road) break;
486 if (c->IsUsed()) continue;
487 if ((c->GetZ()-z)*(c->GetZ()-z)>16.*(seed.GetSigmaZ2()+sz2)) continue;
488 Double_t chi2=seed.GetPredictedChi2(c);
489 if (chi2 > maxchi2) continue;
492 index=krow.GetIndex(i);
498 Float_t l=fSectors->GetPadPitchWidth();
499 Float_t corr=1.; if (i>63) corr=0.67; // new (third) pad response !
500 seed.SetSampledEdx(cl->GetQ()/l*corr,seed.GetNumberOfClusters());
501 seed.Update(cl,maxchi2,index);
511 //_____________________________________________________________________________
512 void AliTPCtracker::MakeSeeds(Int_t i1, Int_t i2) {
513 //-----------------------------------------------------------------
514 // This function creates track seeds.
515 //-----------------------------------------------------------------
516 Double_t x[5], c[15];
518 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
519 Double_t cs=cos(alpha), sn=sin(alpha);
521 Double_t x1 =fSectors->GetX(i1);
522 Double_t xx2=fSectors->GetX(i2);
524 for (Int_t ns=0; ns<fN; ns++) {
525 Int_t nl=fSectors[(ns-1+fN)%fN][i2];
526 Int_t nm=fSectors[ns][i2];
527 Int_t nu=fSectors[(ns+1)%fN][i2];
528 const AliTPCRow& kr1=fSectors[ns][i1];
529 for (Int_t is=0; is < kr1; is++) {
530 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
531 for (Int_t js=0; js < nl+nm+nu; js++) {
532 const AliTPCcluster *kcl;
534 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
537 const AliTPCRow& kr2=fSectors[(ns-1+fN)%fN][i2];
539 y2=kcl->GetY(); z2=kcl->GetZ();
544 const AliTPCRow& kr2=fSectors[ns][i2];
546 x2=xx2; y2=kcl->GetY(); z2=kcl->GetZ();
548 const AliTPCRow& kr2=fSectors[(ns+1)%fN][i2];
550 y2=kcl->GetY(); z2=kcl->GetZ();
555 Double_t zz=z1 - (z1-z3)/(x1-x3)*(x1-x2);
556 if (TMath::Abs(zz-z2)>5.) continue;
558 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
560 Warning("MakeSeeds","Straight seed !");
565 x[4]=f1(x1,y1,x2,y2,x3,y3);
566 if (TMath::Abs(x[4]) >= 0.0066) continue;
567 x[2]=f2(x1,y1,x2,y2,x3,y3);
568 //if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
569 x[3]=f3(x1,y1,x2,y2,z1,z2);
570 if (TMath::Abs(x[3]) > 1.2) continue;
571 Double_t a=asin(x[2]);
572 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
573 if (TMath::Abs(zv-z3)>10.) continue;
575 Double_t sy1=kr1[is]->GetSigmaY2(), sz1=kr1[is]->GetSigmaZ2();
576 Double_t sy2=kcl->GetSigmaY2(), sz2=kcl->GetSigmaZ2();
577 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
578 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
580 Double_t f40=(f1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
581 Double_t f42=(f1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
582 Double_t f43=(f1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
583 Double_t f20=(f2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
584 Double_t f22=(f2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
585 Double_t f23=(f2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
586 Double_t f30=(f3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
587 Double_t f31=(f3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
588 Double_t f32=(f3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
589 Double_t f34=(f3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
593 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
594 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
595 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
596 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
597 c[13]=f30*sy1*f40+f32*sy2*f42;
598 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
600 UInt_t index=kr1.GetIndex(is);
601 AliTPCseed *track=new AliTPCseed(index, x, c, x1, ns*alpha+shift);
602 Float_t l=fSectors->GetPadPitchWidth();
603 track->SetSampledEdx(kr1[is]->GetQ()/l,0);
605 Int_t rc=FollowProlongation(*track, i2);
606 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
607 else fSeeds->AddLast(track);
613 //_____________________________________________________________________________
614 Int_t AliTPCtracker::ReadSeeds(const TFile *inp) {
615 //-----------------------------------------------------------------
616 // This function reades track seeds.
617 //-----------------------------------------------------------------
618 TDirectory *savedir=gDirectory;
620 TFile *in=(TFile*)inp;
622 Error("ReadSeeds","Input file has not been open !");
627 TTree *seedTree=(TTree*)in->Get("Seeds");
629 Error("ReadSeeds","Can't get a tree with track seeds !");
632 AliTPCtrack *seed=new AliTPCtrack;
633 seedTree->SetBranchAddress("tracks",&seed);
635 if (fSeeds==0) fSeeds=new TObjArray(15000);
637 Int_t n=(Int_t)seedTree->GetEntries();
638 for (Int_t i=0; i<n; i++) {
639 seedTree->GetEvent(i);
640 fSeeds->AddLast(new AliTPCseed(*seed,seed->GetAlpha()));
645 delete seedTree; //Thanks to Mariana Bondila
652 //_____________________________________________________________________________
653 Int_t AliTPCtracker::Clusters2Tracks(AliESD *event) {
654 //-----------------------------------------------------------------
655 // This is a track finder.
656 // The clusters must be already loaded !
657 //-----------------------------------------------------------------
660 Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows();
661 Int_t nrows=nlow+nup;
663 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
664 fSectors=fOuterSec; fN=fkNOS;
665 fSeeds=new TObjArray(15000);
666 MakeSeeds(nup-1, nup-1-gap);
667 MakeSeeds(nup-1-shift, nup-1-shift-gap);
671 Int_t nseed=fSeeds->GetEntriesFast();
672 for (Int_t i=0; i<nseed; i++) {
673 //tracking in the outer sectors
674 fSectors=fOuterSec; fN=fkNOS;
676 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
677 if (!FollowProlongation(t)) {
678 delete fSeeds->RemoveAt(i);
682 //tracking in the inner sectors
683 fSectors=fInnerSec; fN=fkNIS;
685 Double_t alpha=t.GetAlpha() - fInnerSec->GetAlphaShift();
686 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
687 if (alpha < 0. ) alpha += 2.*TMath::Pi();
688 Int_t ns=Int_t(alpha/fInnerSec->GetAlpha())%fkNIS;
690 alpha=ns*fInnerSec->GetAlpha()+fInnerSec->GetAlphaShift()-t.GetAlpha();
692 if (t.Rotate(alpha)) {
693 if (FollowProlongation(t)) {
694 if (t.GetNumberOfClusters() >= Int_t(0.4*nrows)) {
696 CookLabel(pt,0.1); //For comparison only
697 pt->PropagateTo(fParam->GetInnerRadiusLow());
699 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
701 event->AddTrack(&iotrack);
707 delete fSeeds->RemoveAt(i);
710 Info("Clusters2Tracks","Number of found tracks : %d",
711 event->GetNumberOfTracks());
713 fSeeds->Clear(); delete fSeeds; fSeeds=0;
718 //_____________________________________________________________________________
719 Int_t AliTPCtracker::RefitInward(AliESD* event) {
721 // The function propagates tracks throught TPC inward
722 // using already associated clusters.
723 // The clusters must be already loaded !
726 Int_t nTracks = event->GetNumberOfTracks();
729 for (Int_t i = 0; i < nTracks; i++) {
730 AliESDtrack* track = event->GetTrack(i);
731 ULong_t status = track->GetStatus();
733 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
734 if ( (status & AliESDtrack::kTPCrefit) != 0 ) continue;
736 AliTPCtrack* tpcTrack = new AliTPCtrack(*track);
737 AliTPCseed* seed = new AliTPCseed(*tpcTrack, tpcTrack->GetAlpha());
739 fSectors = fOuterSec;
741 Int_t res = FollowRefitInward(seed, tpcTrack);
743 Int_t nc = seed->GetNumberOfClusters();
745 fSectors = fInnerSec;
747 res = FollowRefitInward(seed, tpcTrack);
748 UseClusters(seed, nc);
751 seed->PropagateTo(fParam->GetInnerRadiusLow());
752 seed->SetLabel(tpcTrack->GetLabel());
753 seed->SetdEdx(tpcTrack->GetdEdx());
754 track->UpdateTrackParams(seed, AliESDtrack::kTPCrefit);
762 Info("RefitInward","Number of refitted tracks : %d",nRefited);
767 Int_t AliTPCtracker::PropagateBack(AliESD *event) {
768 //-----------------------------------------------------------------
769 // This function propagates tracks back through the TPC.
770 // The clusters must be already loaded !
771 //-----------------------------------------------------------------
772 Int_t nentr=event->GetNumberOfTracks();
773 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
776 for (Int_t i=0; i<nentr; i++) {
777 AliESDtrack *esd=event->GetTrack(i);
778 ULong_t status=esd->GetStatus();
780 if ( (status & AliESDtrack::kTPCin ) == 0 ) continue;
781 if ( (status & AliESDtrack::kTPCout) != 0 ) continue;
783 const AliTPCtrack t(*esd);
784 AliTPCseed s(t,t.GetAlpha());
786 if ( (status & AliESDtrack::kITSout) == 0 ) s.ResetCovariance();
791 Int_t nc=t.GetNumberOfClusters();
792 s.SetNumber(nc); //set number of the cluster to start with
795 fSectors=fInnerSec; fN=fkNIS;
797 Double_t alpha=s.GetAlpha() - fSectors->GetAlphaShift();
798 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
799 if (alpha < 0. ) alpha += 2.*TMath::Pi();
800 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
801 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
804 if (!s.Rotate(alpha)) continue;
805 if (!FollowBackProlongation(s,t)) continue;
810 fSectors=fOuterSec; fN=fkNOS;
812 nc=s.GetNumberOfClusters();
814 alpha=s.GetAlpha() - fSectors->GetAlphaShift();
815 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
816 if (alpha < 0. ) alpha += 2.*TMath::Pi();
817 ns=Int_t(alpha/fSectors->GetAlpha())%fN;
819 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
822 if (!s.Rotate(alpha)) continue;
823 if (!FollowBackProlongation(s,t)) continue;
825 Int_t nrows=fOuterSec->GetNRows()+fInnerSec->GetNRows();
826 if (s.GetNumberOfClusters() < Int_t(0.4*nrows)) continue;
828 s.PropagateTo(fParam->GetOuterRadiusUp());
830 CookLabel(&s,0.1); //For comparison only
832 esd->UpdateTrackParams(&s,AliESDtrack::kTPCout);
835 Info("PropagateBack","Number of back propagated tracks: %d",ntrk);
840 //_________________________________________________________________________
841 AliCluster *AliTPCtracker::GetCluster(Int_t index) const {
842 //--------------------------------------------------------------------
843 // Return pointer to a given cluster
844 //--------------------------------------------------------------------
845 Int_t sec=(index&0xff000000)>>24;
846 Int_t row=(index&0x00ff0000)>>16;
847 Int_t ncl=(index&0x0000ffff)>>00;
849 const AliTPCcluster *cl=0;
852 cl=fInnerSec[sec][row].GetUnsortedCluster(ncl);
855 cl=fOuterSec[sec][row].GetUnsortedCluster(ncl);
858 return (AliCluster*)cl;
861 //__________________________________________________________________________
862 void AliTPCtracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const {
863 //--------------------------------------------------------------------
864 //This function "cooks" a track label. If label<0, this track is fake.
865 //--------------------------------------------------------------------
866 Int_t noc=t->GetNumberOfClusters();
867 Int_t *lb=new Int_t[noc];
868 Int_t *mx=new Int_t[noc];
869 AliCluster **clusters=new AliCluster*[noc];
872 for (i=0; i<noc; i++) {
874 Int_t index=t->GetClusterIndex(i);
875 clusters[i]=GetCluster(index);
879 for (i=0; i<noc; i++) {
880 AliCluster *c=clusters[i];
881 lab=TMath::Abs(c->GetLabel(0));
883 for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
889 for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
891 for (i=0; i<noc; i++) {
892 AliCluster *c=clusters[i];
893 if (TMath::Abs(c->GetLabel(1)) == lab ||
894 TMath::Abs(c->GetLabel(2)) == lab ) max++;
897 if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
900 Int_t tail=Int_t(0.10*noc);
902 for (i=1; i<=tail; i++) {
903 AliCluster *c=clusters[noc-i];
904 if (lab == TMath::Abs(c->GetLabel(0)) ||
905 lab == TMath::Abs(c->GetLabel(1)) ||
906 lab == TMath::Abs(c->GetLabel(2))) max++;
908 if (max < Int_t(0.5*tail)) lab=-lab;
918 //_________________________________________________________________________
919 void AliTPCtracker::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
920 //-----------------------------------------------------------------------
921 // Setup inner sector
922 //-----------------------------------------------------------------------
924 fAlpha=par->GetInnerAngle();
925 fAlphaShift=par->GetInnerAngleShift();
926 fPadPitchWidth=par->GetInnerPadPitchWidth();
927 f1PadPitchLength=par->GetInnerPadPitchLength();
928 f2PadPitchLength=f1PadPitchLength;
930 fN=par->GetNRowLow();
931 fRow=new AliTPCRow[fN];
932 for (Int_t i=0; i<fN; i++) fRow[i].SetX(par->GetPadRowRadiiLow(i));
934 fAlpha=par->GetOuterAngle();
935 fAlphaShift=par->GetOuterAngleShift();
936 fPadPitchWidth=par->GetOuterPadPitchWidth();
937 f1PadPitchLength=par->GetOuter1PadPitchLength();
938 f2PadPitchLength=par->GetOuter2PadPitchLength();
941 fRow=new AliTPCRow[fN];
942 for (Int_t i=0; i<fN; i++){
943 fRow[i].SetX(par->GetPadRowRadiiUp(i));
948 //_________________________________________________________________________
950 AliTPCRow::InsertCluster(const AliTPCcluster* c, Int_t sec, Int_t row) {
951 //-----------------------------------------------------------------------
952 // Insert a cluster into this pad row in accordence with its y-coordinate
953 //-----------------------------------------------------------------------
954 if (fN==kMaxClusterPerRow) {
955 ::Error("InsertCluster","Too many clusters !");
959 Int_t index=(((sec<<8)+row)<<16)+fN;
962 fSize=kMaxClusterPerRow/8;
963 fClusterArray=new AliTPCcluster[fSize];
965 fClusterArray[0]=*c; fClusters[fN++]=fClusterArray;
971 AliTPCcluster *buff=new AliTPCcluster[size];
972 memcpy(buff,fClusterArray,fSize*sizeof(AliTPCcluster));
973 for (Int_t i=0; i<fN; i++)
974 fClusters[i]=buff+(fClusters[i]-fClusterArray);
975 delete[] fClusterArray;
980 Int_t i=Find(c->GetY());
981 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCcluster*));
982 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
984 fClusters[i]=fClusterArray+fN; fClusterArray[fN++]=*c;
987 //___________________________________________________________________
988 Int_t AliTPCtracker::AliTPCRow::Find(Double_t y) const {
989 //-----------------------------------------------------------------------
990 // Return the index of the nearest cluster
991 //-----------------------------------------------------------------------
993 if (y <= fClusters[0]->GetY()) return 0;
994 if (y > fClusters[fN-1]->GetY()) return fN;
995 Int_t b=0, e=fN-1, m=(b+e)/2;
996 for (; b<e; m=(b+e)/2) {
997 if (y > fClusters[m]->GetY()) b=m+1;
1003 //_____________________________________________________________________________
1004 void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
1005 //-----------------------------------------------------------------
1006 // This funtion calculates dE/dX within the "low" and "up" cuts.
1007 //-----------------------------------------------------------------
1009 Int_t nc=GetNumberOfClusters();
1011 Int_t swap;//stupid sorting
1014 for (i=0; i<nc-1; i++) {
1015 if (fdEdxSample[i]<=fdEdxSample[i+1]) continue;
1016 Float_t tmp=fdEdxSample[i];
1017 fdEdxSample[i]=fdEdxSample[i+1]; fdEdxSample[i+1]=tmp;
1022 Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
1024 for (i=nl; i<=nu; i++) dedx += fdEdxSample[i];
1029 Double_t p=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt()));
1031 Double_t log1=TMath::Log(p+0.45), log2=TMath::Log(p+0.12);
1033 if (dedx < 34 + 30/(p+0.45)/(p+0.45) + 24*log1) {SetMass(0.13957); return;}
1034 if (dedx < 34 + 30/(p+0.12)/(p+0.12) + 24*log2) {SetMass(0.49368); return;}
1035 SetMass(0.93827); return;
1039 if (dedx < 34 + 30/(p+0.12)/(p+0.12) + 24*log2) {SetMass(0.13957); return;}
1040 SetMass(0.93827); return;
1043 SetMass(0.13957); return;