]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtracker.cxx
Fixing coding violation problems (Haavard)
[u/mrichter/AliRoot.git] / TPC / AliTPCtracker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //-------------------------------------------------------
19 //          Implementation of the TPC tracker
20 //
21 //   Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
22 //-------------------------------------------------------
23 #include <TObjArray.h>
24 #include <TError.h>
25 #include <TFile.h>
26 #include <TTree.h>
27
28 #include "AliESDEvent.h"
29 #include "AliESDtrack.h"
30
31 #include "AliTPCtracker.h"
32 #include "AliTPCcluster.h"
33 #include "AliTPCParam.h"
34 #include "AliClusters.h"
35
36 ClassImp(AliTPCtracker)
37
38 //_____________________________________________________________________________
39 AliTPCtracker::AliTPCtracker(): 
40   AliTracker(), 
41   fkNIS(0), 
42   fInnerSec(0),         
43   fkNOS(0),
44   fOuterSec(0),
45   fN(0),
46   fSectors(0),
47   fParam(0),
48   fSeeds(0)
49 {
50   //
51   // The default TPC tracker constructor
52   //
53 }
54
55 //_____________________________________________________________________________
56 AliTPCtracker::AliTPCtracker(const AliTPCParam *par): 
57   AliTracker(), 
58   fkNIS(par->GetNInnerSector()/2), 
59   fInnerSec(new AliTPCSector[fkNIS]),         
60   fkNOS(par->GetNOuterSector()/2),
61   fOuterSec(new AliTPCSector[fkNOS]),
62   fN(0),
63   fSectors(0),
64   fParam((AliTPCParam*) par),
65   fSeeds(0)
66 {
67   //---------------------------------------------------------------------
68   // The main TPC tracker constructor
69   //---------------------------------------------------------------------
70
71   Int_t i;
72   for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
73   for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
74
75 }
76
77 //_____________________________________________________________________________
78 AliTPCtracker::~AliTPCtracker() {
79   //------------------------------------------------------------------
80   // TPC tracker destructor
81   //------------------------------------------------------------------
82   delete[] fInnerSec;
83   delete[] fOuterSec;
84   if (fSeeds) {
85     fSeeds->Delete(); 
86     delete fSeeds;
87   }
88 }
89
90 //_____________________________________________________________________________
91 Double_t f1(Double_t x1,Double_t y1,
92                    Double_t x2,Double_t y2,
93                    Double_t x3,Double_t y3) 
94 {
95   //-----------------------------------------------------------------
96   // Initial approximation of the track curvature
97   //-----------------------------------------------------------------
98   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
99   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
100                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
101   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
102                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
103
104   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
105
106   return -xr*yr/sqrt(xr*xr+yr*yr); 
107 }
108
109
110 //_____________________________________________________________________________
111 Double_t f2(Double_t x1,Double_t y1,
112                    Double_t x2,Double_t y2,
113                    Double_t x3,Double_t y3) 
114 {
115   //-----------------------------------------------------------------
116   // Initial approximation of the track curvature times center of curvature
117   //-----------------------------------------------------------------
118   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
119   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
120                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
121   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
122                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
123
124   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
125   
126   return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
127 }
128
129 //_____________________________________________________________________________
130 Double_t f3(Double_t x1,Double_t y1, 
131                    Double_t x2,Double_t y2,
132                    Double_t z1,Double_t z2) 
133 {
134   //-----------------------------------------------------------------
135   // Initial approximation of the tangent of the track dip angle
136   //-----------------------------------------------------------------
137   return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
138 }
139
140 //_____________________________________________________________________________
141 Int_t AliTPCtracker::LoadClusters(TTree *cTree) {
142   //-----------------------------------------------------------------
143   // This function loads TPC clusters.
144   //-----------------------------------------------------------------
145   TBranch *branch=cTree->GetBranch("Segment");
146   if (!branch) {
147      Error("LoadClusters","Can't get the branch !");
148      return 1;
149   }
150
151   AliClusters carray, *addr=&carray;
152   carray.SetClass("AliTPCcluster");
153   carray.SetArray(0);
154   branch->SetAddress(&addr);
155
156   Int_t nentr=(Int_t)cTree->GetEntries();
157
158   for (Int_t i=0; i<nentr; i++) {
159       cTree->GetEvent(i);
160
161       Int_t ncl=carray.GetArray()->GetEntriesFast();
162
163       Int_t nir=fInnerSec->GetNRows(), nor=fOuterSec->GetNRows();
164       Int_t id=carray.GetID();
165       if ((id<0) || (id>2*(fkNIS*nir + fkNOS*nor))) {
166          Fatal("LoadClusters","Wrong index !");
167       }
168       Int_t outindex = 2*fkNIS*nir;
169       if (id<outindex) {
170          Int_t sec = id/nir;
171          Int_t row = id - sec*nir;
172          sec %= fkNIS;
173          AliTPCRow &padrow=fInnerSec[sec][row];
174          while (ncl--) {
175            AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
176            padrow.InsertCluster(c,sec,row);
177          }
178       } else {
179          id -= outindex;
180          Int_t sec = id/nor;
181          Int_t row = id - sec*nor;
182          sec %= fkNOS;
183          AliTPCRow &padrow=fOuterSec[sec][row];
184          while (ncl--) {
185            AliTPCcluster *c=(AliTPCcluster*)carray[ncl];
186            padrow.InsertCluster(c,sec+fkNIS,row);
187          }
188       }
189       carray.GetArray()->Clear();
190   }
191
192   return 0;
193 }
194
195 //_____________________________________________________________________________
196 void AliTPCtracker::UnloadClusters() {
197   //-----------------------------------------------------------------
198   // This function unloads TPC clusters.
199   //-----------------------------------------------------------------
200   Int_t i;
201   for (i=0; i<fkNIS; i++) {
202     Int_t nr=fInnerSec->GetNRows();
203     for (Int_t n=0; n<nr; n++) fInnerSec[i][n].ResetClusters();
204   }
205   for (i=0; i<fkNOS; i++) {
206     Int_t nr=fOuterSec->GetNRows();
207     for (Int_t n=0; n<nr; n++) fOuterSec[i][n].ResetClusters();
208   }
209 }
210
211 //_____________________________________________________________________________
212 Int_t AliTPCtracker::FollowProlongation(AliTPCseed& t, Int_t rf) {
213   //-----------------------------------------------------------------
214   // This function tries to find a track prolongation.
215   //-----------------------------------------------------------------
216   Double_t xt=t.GetX();
217   const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? kRowsToSkip : 
218                     Int_t(0.5*fSectors->GetNRows());
219   Int_t tryAgain=kSKIP;
220
221   Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
222   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();  
223   if (alpha < 0.            ) alpha += 2.*TMath::Pi();  
224   Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
225
226   Int_t nrows=fSectors->GetRowNumber(xt)-1;
227   for (Int_t nr=nrows; nr>=rf; nr--) {
228     Double_t x=fSectors->GetX(nr), ymax=fSectors->GetMaxY(nr);
229     if (!t.PropagateTo(x)) return 0;
230
231     AliTPCcluster *cl=0;
232     UInt_t index=0;
233     Double_t maxchi2=kMaxCHI2;
234     const AliTPCRow &krow=fSectors[s][nr];
235     Double_t pt=t.GetSignedPt();
236     Double_t sy2=AliTPCcluster::SigmaY2(t.GetX(),t.GetTgl(),pt);
237     Double_t sz2=AliTPCcluster::SigmaZ2(t.GetX(),t.GetTgl());
238     Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
239
240     if (road>kMaxROAD) {
241       if (t.GetNumberOfClusters()>4) 
242          Warning("FindProlongation","Too broad road !"); 
243       return 0;
244     }
245
246     if (krow) {
247       for (Int_t i=krow.Find(y-road); i<krow; i++) {
248         AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
249         if (c->GetY() > y+road) break;
250         if (c->IsUsed()) continue;
251         if ((c->GetZ()-z)*(c->GetZ()-z) > 16.*(t.GetSigmaZ2()+sz2)) continue;
252         Double_t chi2=t.GetPredictedChi2(c);
253         if (chi2 > maxchi2) continue;
254         maxchi2=chi2;
255         cl=c;
256         index=krow.GetIndex(i);       
257       }
258     }
259     if (cl) {
260       Float_t l=fSectors->GetPadPitchWidth();
261       Float_t corr=1.; if (nr>63) corr=0.67; // new (third) pad response !
262       t.SetSampledEdx(cl->GetQ()/l*corr,t.GetNumberOfClusters());
263       if (!t.Update(cl,maxchi2,index)) {
264          if (!tryAgain--) return 0;
265       } else tryAgain=kSKIP;
266     } else {
267       if (tryAgain==0) break;
268       if (y > ymax) {
269          s = (s+1) % fN;
270          if (!t.Rotate(fSectors->GetAlpha())) return 0;
271       } else if (y <-ymax) {
272          s = (s-1+fN) % fN;
273          if (!t.Rotate(-fSectors->GetAlpha())) return 0;
274       }
275       tryAgain--;
276     }
277   }
278
279   return 1;
280 }
281 //_____________________________________________________________________________
282
283 Int_t AliTPCtracker::FollowRefitInward(AliTPCseed *seed, AliTPCtrack *track) {
284   //
285   // This function propagates seed inward TPC using old clusters
286   // from the track.
287   // 
288   // Sylwester Radomski, GSI
289   // 26.02.2003
290   //
291
292   // loop over rows
293
294   Int_t nRows = fSectors->GetNRows();
295   for (Int_t iRow = nRows-1; iRow >= 0; iRow--) {
296
297     Double_t x = fSectors->GetX(iRow);
298     if (!seed->PropagateTo(x)) return 0;
299
300     // try to find an assigned cluster in this row
301
302     AliTPCcluster* cluster = NULL;
303     Int_t idx = -1;
304     Int_t sec = -1;
305     for (Int_t iCluster = 0; iCluster < track->GetNumberOfClusters(); iCluster++) {
306       idx = track->GetClusterIndex(iCluster); 
307       sec = (idx&0xff000000)>>24; 
308       Int_t row = (idx&0x00ff0000)>>16;
309       if (((fSectors == fInnerSec) && (sec >= fkNIS)) ||
310           ((fSectors == fOuterSec) && (sec < fkNIS))) continue;
311       if (row == iRow) {
312         cluster = (AliTPCcluster*) GetCluster(idx);
313         break;
314       }
315     }
316
317     // update the track seed with the found cluster
318
319     if (cluster) {
320       Double_t dAlpha = fParam->GetAngle(sec) - seed->GetAlpha();
321       if (TMath::Abs(dAlpha) > 0.0001) {
322         if (!seed->Rotate(dAlpha)) return 0;
323         if (!seed->PropagateTo(x)) return 0;
324       }
325
326       seed->Update(cluster, seed->GetPredictedChi2(cluster), idx);
327     }
328   }
329
330   return 1;
331 }
332
333 //_____________________________________________________________________________
334 Int_t AliTPCtracker::FollowBackProlongation
335 (AliTPCseed& seed, const AliTPCtrack &track) {
336   //-----------------------------------------------------------------
337   // This function propagates tracks back through the TPC
338   //-----------------------------------------------------------------
339   Double_t alpha=seed.GetAlpha() - fSectors->GetAlphaShift();
340   if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();  
341   if (alpha < 0.            ) alpha += 2.*TMath::Pi();  
342   Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN;
343
344   Int_t idx=-1, sec=-1, row=-1;
345   Int_t nc=track.GetNumberOfClusters();
346
347   if (nc--) {
348      idx=track.GetClusterIndex(nc);
349      sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
350   }
351   if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; } 
352   else { if (sec <  fkNIS) row=-1; }   
353
354   Int_t nr=fSectors->GetNRows();
355   for (Int_t i=0; i<nr; i++) {
356     Double_t x=fSectors->GetX(i), ymax=fSectors->GetMaxY(i);
357     Double_t y;
358     if (!seed.GetYAt(x,GetBz(),y)) return 0;
359  
360     if (y > ymax) {
361        s = (s+1) % fN;
362        if (!seed.Rotate(fSectors->GetAlpha())) return 0;
363     } else if (y <-ymax) {
364        s = (s-1+fN) % fN;
365        if (!seed.Rotate(-fSectors->GetAlpha())) return 0;
366     }
367
368     if (!seed.PropagateTo(x)) return 0;
369
370     AliTPCcluster *cl=0;
371     Int_t index=0;
372     Double_t maxchi2=kMaxCHI2;
373     Double_t pt=seed.GetSignedPt();
374     Double_t sy2=AliTPCcluster::SigmaY2(seed.GetX(),seed.GetTgl(),pt);
375     Double_t sz2=AliTPCcluster::SigmaZ2(seed.GetX(),seed.GetTgl());
376     Double_t road=4.*sqrt(seed.GetSigmaY2() + sy2), z=seed.GetZ();
377     if (road>kMaxROAD) {
378       Warning("FollowBackProlongation","Too broad road !"); 
379       return 0;
380     }
381
382     Int_t accepted=seed.GetNumberOfClusters();
383     if (row==i) {
384        //try to accept already found cluster
385        AliTPCcluster *c=(AliTPCcluster*)GetCluster(idx);
386        Double_t chi2;
387        if ((chi2=seed.GetPredictedChi2(c))<maxchi2 || accepted<27) { 
388           index=idx; cl=c; maxchi2=chi2;
389        } //else cerr<<"AliTPCtracker::FollowBackProlongation: oulier !\n";
390           
391        if (nc--) {
392           idx=track.GetClusterIndex(nc); 
393           sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16;
394        } 
395        if (fSectors==fInnerSec) { if (sec >= fkNIS) row=-1; }
396        else { if (sec < fkNIS) row=-1; }   
397
398     }
399     if (!cl) {
400        //try to fill the gap
401        const AliTPCRow &krow=fSectors[s][i];
402        if (accepted>27)
403        if (krow) {
404           for (Int_t i=krow.Find(y-road); i<krow; i++) {
405             AliTPCcluster *c=(AliTPCcluster*)(krow[i]);
406             if (c->GetY() > y+road) break;
407             if (c->IsUsed()) continue;
408          if ((c->GetZ()-z)*(c->GetZ()-z)>16.*(seed.GetSigmaZ2()+sz2)) continue;
409             Double_t chi2=seed.GetPredictedChi2(c);
410             if (chi2 > maxchi2) continue;
411             maxchi2=chi2;
412             cl=c;
413             index=krow.GetIndex(i);
414           }
415        }
416     }
417
418     if (cl) {
419       Float_t l=fSectors->GetPadPitchWidth();
420       Float_t corr=1.; if (i>63) corr=0.67; // new (third) pad response !
421       seed.SetSampledEdx(cl->GetQ()/l*corr,seed.GetNumberOfClusters());
422       seed.Update(cl,maxchi2,index);
423     }
424
425   }
426
427   return 1;
428 }
429
430 //_____________________________________________________________________________
431 void AliTPCtracker::MakeSeeds(Int_t i1, Int_t i2) {
432   //-----------------------------------------------------------------
433   // This function creates track seeds.
434   //-----------------------------------------------------------------
435   Double_t x[5], c[15];
436
437   Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
438   Double_t cs=cos(alpha), sn=sin(alpha);
439
440   Double_t x1 =fSectors->GetX(i1);
441   Double_t xx2=fSectors->GetX(i2);
442
443   for (Int_t ns=0; ns<fN; ns++) {
444     Int_t nl=fSectors[(ns-1+fN)%fN][i2];
445     Int_t nm=fSectors[ns][i2];
446     Int_t nu=fSectors[(ns+1)%fN][i2];
447     const AliTPCRow& kr1=fSectors[ns][i1];
448     for (Int_t is=0; is < kr1; is++) {
449       Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
450       for (Int_t js=0; js < nl+nm+nu; js++) {
451         const AliTPCcluster *kcl;
452         Double_t x2,   y2,   z2;
453         Double_t x3=GetX(), y3=GetY(), z3=GetZ();
454
455         if (js<nl) {
456           const AliTPCRow& kr2=fSectors[(ns-1+fN)%fN][i2];
457           kcl=kr2[js];
458           y2=kcl->GetY(); z2=kcl->GetZ();
459           x2= xx2*cs+y2*sn;
460           y2=-xx2*sn+y2*cs;
461         } else 
462           if (js<nl+nm) {
463             const AliTPCRow& kr2=fSectors[ns][i2];
464             kcl=kr2[js-nl];
465             x2=xx2; y2=kcl->GetY(); z2=kcl->GetZ();
466           } else {
467             const AliTPCRow& kr2=fSectors[(ns+1)%fN][i2];
468             kcl=kr2[js-nl-nm];
469             y2=kcl->GetY(); z2=kcl->GetZ();
470             x2=xx2*cs-y2*sn;
471             y2=xx2*sn+y2*cs;
472           }
473
474         Double_t zz=z1 - (z1-z3)/(x1-x3)*(x1-x2); 
475         if (TMath::Abs(zz-z2)>5.) continue;
476
477         Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
478         if (d==0.) {
479            Warning("MakeSeeds","Straight seed !"); 
480            continue;
481         }
482         x[0]=y1;
483         x[1]=z1;
484         x[4]=f1(x1,y1,x2,y2,x3,y3);
485         if (TMath::Abs(x[4]) >= 0.0066) continue;
486         x[2]=f2(x1,y1,x2,y2,x3,y3);
487         //if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
488         x[3]=f3(x1,y1,x2,y2,z1,z2);
489         if (TMath::Abs(x[3]) > 1.2) continue;
490         Double_t a=asin(x[2]);
491         Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
492         if (TMath::Abs(zv-z3)>10.) continue; 
493
494         Double_t sy1=kr1[is]->GetSigmaY2(), sz1=kr1[is]->GetSigmaZ2();
495         Double_t sy2=kcl->GetSigmaY2(),     sz2=kcl->GetSigmaZ2();
496         //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
497         Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
498
499         Double_t f40=(f1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
500         Double_t f42=(f1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
501         Double_t f43=(f1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
502         Double_t f20=(f2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
503         Double_t f22=(f2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
504         Double_t f23=(f2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
505         Double_t f30=(f3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
506         Double_t f31=(f3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
507         Double_t f32=(f3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
508         Double_t f34=(f3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
509
510         c[0]=sy1;
511         c[1]=0.;       c[2]=sz1;
512         c[3]=f20*sy1;  c[4]=0.;       c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
513         c[6]=f30*sy1;  c[7]=f31*sz1;  c[8]=f30*sy1*f20+f32*sy2*f22;
514                        c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
515         c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
516         c[13]=f30*sy1*f40+f32*sy2*f42;
517         c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
518
519         Int_t index=kr1.GetIndex(is);
520         AliTPCseed *track=new AliTPCseed(x1, ns*alpha+shift, x, c, index);
521         Float_t l=fSectors->GetPadPitchWidth();
522         track->SetSampledEdx(kr1[is]->GetQ()/l,0);
523
524         Int_t rc=FollowProlongation(*track, i2);
525         if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
526         else fSeeds->AddLast(track); 
527       }
528     }
529   }
530 }
531
532 //_____________________________________________________________________________
533 Int_t AliTPCtracker::ReadSeeds(const TFile *inp) {
534   //-----------------------------------------------------------------
535   // This function reades track seeds.
536   //-----------------------------------------------------------------
537   TDirectory *savedir=gDirectory; 
538
539   TFile *in=(TFile*)inp;
540   if (!in->IsOpen()) {
541      Error("ReadSeeds","Input file has not been open !");
542      return 1;
543   }
544
545   in->cd();
546   TTree *seedTree=(TTree*)in->Get("Seeds");
547   if (!seedTree) {
548      Error("ReadSeeds","Can't get a tree with track seeds !");
549      return 2;
550   }
551   AliTPCtrack *seed=new AliTPCtrack; 
552   seedTree->SetBranchAddress("tracks",&seed);
553   
554   if (fSeeds==0) fSeeds=new TObjArray(15000);
555
556   Int_t n=(Int_t)seedTree->GetEntries();
557   for (Int_t i=0; i<n; i++) {
558      seedTree->GetEvent(i);
559      seed->ResetClusters();
560      fSeeds->AddLast(new AliTPCseed(*seed));
561   }
562   
563   delete seed;
564
565   delete seedTree; //Thanks to Mariana Bondila
566
567   savedir->cd();
568
569   return 0;
570 }
571
572 //_____________________________________________________________________________
573 Int_t AliTPCtracker::Clusters2Tracks(AliESDEvent *event) {
574   //-----------------------------------------------------------------
575   // This is a track finder.
576   // The clusters must be already loaded ! 
577   //-----------------------------------------------------------------
578
579   //find track seeds
580   Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows();
581   Int_t nrows=nlow+nup;
582   if (fSeeds==0) {
583      Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
584      fSectors=fOuterSec; fN=fkNOS;
585      fSeeds=new TObjArray(15000);
586      MakeSeeds(nup-1, nup-1-gap);
587      MakeSeeds(nup-1-shift, nup-1-shift-gap);
588   }
589   fSeeds->Sort();
590
591   Int_t nseed=fSeeds->GetEntriesFast();
592   for (Int_t i=0; i<nseed; i++) {
593     //tracking in the outer sectors
594     fSectors=fOuterSec; fN=fkNOS;
595
596     AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
597     if (!FollowProlongation(t)) {
598        delete fSeeds->RemoveAt(i);
599        continue;
600     }
601
602     //tracking in the inner sectors
603     fSectors=fInnerSec; fN=fkNIS;
604
605     Double_t alpha=t.GetAlpha() - fInnerSec->GetAlphaShift();
606     if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
607     if (alpha < 0.            ) alpha += 2.*TMath::Pi();
608     Int_t ns=Int_t(alpha/fInnerSec->GetAlpha())%fkNIS;
609
610     alpha=ns*fInnerSec->GetAlpha()+fInnerSec->GetAlphaShift()-t.GetAlpha();
611
612     if (t.Rotate(alpha)) {
613       if (FollowProlongation(t)) {
614         if (t.GetNumberOfClusters() >= Int_t(0.4*nrows)) {
615           t.CookdEdx();
616           CookLabel(pt,0.1); //For comparison only
617           pt->PropagateTo(fParam->GetInnerRadiusLow());
618           AliESDtrack iotrack;
619           iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
620
621           event->AddTrack(&iotrack);
622
623           UseClusters(&t);
624         }
625       }
626     }
627     delete fSeeds->RemoveAt(i);
628   }
629
630   Info("Clusters2Tracks","Number of found tracks : %d",
631        event->GetNumberOfTracks());
632
633   fSeeds->Clear(); delete fSeeds; fSeeds=0;
634
635   return 0;
636 }
637
638 //_____________________________________________________________________________
639 Int_t AliTPCtracker::RefitInward(AliESDEvent* event) {
640   //
641   // The function propagates tracks throught TPC inward
642   // using already associated clusters.
643   // The clusters must be already loaded !
644   //
645
646   Int_t nTracks = event->GetNumberOfTracks();
647   Int_t nRefited = 0;
648
649   for (Int_t i = 0; i < nTracks; i++) {
650     AliESDtrack* track = event->GetTrack(i);
651     ULong_t status = track->GetStatus();
652
653     if ( (status & AliESDtrack::kTPCrefit) != 0 ) continue;    
654     if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
655
656     if ( (status & AliESDtrack::kTRDout ) != 0 ) 
657       if ( (status & AliESDtrack::kTRDrefit ) == 0 ) continue;
658
659     AliTPCtrack* tpcTrack = new AliTPCtrack(*track);
660     AliTPCseed* seed=new AliTPCseed(*tpcTrack); seed->ResetClusters(); 
661
662     if ( (status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance(10.);
663
664     fSectors = fOuterSec;
665
666     Int_t res = FollowRefitInward(seed, tpcTrack);
667     UseClusters(seed);
668     Int_t nc = seed->GetNumberOfClusters();
669
670     fSectors = fInnerSec;
671
672     res = FollowRefitInward(seed, tpcTrack);
673     UseClusters(seed, nc);
674
675     if (res) {
676       seed->PropagateTo(fParam->GetInnerRadiusLow());
677       seed->SetLabel(tpcTrack->GetLabel());
678       seed->SetdEdx(tpcTrack->GetdEdx());
679       track->UpdateTrackParams(seed, AliESDtrack::kTPCrefit);
680       nRefited++;
681     }
682
683     delete seed;
684     delete tpcTrack;
685   }
686
687   Info("RefitInward","Number of refitted tracks : %d",nRefited);
688
689   return 0;
690 }
691
692 Int_t AliTPCtracker::PropagateBack(AliESDEvent *event) {
693   //-----------------------------------------------------------------
694   // This function propagates tracks back through the TPC.
695   // The clusters must be already loaded !
696   //-----------------------------------------------------------------
697   Int_t nentr=event->GetNumberOfTracks();
698   Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
699
700   Int_t ntrk=0;
701   for (Int_t i=0; i<nentr; i++) {
702     AliESDtrack *esd=event->GetTrack(i);
703     ULong_t status=esd->GetStatus();
704
705     if ( (status & AliESDtrack::kTPCin ) == 0 ) continue;
706     if ( (status & AliESDtrack::kTPCout) != 0 ) continue;
707     if ( (status & AliESDtrack::kITSin) != 0 )
708        if ( (status & AliESDtrack::kITSout) == 0 ) continue;
709
710     const AliTPCtrack t(*esd);
711     AliTPCseed s(t); s.ResetClusters();
712
713     if ( (status & AliESDtrack::kITSout) == 0 ) s.ResetCovariance(10.);
714
715     //inner sectors
716     fSectors=fInnerSec; fN=fkNIS;
717
718     Double_t alpha=s.GetAlpha() - fSectors->GetAlphaShift();
719     if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
720     if (alpha < 0.            ) alpha += 2.*TMath::Pi();
721     Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
722     alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
723     alpha-=s.GetAlpha();
724
725     if (!s.Rotate(alpha)) continue;
726     if (!FollowBackProlongation(s,t)) continue;
727
728     UseClusters(&s);
729
730     //outer sectors
731     fSectors=fOuterSec; fN=fkNOS;
732
733     Int_t nc=s.GetNumberOfClusters();
734
735     alpha=s.GetAlpha() - fSectors->GetAlphaShift();
736     if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
737     if (alpha < 0.            ) alpha += 2.*TMath::Pi();
738     ns=Int_t(alpha/fSectors->GetAlpha())%fN;
739
740     alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
741     alpha-=s.GetAlpha();
742
743     if (!s.Rotate(alpha)) continue;
744     if (!FollowBackProlongation(s,t)) continue;
745     {
746     Int_t nrows=fOuterSec->GetNRows()+fInnerSec->GetNRows();
747     if (s.GetNumberOfClusters() < Int_t(0.4*nrows)) continue;
748     }
749     s.PropagateTo(fParam->GetOuterRadiusUp());
750     s.CookdEdx();
751     CookLabel(&s,0.1); //For comparison only
752     UseClusters(&s,nc);
753     esd->UpdateTrackParams(&s,AliESDtrack::kTPCout);
754     ntrk++;
755   }
756   Info("PropagateBack","Number of back propagated tracks: %d",ntrk);
757
758   return 0;
759 }
760
761 //_________________________________________________________________________
762 AliCluster *AliTPCtracker::GetCluster(Int_t index) const {
763   //--------------------------------------------------------------------
764   //       Return pointer to a given cluster
765   //--------------------------------------------------------------------
766   Int_t sec=(index&0xff000000)>>24; 
767   Int_t row=(index&0x00ff0000)>>16; 
768   Int_t ncl=(index&0x0000ffff)>>00;
769
770   const AliTPCcluster *cl=0;
771
772   if (sec<fkNIS) {
773     cl=fInnerSec[sec][row].GetUnsortedCluster(ncl); 
774   } else {
775     sec-=fkNIS;
776     cl=fOuterSec[sec][row].GetUnsortedCluster(ncl);
777   }
778
779   return (AliCluster*)cl;      
780 }
781
782 //__________________________________________________________________________
783 void AliTPCtracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const {
784   //--------------------------------------------------------------------
785   //This function "cooks" a track label. If label<0, this track is fake.
786   //--------------------------------------------------------------------
787   Int_t noc=t->GetNumberOfClusters();
788   Int_t *lb=new Int_t[noc];
789   Int_t *mx=new Int_t[noc];
790   AliCluster **clusters=new AliCluster*[noc];
791
792   Int_t i;
793   for (i=0; i<noc; i++) {
794      lb[i]=mx[i]=0;
795      Int_t index=t->GetClusterIndex(i);
796      clusters[i]=GetCluster(index);
797   }
798
799   Int_t lab=123456789;
800   for (i=0; i<noc; i++) {
801     AliCluster *c=clusters[i];
802     lab=TMath::Abs(c->GetLabel(0));
803     Int_t j;
804     for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
805     lb[j]=lab;
806     (mx[j])++;
807   }
808
809   Int_t max=0;
810   for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
811     
812   for (i=0; i<noc; i++) {
813     AliCluster *c=clusters[i];
814     if (TMath::Abs(c->GetLabel(1)) == lab ||
815         TMath::Abs(c->GetLabel(2)) == lab ) max++;
816   }
817
818   if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
819
820   else {
821      Int_t tail=Int_t(0.10*noc);
822      max=0;
823      for (i=1; i<=tail; i++) {
824        AliCluster *c=clusters[noc-i];
825        if (lab == TMath::Abs(c->GetLabel(0)) ||
826            lab == TMath::Abs(c->GetLabel(1)) ||
827            lab == TMath::Abs(c->GetLabel(2))) max++;
828      }
829      if (max < Int_t(0.5*tail)) lab=-lab;
830   }
831
832   t->SetLabel(lab);
833
834   delete[] lb;
835   delete[] mx;
836   delete[] clusters;
837 }
838
839 //_________________________________________________________________________
840 void AliTPCtracker::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
841   //-----------------------------------------------------------------------
842   // Setup inner sector
843   //-----------------------------------------------------------------------
844   if (f==0) {
845      fAlpha=par->GetInnerAngle();
846      fAlphaShift=par->GetInnerAngleShift();
847      fPadPitchWidth=par->GetInnerPadPitchWidth();
848      f1PadPitchLength=par->GetInnerPadPitchLength();
849      f2PadPitchLength=f1PadPitchLength;
850
851      fN=par->GetNRowLow();
852      fRow=new AliTPCRow[fN];
853      for (Int_t i=0; i<fN; i++) fRow[i].SetX(par->GetPadRowRadiiLow(i));
854   } else {
855      fAlpha=par->GetOuterAngle();
856      fAlphaShift=par->GetOuterAngleShift();
857      fPadPitchWidth=par->GetOuterPadPitchWidth();
858      f1PadPitchLength=par->GetOuter1PadPitchLength();
859      f2PadPitchLength=par->GetOuter2PadPitchLength();
860
861      fN=par->GetNRowUp();
862      fRow=new AliTPCRow[fN];
863      for (Int_t i=0; i<fN; i++){ 
864      fRow[i].SetX(par->GetPadRowRadiiUp(i));
865      }
866   } 
867 }
868
869 //_________________________________________________________________________
870 void AliTPCtracker::
871 AliTPCRow::InsertCluster(const AliTPCcluster* c, Int_t sec, Int_t row) {
872   //-----------------------------------------------------------------------
873   // Insert a cluster into this pad row in accordence with its y-coordinate
874   //-----------------------------------------------------------------------
875   if (fN==kMaxClusterPerRow) {
876      ::Error("InsertCluster","Too many clusters !"); 
877      return;
878   }
879
880   Int_t index=(((sec<<8)+row)<<16)+fN;
881
882   if (fN==0) {
883      fSize=kMaxClusterPerRow/8;
884      fClusterArray=new AliTPCcluster[fSize];
885      fIndex[0]=index;
886      fClusterArray[0]=*c; fClusters[fN++]=fClusterArray; 
887      return;
888   }
889
890   if (fN==fSize) {
891      Int_t size=fSize*2;
892      AliTPCcluster *buff=new AliTPCcluster[size];
893      memcpy(buff,fClusterArray,fSize*sizeof(AliTPCcluster));
894      for (Int_t i=0; i<fN; i++) 
895         fClusters[i]=buff+(fClusters[i]-fClusterArray);
896      delete[] fClusterArray;
897      fClusterArray=buff;
898      fSize=size;
899   }
900
901   Int_t i=Find(c->GetY());
902   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCcluster*));
903   memmove(fIndex   +i+1 ,fIndex   +i,(fN-i)*sizeof(UInt_t));
904   fIndex[i]=index; 
905   fClusters[i]=fClusterArray+fN; fClusterArray[fN++]=*c;
906 }
907
908 //___________________________________________________________________
909 Int_t AliTPCtracker::AliTPCRow::Find(Double_t y) const {
910   //-----------------------------------------------------------------------
911   // Return the index of the nearest cluster 
912   //-----------------------------------------------------------------------
913   if(fN<=0) return 0;
914   if (y <= fClusters[0]->GetY()) return 0;
915   if (y > fClusters[fN-1]->GetY()) return fN;
916   Int_t b=0, e=fN-1, m=(b+e)/2;
917   for (; b<e; m=(b+e)/2) {
918     if (y > fClusters[m]->GetY()) b=m+1;
919     else e=m; 
920   }
921   return m;
922 }
923
924 //_____________________________________________________________________________
925 void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
926   //-----------------------------------------------------------------
927   // This funtion calculates dE/dX within the "low" and "up" cuts.
928   //-----------------------------------------------------------------
929   Int_t i;
930   Int_t nc=GetNumberOfClusters();
931
932   Int_t swap;//stupid sorting
933   do {
934     swap=0;
935     for (i=0; i<nc-1; i++) {
936       if (fdEdxSample[i]<=fdEdxSample[i+1]) continue;
937       Float_t tmp=fdEdxSample[i];
938       fdEdxSample[i]=fdEdxSample[i+1]; fdEdxSample[i+1]=tmp;
939       swap++;
940     }
941   } while (swap);
942
943   Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
944   Float_t dedx=0;
945   for (i=nl; i<=nu; i++) dedx += fdEdxSample[i];
946   dedx /= (nu-nl+1);
947   SetdEdx(dedx);
948
949   return;
950 }
951
952
953