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