]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUTrackerCooked.cxx
c0f03d45cedeb78f229cc96af46b6053b52238f2
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUTrackerCooked.cxx
1 //-------------------------------------------------------------------------
2 //               Implementation of the ITS tracker class
3 //    The pattern recongintion based on the "cooked covariance" approach
4 //-------------------------------------------------------------------------
5
6 #include <TTree.h>
7 #include <TClonesArray.h>
8
9 #include "AliLog.h"
10 #include "AliESDEvent.h"
11 #include "AliITSUClusterPix.h"
12 #include "AliITSUGeomTGeo.h"
13 #include "AliITSUTrackerCooked.h"
14 #include "AliITSUTrackCooked.h" 
15
16 ClassImp(AliITSUTrackerCooked)
17
18 //************************************************
19 // Constants hardcoded for the moment:
20 //************************************************
21 // radial positions of layers: default contructor
22 const 
23 Double_t klRadius[7]={2.34, 3.15, 3.93, 19.61, 24.55, 34.39, 39.34}; //tdr6
24 // seed "windows" in z and phi: MakeSeeds
25 const Double_t kzWin=0.33, kpWin=3.14/4;
26 // Maximal accepted impact parameters for the seeds 
27 const Double_t kmaxDCAxy=3.;
28 const Double_t kmaxDCAz= 3.;
29 // Layers for the seeding
30 const Int_t kSeedingLayer1=6, kSeedingLayer2=4, kSeedingLayer3=5;
31 // Space point resolution
32 const Double_t kSigma2=0.0005*0.0005;
33 // Max accepted chi2 per cluster
34 const Double_t kmaxChi2PerCluster=10.;
35 // Tracking "road" from layer to layer
36 const Double_t kRoadY=0.7;
37 const Double_t kRoadZ=0.7;
38 // Minimal number of attached clusters
39 const Int_t kminNumberOfClusters=3;
40
41 //************************************************
42 // TODO:
43 //************************************************
44 // Seeding:
45 // Precalculate cylidnrical (r,phi) for the clusters;
46 // use exact r's for the clusters
47
48
49 AliITSUTrackerCooked::AliITSUlayer
50               AliITSUTrackerCooked::fgLayers[AliITSUTrackerCooked::kNLayers];
51
52 AliITSUTrackerCooked::AliITSUTrackerCooked(): 
53 AliTracker(),
54 fSeeds(0),
55 fI(kNLayers-1),
56 fBestTrack(0), 
57 fTrackToFollow(0) 
58 {
59   //--------------------------------------------------------------------
60   // This default constructor needs to be provided
61   //--------------------------------------------------------------------
62
63   AliITSUGeomTGeo *gm  = new AliITSUGeomTGeo(kTRUE,kTRUE);
64   AliITSUClusterPix::SetGeom(gm);
65
66   for (Int_t i=0; i<kNLayers; i++) fgLayers[i].SetR(klRadius[i]);
67
68   // Some default primary vertex
69   Double_t xyz[]={0.,0.,0.};
70   Double_t ers[]={2.,2.,2.};
71
72   SetVertex(xyz,ers);
73
74 }
75
76 void AliITSUTrackerCooked::ResetTrackToFollow(const AliITSUTrackCooked &t) {
77   //--------------------------------------------------------------------
78   // Prepare to follow a new track seed
79   //--------------------------------------------------------------------
80      delete fTrackToFollow;
81      fTrackToFollow = new AliITSUTrackCooked(t);
82 }
83   
84 void AliITSUTrackerCooked::ResetBestTrack() {
85   //--------------------------------------------------------------------
86   // Replace the best track branch
87   //--------------------------------------------------------------------
88      delete fBestTrack;
89      fBestTrack = new AliITSUTrackCooked(*fTrackToFollow);
90 }
91   
92 AliITSUTrackerCooked::~AliITSUTrackerCooked() 
93 {
94   //--------------------------------------------------------------------
95   // Virtual destructor
96   //--------------------------------------------------------------------
97
98   if (fSeeds) fSeeds->Delete(); delete fSeeds; 
99   delete fBestTrack;
100   delete fTrackToFollow;
101
102 }
103
104 static Double_t 
105 f1(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x3, Double_t y3)
106 {
107     //-----------------------------------------------------------------
108     // Initial approximation of the track curvature
109     //-----------------------------------------------------------------
110     Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
111     Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
112                     (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
113     Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
114                     (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
115     
116     Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
117     
118     return -xr*yr/sqrt(xr*xr+yr*yr);
119 }
120
121 static Double_t 
122 f2(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x3, Double_t y3)
123 {
124     //-----------------------------------------------------------------
125     // Initial approximation of the track curvature times center of curvature
126     //-----------------------------------------------------------------
127     Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
128     Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
129                     (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
130     Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
131                     (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
132     
133     Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
134     
135     return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
136 }
137
138 static Double_t 
139 f3(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t z1, Double_t z2)
140 {
141     //-----------------------------------------------------------------
142     // Initial approximation of the tangent of the track dip angle
143     //-----------------------------------------------------------------
144     return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
145 }
146
147 Bool_t AliITSUTrackerCooked::
148 AddCookedSeed(const Float_t r1[3], Int_t l1, Int_t i1, 
149               const Float_t r2[3], Int_t l2, Int_t i2,
150               const AliCluster *c3,Int_t l3, Int_t i3) 
151 {
152     //--------------------------------------------------------------------
153     // This is the main cooking function.
154     // Creates seed parameters out of provided clusters.
155     //--------------------------------------------------------------------
156     Float_t x,a;
157     if (!c3->GetXAlphaRefPlane(x,a)) return kFALSE;
158
159     Double_t ca=TMath::Cos(a), sa=TMath::Sin(a);
160     Double_t x1 = r1[0]*ca + r1[1]*sa,
161              y1 =-r1[0]*sa + r1[1]*ca, z1 = r1[2];
162     Double_t x2 = r2[0]*ca + r2[1]*sa,
163              y2 =-r2[0]*sa + r2[1]*ca, z2 = r2[2];
164     Double_t x3 = x,  y3 = c3->GetY(), z3 = c3->GetZ();
165
166     Double_t par[5];
167     par[0]=y3;
168     par[1]=z3;
169     Double_t crv=f1(x1, y1, x2, y2, x3, y3); //curvature
170     Double_t cx0=f2(x1, y1, x2, y2, x3, y3); //curvature*x0
171     Double_t tgl12=f3(x1, y1, x2, y2, z1, z2);
172     Double_t tgl23=f3(x2, y2, x3, y3, z2, z3);
173
174     Double_t sf=x*crv - cx0;
175     if (TMath::Abs(sf) >= kAlmost1) return kFALSE;
176     par[2]=sf;
177
178     par[3]=0.5*(tgl12 + tgl23);
179     Double_t bz=GetBz();
180     par[4]=(TMath::Abs(bz) < kAlmost0Field) ? kAlmost0 : crv/(bz*kB2C);
181
182     Double_t cov[15];
183     for (Int_t i=0; i<15; i++) cov[i]=0.;
184     cov[0] =kSigma2*10;
185     cov[2] =kSigma2*10;
186     cov[5] =0.007*0.007*10;   //FIXME all these lines
187     cov[9] =0.007*0.007*10;
188     cov[14]=0.1*0.1*10;
189
190     AliITSUTrackCooked *seed=new AliITSUTrackCooked();
191     seed->Set(Double_t(x), Double_t(a), par, cov);
192
193     Float_t dz[2]; 
194     seed->GetDZ(GetX(),GetY(),GetZ(),GetBz(),dz);
195     if (TMath::Abs(dz[0]) > kmaxDCAxy) {delete seed; return kFALSE;} 
196     if (TMath::Abs(dz[1]) > kmaxDCAz ) {delete seed; return kFALSE;} 
197
198     seed->SetClusterIndex(l1,i1);
199     seed->SetClusterIndex(l2,i2);
200     seed->SetClusterIndex(l3,i3);
201
202     fSeeds->AddLast(seed);
203
204     return kTRUE;
205 }
206
207 Int_t AliITSUTrackerCooked::MakeSeeds() {
208   //--------------------------------------------------------------------
209   // This is the main pattern recongition function.
210   // Creates seeds out of two clusters and another point.
211   //--------------------------------------------------------------------
212    if (fSeeds) {fSeeds->Delete(); delete fSeeds;}
213    fSeeds=new TObjArray(77777);
214
215    Double_t zv=GetZ();
216
217    AliITSUlayer &layer1=fgLayers[kSeedingLayer1];
218    AliITSUlayer &layer2=fgLayers[kSeedingLayer2];
219    AliITSUlayer &layer3=fgLayers[kSeedingLayer3];
220    Double_t r1=layer1.GetR();
221    Double_t r2=layer2.GetR();
222    Double_t r3=layer3.GetR();
223
224    Int_t nClusters1=layer1.GetNumberOfClusters();
225    Int_t nClusters2=layer2.GetNumberOfClusters();
226    Int_t nClusters3=layer3.GetNumberOfClusters();
227    for (Int_t n1=0; n1<nClusters1; n1++) {
228      AliCluster *c1=layer1.GetCluster(n1);
229      //
230      //Int_t lab=c1->GetLabel(0);
231      //
232      Double_t z1=c1->GetZ();
233      Float_t xyz1[3]; c1->GetGlobalXYZ(xyz1);
234      Double_t phi1=TMath::ATan2(xyz1[1],xyz1[0]);
235      Double_t zr2=zv + r2/r1*(z1-zv);
236      Int_t start2=layer2.FindClusterIndex(zr2-kzWin);
237      for (Int_t n2=start2; n2<nClusters2; n2++) {
238          AliCluster *c2=layer2.GetCluster(n2);
239          //
240          //if (c2->GetLabel(0)!=lab) continue;
241          //
242          Double_t z2=c2->GetZ();
243          
244          if (z2 > (zr2+kzWin)) break;  //check in Z
245          Float_t xyz2[3]; c2->GetGlobalXYZ(xyz2);
246          Double_t phi2=TMath::ATan2(xyz2[1],xyz2[0]);
247          if (TMath::Abs(phi2-phi1) > kpWin) continue;  //check in Phi
248  
249          Double_t zr3=z1 + (r3-r1)/(r2-r1)*(z2-z1);
250          Double_t d13=r1-r3, d32=r3-r2, d=r1-r2;
251          Double_t phir3=d32/d*phi1 + d13/d*phi2;  // FIXME
252          Int_t start3=layer3.FindClusterIndex(zr3-kzWin/2);
253          for (Int_t n3=start3; n3<nClusters3; n3++) {
254              AliCluster *c3=layer3.GetCluster(n3);
255              //
256              //if (c3->GetLabel(0)!=lab) continue;
257              //
258              Double_t z3=c3->GetZ();
259          
260              if (z3 > (zr3+kzWin/2)) break;  //check in Z
261              Float_t xyz3[3]; c3->GetGlobalXYZ(xyz3);
262              Double_t phi3=TMath::ATan2(xyz3[1],xyz3[0]);
263              if (TMath::Abs(phir3-phi3) > kpWin/100) continue;  //check in Phi
264
265              AliITSUClusterPix cc(*((AliITSUClusterPix*)c2));
266              cc.GoToFrameTrk();
267              AddCookedSeed(xyz1, kSeedingLayer1, n1,
268                            xyz3, kSeedingLayer3, n3, 
269                            &cc,  kSeedingLayer2, n2);
270
271          }
272      }
273    }
274
275    for (Int_t n1=0; n1<nClusters1; n1++) {
276      AliCluster *c1=layer1.GetCluster(n1);
277      ((AliITSUClusterPix*)c1)->GoToFrameTrk();
278    }
279    for (Int_t n2=0; n2<nClusters2; n2++) {
280      AliCluster *c2=layer2.GetCluster(n2);
281      ((AliITSUClusterPix*)c2)->GoToFrameTrk();
282    }
283    for (Int_t n3=0; n3<nClusters3; n3++) {
284      AliCluster *c3=layer3.GetCluster(n3);
285      ((AliITSUClusterPix*)c3)->GoToFrameTrk();
286    }
287
288    fSeeds->Sort();
289    return fSeeds->GetEntriesFast();
290 }
291
292 Int_t AliITSUTrackerCooked::Clusters2Tracks(AliESDEvent *event) {
293   //--------------------------------------------------------------------
294   // This is the main tracking function
295   // The clusters must already be loaded
296   //--------------------------------------------------------------------
297
298   Int_t nSeeds=MakeSeeds();
299
300   // Possibly, icrement the seeds with additional clusters (Kalman)
301
302   // Possibly, (re)fit the found tracks 
303
304   Int_t ngood=0;
305   for (Int_t s=0; s<nSeeds; s++) {
306       const AliITSUTrackCooked *track=(AliITSUTrackCooked*)fSeeds->At(s);
307       ResetTrackToFollow(*track);
308       ResetBestTrack();
309       fI=kSeedingLayer2;
310       fgLayers[fI].ResetTrack(*track);
311
312       for (FollowProlongation(); fI<kSeedingLayer2; fI++) {
313           while (TakeNextProlongation()) FollowProlongation();
314       }
315
316       if (fBestTrack->GetNumberOfClusters() < kminNumberOfClusters) continue;
317
318       CookLabel(fBestTrack,0.); //For comparison only
319       Int_t label=fBestTrack->GetLabel();
320       if (label>0) ngood++;
321
322       AliESDtrack iotrack;
323       iotrack.UpdateTrackParams(fBestTrack,AliESDtrack::kITSin);
324       iotrack.SetLabel(label);
325       event->AddTrack(&iotrack);
326       UseClusters(fBestTrack);
327   }
328
329   Info("Clusters2Tracks","Seeds: %d",nSeeds);
330   if (nSeeds)
331   Info("Clusters2Tracks","Good tracks/seeds: %f",Float_t(ngood)/nSeeds);
332
333   if (fSeeds) {fSeeds->Delete(); delete fSeeds;}
334   fSeeds=0;
335     
336   return 0;
337 }
338
339 void AliITSUTrackerCooked::FollowProlongation() {
340   //--------------------------------------------------------------------
341   // Push this track tree branch towards the primary vertex
342   //--------------------------------------------------------------------
343   while (fI) {
344     fI--;
345     AliITSUlayer &layer = fgLayers[fI]; //fI is the number of the next layer
346     Double_t r=layer.GetR();
347
348     //Find intersection (fTrackToFollow is still at the previous layer)
349     Double_t phi,z;  
350     if (!fTrackToFollow->GetPhiZat(r,phi,z)) {
351       //Warning("FollowProlongation","failed to estimate track !\n");
352       return;
353     }
354
355     //if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
356     Double_t zMin = z - kRoadZ; 
357     Double_t zMax = z + kRoadZ;
358     Double_t phiMin = phi - kRoadY/r;
359     Double_t phiMax = phi + kRoadY/r;
360     if (layer.SelectClusters(zMin, zMax, phiMin, phiMax)==0) return;  
361
362     if (!TakeNextProlongation()) return;
363
364   } 
365
366   //deal with the best track
367   Int_t ncl=fTrackToFollow->GetNumberOfClusters();
368   Int_t nclb=fBestTrack->GetNumberOfClusters();
369   if (ncl)
370   if (ncl >= nclb) {
371      Double_t chi2=fTrackToFollow->GetChi2();
372      if (chi2/ncl < kmaxChi2PerCluster) {        
373         if (ncl > nclb || chi2 < fBestTrack->GetChi2()) {
374            ResetBestTrack();
375         }
376      }
377   }
378
379 }
380
381 Int_t AliITSUTrackerCooked::TakeNextProlongation() {
382   //--------------------------------------------------------------------
383   // Switch to the next track tree branch
384   //--------------------------------------------------------------------
385   AliITSUlayer &layer=fgLayers[fI];
386
387   const AliCluster *c=0; Int_t ci=-1;
388   const AliCluster *cc=0; Int_t cci=-1;
389   UShort_t volId=-1;
390   Double_t x=0., alpha=0.;
391   Double_t z=0., dz=0., y=0., dy=0., chi2=0.; 
392   while ((c=layer.GetNextCluster(ci))!=0) {
393     if (c->IsClusterUsed()) continue;
394     Int_t id=c->GetVolumeId();
395     if (id != volId) {
396        volId=id;
397        Float_t xr,ar; c->GetXAlphaRefPlane(xr, ar);
398        x=xr; alpha=ar;
399        const AliITSUTrackCooked *t = fgLayers[fI+1].GetTrack();
400        ResetTrackToFollow(*t);
401        if (!fTrackToFollow->Propagate(alpha, x, GetBz())) {
402          //Warning("TakeNextProlongation","propagation failed !\n");
403           continue;
404        }
405        dz=7*TMath::Sqrt(fTrackToFollow->GetSigmaZ2() + kSigma2);
406        dy=7*TMath::Sqrt(fTrackToFollow->GetSigmaY2() + kSigma2);
407        z=fTrackToFollow->GetZ();
408        y=fTrackToFollow->GetY();
409     }
410
411     //if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
412
413     if (TMath::Abs(z - c->GetZ()) > dz) continue;
414     if (TMath::Abs(y - c->GetY()) > dy) continue;
415
416     Double_t ch2=fTrackToFollow->GetPredictedChi2(c); 
417     if (ch2 > kmaxChi2PerCluster) continue;
418     chi2=ch2;
419     cc=c; cci=ci;
420     break;
421   }
422
423   if (!cc) return 0;
424
425   if (!fTrackToFollow->Update(cc,chi2,(fI<<28)+cci)) {
426      //Warning("TakeNextProlongation","filtering failed !\n");
427      return 0;
428   }
429   Double_t xx0 = (fI > 2) ? 0.008 : 0.003;  // Rough layer thickness
430   Double_t x0  = 9.36; // Radiation length of Si [cm]
431   Double_t rho = 2.33; // Density of Si [g/cm^3] 
432   Double_t mass = fTrackToFollow->GetMass();
433   fTrackToFollow->CorrectForMeanMaterial(xx0, xx0*x0*rho, mass, kTRUE);
434   layer.ResetTrack(*fTrackToFollow); 
435
436   return 1;
437 }
438
439 Int_t AliITSUTrackerCooked::PropagateBack(AliESDEvent *event) {
440   //--------------------------------------------------------------------
441   // Here, we implement the Kalman smoother ?
442   // The clusters must already be loaded
443   //--------------------------------------------------------------------
444   Int_t n=event->GetNumberOfTracks();
445   Int_t ntrk=0;
446   Int_t ngood=0;
447   for (Int_t i=0; i<n; i++) {
448       AliESDtrack *esdTrack=event->GetTrack(i);
449
450       if ((esdTrack->GetStatus()&AliESDtrack::kITSin)==0) continue;
451
452       AliITSUTrackCooked track(*esdTrack);
453
454       ResetTrackToFollow(track);
455
456       fTrackToFollow->ResetCovariance(10.); fTrackToFollow->ResetClusters();
457       if (RefitAt(40., fTrackToFollow, &track)) {
458
459          CookLabel(fTrackToFollow, 0.); //For comparison only
460          Int_t label=fTrackToFollow->GetLabel();
461          if (label>0) ngood++;
462
463          esdTrack->UpdateTrackParams(fTrackToFollow,AliESDtrack::kITSout);
464          //UseClusters(fTrackToFollow);
465          ntrk++;
466       }
467   }
468
469   Info("PropagateBack","Back propagated tracks: %d",ntrk);
470   if (ntrk)
471   Info("PropagateBack","Good tracks/back propagated: %f",Float_t(ngood)/ntrk);
472     
473   return 0;
474 }
475
476 Bool_t AliITSUTrackerCooked::
477 RefitAt(Double_t xx, AliITSUTrackCooked *t, const AliITSUTrackCooked *c) {
478   //--------------------------------------------------------------------
479   // This function refits the track "t" at the position "x" using
480   // the clusters from "c"
481   //--------------------------------------------------------------------
482   Int_t index[kNLayers];
483   Int_t k;
484   for (k=0; k<kNLayers; k++) index[k]=-1;
485   Int_t nc=c->GetNumberOfClusters();
486   for (k=0; k<nc; k++) {
487     Int_t idx=c->GetClusterIndex(k), nl=(idx&0xf0000000)>>28;
488     index[nl]=idx;
489   }
490
491   Int_t from, to, step;
492   if (xx > t->GetX()) {
493       from=0; to=kNLayers;
494       step=+1;
495   } else {
496       from=kNLayers-1; to=-1;
497       step=-1;
498   }
499
500   for (Int_t i=from; i != to; i += step) {
501      Int_t idx=index[i];
502      if (idx>=0) {
503         const AliCluster *cl=GetCluster(idx);
504         Float_t xr,ar; cl->GetXAlphaRefPlane(xr, ar);
505         if (!t->Propagate(Double_t(ar), Double_t(xr), GetBz())) {
506            //Warning("RefitAt","propagation failed !\n");
507            return kFALSE;
508         }
509         Double_t chi2=t->GetPredictedChi2(cl);
510         if (chi2 < kmaxChi2PerCluster) t->Update(cl, chi2, idx);
511      } else {
512         Double_t r=fgLayers[i].GetR();
513         Double_t phi,z;
514         if (!t->GetPhiZat(r,phi,z)) {
515            //Warning("RefitAt","failed to estimate track !\n");
516            return kFALSE;
517         }
518         if (!t->Propagate(phi, r, GetBz())) {
519            //Warning("RefitAt","propagation failed !\n");
520            return kFALSE;
521         }
522      }
523      Double_t xx0 = (i > 2) ? 0.008 : 0.003;  // Rough layer thickness
524      Double_t x0  = 9.36; // Radiation length of Si [cm]
525      Double_t rho = 2.33; // Density of Si [g/cm^3]
526      Double_t mass = t->GetMass();
527      t->CorrectForMeanMaterial(xx0, -step*xx0*x0*rho, mass, kTRUE);
528   }
529
530   if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
531   return kTRUE;
532 }
533
534 Int_t AliITSUTrackerCooked::RefitInward(AliESDEvent *event) {
535   //--------------------------------------------------------------------
536   // Some final refit, after the outliers get removed by the smoother ?  
537   // The clusters must be loaded
538   //--------------------------------------------------------------------
539   Int_t n=event->GetNumberOfTracks();
540   Int_t ntrk=0;
541   Int_t ngood=0;
542   for (Int_t i=0; i<n; i++) {
543       AliESDtrack *esdTrack=event->GetTrack(i);
544
545       if ((esdTrack->GetStatus()&AliESDtrack::kITSout) == 0) continue;
546
547       AliITSUTrackCooked track(*esdTrack);
548       ResetTrackToFollow(track);
549
550       fTrackToFollow->ResetCovariance(10.); fTrackToFollow->ResetClusters();
551       if (!RefitAt(2.1, fTrackToFollow, &track)) continue;
552       //Cross the beam pipe
553       if (!fTrackToFollow->PropagateTo(1.8, 2.27e-3, 35.28*1.848)) continue;
554
555       CookLabel(fTrackToFollow, 0.); //For comparison only
556       Int_t label=fTrackToFollow->GetLabel();
557       if (label>0) ngood++;
558
559       esdTrack->UpdateTrackParams(fTrackToFollow,AliESDtrack::kITSrefit);
560       //esdTrack->RelateToVertex(event->GetVertex(),GetBz(),33.);
561       //UseClusters(fTrackToFollow);
562       ntrk++;
563   }
564
565   Info("RefitInward","Refitted tracks: %d",ntrk);
566   if (ntrk)
567   Info("RefitInward","Good tracks/refitted: %f",Float_t(ngood)/ntrk);
568     
569   return 0;
570 }
571
572 Int_t AliITSUTrackerCooked::LoadClusters(TTree *cTree) {
573   //--------------------------------------------------------------------
574   // This function reads the ITSU clusters from the tree,
575   // sort them, distribute over the internal tracker arrays, etc
576   //--------------------------------------------------------------------
577   if (!cTree) {
578      AliFatal("No cluster tree !");
579      return 1;
580   }
581
582   //This TClonesArray is not the owner of the clusters
583   TClonesArray dummy("AliITSUClusterPix",kMaxClusterPerLayer), *clusters=&dummy;
584
585   for (Int_t i=0; i<kNLayers; i++) {
586       TBranch *br = cTree->GetBranch(Form("ITSRecPoints%d",i));
587       if (!br) AliFatal(Form("No cluster branch for layer %d",i));
588       br->SetAddress(&clusters);
589       br->GetEvent(0);
590       switch (i) {
591       case kSeedingLayer1: 
592       case kSeedingLayer2: 
593       case kSeedingLayer3: 
594          fgLayers[i].InsertClusters(clusters,kTRUE);
595          break;
596       default:
597          fgLayers[i].InsertClusters(clusters,kFALSE);
598          break;
599       }
600       clusters->Delete();
601   }
602
603   return 0;
604 }
605
606 void AliITSUTrackerCooked::UnloadClusters() {
607   //--------------------------------------------------------------------
608   // This function unloads ITSU clusters from the RAM
609   //--------------------------------------------------------------------
610   for (Int_t i=0; i<kNLayers; i++) fgLayers[i].DeleteClusters();
611 }
612
613 AliCluster *AliITSUTrackerCooked::GetCluster(Int_t index) const {
614   //--------------------------------------------------------------------
615   //       Return pointer to a given cluster
616   //--------------------------------------------------------------------
617     Int_t l=(index & 0xf0000000) >> 28;
618     Int_t c=(index & 0x0fffffff) >> 00;
619     return fgLayers[l].GetCluster(c);
620 }
621
622 AliITSUTrackerCooked::AliITSUlayer::AliITSUlayer():
623   fR(0),
624   fN(0),
625   fNsel(0),
626   fTrack(0) 
627 {
628   //--------------------------------------------------------------------
629   // This default constructor needs to be provided
630   //--------------------------------------------------------------------
631   for (Int_t i=0; i<kMaxClusterPerLayer; i++) fClusters[i]=0;
632   for (Int_t i=0; i<kMaxSelected; i++) fIndex[i]=-1;
633 }
634
635 AliITSUTrackerCooked::AliITSUlayer::~AliITSUlayer()
636 {
637   //--------------------------------------------------------------------
638   // Simple destructor
639   //--------------------------------------------------------------------
640   DeleteClusters();
641   delete fTrack;
642 }
643
644 void 
645 AliITSUTrackerCooked::AliITSUlayer::ResetTrack(const AliITSUTrackCooked &t) {
646   //--------------------------------------------------------------------
647   // Replace the track estimate at this layer
648   //--------------------------------------------------------------------
649    delete fTrack;
650    fTrack=new AliITSUTrackCooked(t);
651 }
652
653 void AliITSUTrackerCooked::
654   AliITSUlayer::InsertClusters(TClonesArray *clusters, Bool_t seedingLayer)
655 {
656   //--------------------------------------------------------------------
657   // Load clusters to this layer
658   //--------------------------------------------------------------------
659   Int_t ncl=clusters->GetEntriesFast();
660  
661   while (ncl--) {
662      AliITSUClusterPix *c=(AliITSUClusterPix*)clusters->UncheckedAt(ncl);
663      (seedingLayer) ? c->GoToFrameGlo() : c->GoToFrameTrk();
664      //if (!c->Misalign()) AliWarning("Can't misalign this cluster !");
665      InsertCluster(new AliITSUClusterPix(*c));
666   }
667
668 }
669
670 void AliITSUTrackerCooked::AliITSUlayer::DeleteClusters()
671 {
672   //--------------------------------------------------------------------
673   // Load clusters to this layer
674   //--------------------------------------------------------------------
675   for (Int_t i=0; i<fN; i++) {delete fClusters[i]; fClusters[i]=0;}
676   fN=0;
677 }
678
679 Int_t 
680 AliITSUTrackerCooked::AliITSUlayer::InsertCluster(AliCluster *c) {
681   //--------------------------------------------------------------------
682   // This function inserts a cluster to this layer in increasing
683   // order of the cluster's fZ
684   //--------------------------------------------------------------------
685   if (fN>=kMaxClusterPerLayer) {
686      ::Error("InsertCluster","Too many clusters !\n");
687      return 1;
688   }
689   if (fN==0) fClusters[0]=c;
690   else {
691      Int_t i=FindClusterIndex(c->GetZ());
692      Int_t k=fN-i;
693      memmove(fClusters+i+1 ,fClusters+i,k*sizeof(AliCluster*));
694      fClusters[i]=c;
695   }
696   fN++;
697   return 0;
698 }
699
700 Int_t 
701 AliITSUTrackerCooked::AliITSUlayer::FindClusterIndex(Double_t z) const {
702   //--------------------------------------------------------------------
703   // This function returns the index of the first 
704   // with its fZ >= "z". 
705   //--------------------------------------------------------------------
706   if (fN==0) return 0;
707
708   Int_t b=0;
709   if (z <= fClusters[b]->GetZ()) return b;
710
711   Int_t e=b+fN-1;
712   if (z > fClusters[e]->GetZ()) return e+1;
713
714   Int_t m=(b+e)/2;
715   for (; b<e; m=(b+e)/2) {
716     if (z > fClusters[m]->GetZ()) b=m+1;
717     else e=m; 
718   }
719   return m;
720 }
721
722 Int_t AliITSUTrackerCooked::AliITSUlayer::
723 SelectClusters(Float_t zMin,Float_t zMax,Float_t phiMin, Float_t phiMax) {
724   //--------------------------------------------------------------------
725   // This function selects clusters within the "road"
726   //--------------------------------------------------------------------
727   UShort_t volId=-1;
728   Float_t x=0., alpha=0.;
729   for (Int_t i=FindClusterIndex(zMin); i<fN; i++) {
730       AliCluster *c=fClusters[i];
731       if (c->GetZ() > zMax) break;
732       if (c->IsClusterUsed()) continue;
733       UShort_t id=c->GetVolumeId();
734       if (id != volId) {
735         volId=id;
736         c->GetXAlphaRefPlane(x,alpha); //FIXME
737       }
738       Double_t cPhi=alpha + c->GetY()/fR;
739       if (cPhi<0.) cPhi+=2.*TMath::Pi();
740       else if (cPhi >= 2.*TMath::Pi()) cPhi-=2.*TMath::Pi();
741       if (cPhi <= phiMin) continue;
742       if (cPhi >  phiMax) continue;
743       fIndex[fNsel++]=i;
744       if (fNsel==kMaxSelected) break;
745   }
746   return fNsel;
747 }
748
749 const AliCluster *AliITSUTrackerCooked::AliITSUlayer::GetNextCluster(Int_t &ci){
750   //--------------------------------------------------------------------
751   // This function returns clusters within the "road" 
752   //--------------------------------------------------------------------
753   AliCluster *c=0;
754   ci=-1;
755   if (fNsel) {
756      fNsel--;
757      ci=fIndex[fNsel]; 
758      c=fClusters[ci];
759   }
760   return c; 
761 }
762