]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
Bug fix - chack the abs values
[u/mrichter/AliRoot.git] / ITS / AliITStrackerV2.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 /* $Id$ */
16 //-------------------------------------------------------------------------
17 //               Implementation of the ITS tracker class
18 //    It reads AliITSRecPoint clusters and creates AliITStrackV2 tracks
19 //                   and fills with them the ESD
20 //          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
21 //     dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
22 //-------------------------------------------------------------------------
23
24 #include <new>
25
26 #include <TError.h>
27 #include <TFile.h>
28 #include <TTree.h>
29 #include <TRandom.h>
30 #include <TGeoMatrix.h>
31
32 #include "AliITSgeomTGeo.h"
33 #include "AliAlignObj.h"
34 #include "AliITSRecPoint.h"
35 #include "AliESDEvent.h"
36 #include "AliESDtrack.h"
37 #include "AliITSRecPoint.h"
38 #include "AliITSReconstructor.h"
39 #include "AliITStrackerV2.h"
40
41 ClassImp(AliITStrackerV2)
42
43 AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[AliITSgeomTGeo::kNLayers]; //ITS layers
44
45 AliITStrackerV2::AliITStrackerV2(): 
46   AliTracker(), 
47   fI(AliITSgeomTGeo::GetNLayers()),
48   fBestTrack(),
49   fTrackToFollow(),
50   fPass(0),
51   fLastLayerToTrackTo(AliITSRecoParam::GetLastLayerToTrackTo())
52 {
53   //--------------------------------------------------------------------
54   //This is the AliITStrackerV2 default constructor
55   //--------------------------------------------------------------------
56
57   for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) new(fgLayers+i-1) AliITSlayer();
58
59   fConstraint[0]=1; fConstraint[1]=0;
60
61   Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
62                   AliITSReconstructor::GetRecoParam()->GetYVdef(),
63                   AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
64   Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
65                   AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
66                   AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
67   SetVertex(xyz,ers);
68
69   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=AliITSRecoParam::GetLayersNotToSkip(i);
70
71 }
72
73 AliITStrackerV2::AliITStrackerV2(const AliITStrackerV2 &t): 
74   AliTracker(t), 
75   fI(t.fI),
76   fBestTrack(t.fBestTrack),
77   fTrackToFollow(t.fTrackToFollow),
78   fPass(t.fPass),
79   fLastLayerToTrackTo(t.fLastLayerToTrackTo)
80 {
81   //--------------------------------------------------------------------
82   //This is the AliITStrackerV2 copy constructor
83   //--------------------------------------------------------------------
84
85   //for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) new(fgLayers+i-1) AliITSlayer();
86
87   fConstraint[0]=t.fConstraint[0]; fConstraint[1]=t.fConstraint[1];
88
89   Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
90                   AliITSReconstructor::GetRecoParam()->GetYVdef(),
91                   AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
92   Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
93                   AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
94                   AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
95   xyz[0]=t.GetX(); xyz[1]=t.GetY(); xyz[2]=t.GetZ(); 
96   ers[0]=t.GetSigmaX(); ers[1]=t.GetSigmaY(); ers[2]=t.GetSigmaZ(); 
97   SetVertex(xyz,ers);
98
99   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=t.fLayersNotToSkip[i];
100
101 }
102
103 AliITStrackerV2::AliITStrackerV2(const Char_t *geom) : 
104   AliTracker(), 
105   fI(AliITSgeomTGeo::GetNLayers()),
106   fBestTrack(),
107   fTrackToFollow(),
108   fPass(0),
109   fLastLayerToTrackTo(AliITSRecoParam::GetLastLayerToTrackTo())
110 {
111   //--------------------------------------------------------------------
112   //This is the AliITStrackerV2 constructor
113   //--------------------------------------------------------------------
114   if (geom) {
115     AliWarning("\"geom\" is actually a dummy argument !");
116   }
117
118   for (Int_t i=1; i<AliITSgeomTGeo::GetNLayers()+1; i++) {
119     Int_t nlad=AliITSgeomTGeo::GetNLadders(i);
120     Int_t ndet=AliITSgeomTGeo::GetNDetectors(i);
121
122     Double_t xyz[3], &x=xyz[0], &y=xyz[1], &z=xyz[2];
123     AliITSgeomTGeo::GetOrigTranslation(i,1,1,xyz); 
124     Double_t poff=TMath::ATan2(y,x);
125     Double_t zoff=z;
126     Double_t r=TMath::Sqrt(x*x + y*y);
127
128     AliITSgeomTGeo::GetOrigTranslation(i,1,2,xyz);
129     r += TMath::Sqrt(x*x + y*y);
130     AliITSgeomTGeo::GetOrigTranslation(i,2,1,xyz);
131     r += TMath::Sqrt(x*x + y*y);
132     AliITSgeomTGeo::GetOrigTranslation(i,2,2,xyz);
133     r += TMath::Sqrt(x*x + y*y);
134     r*=0.25;
135
136     new (fgLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
137
138     for (Int_t j=1; j<nlad+1; j++) {
139       for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
140         TGeoHMatrix m; AliITSgeomTGeo::GetOrigMatrix(i,j,k,m);
141         const TGeoHMatrix *tm=AliITSgeomTGeo::GetTracking2LocalMatrix(i,j,k);
142         m.Multiply(tm);
143         Double_t txyz[3]={0.}; 
144         xyz[0]=0.; xyz[1]=0.; xyz[2]=0.;
145         m.LocalToMaster(txyz,xyz);
146         r=TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]);
147         Double_t phi=TMath::ATan2(xyz[1],xyz[0]);
148
149         if (phi<0) phi+=TMath::TwoPi();
150         else if (phi>=TMath::TwoPi()) phi-=TMath::TwoPi();
151
152         AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1); 
153         new(&det) AliITSdetector(r,phi); 
154       } 
155     }  
156
157   }
158
159   fConstraint[0]=1; fConstraint[1]=0;
160
161   Double_t xyz[]={AliITSReconstructor::GetRecoParam()->GetXVdef(),
162                   AliITSReconstructor::GetRecoParam()->GetYVdef(),
163                   AliITSReconstructor::GetRecoParam()->GetZVdef()}; 
164   Double_t ers[]={AliITSReconstructor::GetRecoParam()->GetSigmaXVdef(),
165                   AliITSReconstructor::GetRecoParam()->GetSigmaYVdef(),
166                   AliITSReconstructor::GetRecoParam()->GetSigmaZVdef()}; 
167   SetVertex(xyz,ers);
168
169   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=AliITSRecoParam::GetLayersNotToSkip(i);
170
171 }
172
173 void AliITStrackerV2::SetLayersNotToSkip(Int_t *l) {
174   //--------------------------------------------------------------------
175   //This function set masks of the layers which must be not skipped
176   //--------------------------------------------------------------------
177   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fLayersNotToSkip[i]=l[i];
178 }
179
180 Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
181   //--------------------------------------------------------------------
182   //This function loads ITS clusters
183   //--------------------------------------------------------------------
184   TBranch *branch=cTree->GetBranch("ITSRecPoints");
185   if (!branch) { 
186     Error("LoadClusters"," can't get the branch !\n");
187     return 1;
188   }
189
190   TClonesArray dummy("AliITSRecPoint",10000), *clusters=&dummy;
191   branch->SetAddress(&clusters);
192
193   Int_t j=0;
194   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) {
195     Int_t ndet=fgLayers[i].GetNdetectors();
196     Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
197
198     Double_t r=fgLayers[i].GetR();
199     Double_t circ=TMath::TwoPi()*r;
200
201     for (; j<jmax; j++) {           
202       if (!cTree->GetEvent(j)) continue;
203       Int_t ncl=clusters->GetEntriesFast();
204  
205       while (ncl--) {
206         AliITSRecPoint *c=(AliITSRecPoint*)clusters->UncheckedAt(ncl);
207
208         if (!c->Misalign()) AliWarning("Can't misalign this cluster !");
209
210         Int_t idx=c->GetDetectorIndex();
211         AliITSdetector &det=fgLayers[i].GetDetector(idx);
212    
213         Double_t y=r*det.GetPhi()+c->GetY();
214         if (y>circ) y-=circ; else if (y<0) y+=circ;
215         c->SetPhiR(y);
216
217         fgLayers[i].InsertCluster(new AliITSRecPoint(*c));
218       }
219       clusters->Delete();
220     }
221     fgLayers[i].ResetRoad(); //road defined by the cluster density
222   }
223
224   return 0;
225 }
226
227 void AliITStrackerV2::UnloadClusters() {
228   //--------------------------------------------------------------------
229   //This function unloads ITS clusters
230   //--------------------------------------------------------------------
231   for (Int_t i=0; i<AliITSgeomTGeo::GetNLayers(); i++) fgLayers[i].ResetClusters();
232 }
233
234 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
235   //--------------------------------------------------------------------
236   // Correction for the material between the TPC and the ITS
237   // (should it belong to the TPC code ?)
238   //--------------------------------------------------------------------
239   Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ? 
240   Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
241   Double_t yr=12.8, dr=0.03; // rods ?
242   Double_t zm=0.2, dm=0.40;  // membrane
243   //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
244   Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
245
246   if (t->GetX() > riw) {
247      if (!t->PropagateTo(riw,diw,x0iw)) return 1;
248      if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
249      if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
250      if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
251      //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
252      //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r); 
253      if (!t->PropagateTo(rs,ds)) return 1;
254   } else if (t->GetX() < rs) {
255      if (!t->PropagateTo(rs,-ds)) return 1;
256      //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
257      //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r); 
258      if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
259      if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
260   } else {
261   ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
262     return 1;
263   }
264   
265   return 0;
266 }
267
268 Int_t AliITStrackerV2::Clusters2Tracks(AliESDEvent *event) {
269   //--------------------------------------------------------------------
270   // This functions reconstructs ITS tracks
271   // The clusters must be already loaded !
272   //--------------------------------------------------------------------
273   TObjArray itsTracks(15000);
274
275   {/* Read ESD tracks */
276     Int_t nentr=event->GetNumberOfTracks();
277     Info("Clusters2Tracks", "Number of ESD tracks: %d\n", nentr);
278     while (nentr--) {
279       AliESDtrack *esd=event->GetTrack(nentr);
280
281       if ((esd->GetStatus()&AliESDtrack::kTPCin)==0) continue;
282       if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
283       if (esd->GetStatus()&AliESDtrack::kITSin) continue;
284
285       AliITStrackV2 *t = new AliITStrackV2(*esd);
286
287       if (TMath::Abs(t->GetD(GetX(),GetY()))>4) {
288         delete t;
289         continue;
290       }
291
292       if (CorrectForDeadZoneMaterial(t)!=0) {
293          Warning("Clusters2Tracks",
294                  "failed to correct for the material in the dead zone !\n");
295          delete t;
296          continue;
297       }
298       itsTracks.AddLast(t);
299     }
300   } /* End Read ESD tracks */
301
302   itsTracks.Sort();
303   Int_t nentr=itsTracks.GetEntriesFast();
304
305   Int_t ntrk=0;
306   for (fPass=0; fPass<2; fPass++) {
307      Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
308      for (Int_t i=0; i<nentr; i++) {
309        AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
310        if (t==0) continue;           //this track has been already tracked
311        Int_t tpcLabel=t->GetLabel(); //save the TPC track label
312
313        ResetTrackToFollow(*t);
314        ResetBestTrack();
315
316        for (FollowProlongation(); fI<AliITSgeomTGeo::GetNLayers(); fI++) {
317           while (TakeNextProlongation()) FollowProlongation();
318        }
319
320        if (fBestTrack.GetNumberOfClusters() == 0) continue;
321
322        if (fConstraint[fPass]) {
323           ResetTrackToFollow(*t);
324           if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
325           ResetBestTrack();
326        }
327
328        fBestTrack.SetLabel(tpcLabel);
329        fBestTrack.CookdEdx();
330        CookLabel(&fBestTrack,0.); //For comparison only
331        fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
332        UseClusters(&fBestTrack);
333        delete itsTracks.RemoveAt(i);
334        ntrk++;
335      }
336   }
337
338   itsTracks.Delete();
339
340   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
341
342   return 0;
343 }
344
345 Int_t AliITStrackerV2::PropagateBack(AliESDEvent *event) {
346   //--------------------------------------------------------------------
347   // This functions propagates reconstructed ITS tracks back
348   // The clusters must be loaded !
349   //--------------------------------------------------------------------
350   Int_t nentr=event->GetNumberOfTracks();
351   Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
352
353   Int_t ntrk=0;
354   for (Int_t i=0; i<nentr; i++) {
355      AliESDtrack *esd=event->GetTrack(i);
356
357      if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
358      if (esd->GetStatus()&AliESDtrack::kITSout) continue;
359
360      AliITStrackV2 *t = new AliITStrackV2(*esd);
361
362      ResetTrackToFollow(*t);
363
364      // propagete to vertex [SR, GSI 17.02.2003]
365      // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
366      if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
367        if (fTrackToFollow.PropagateToVertex(event->GetVertex())) {
368           fTrackToFollow.StartTimeIntegral();
369        }
370        fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
371      }
372
373      fTrackToFollow.ResetCovariance(10.); fTrackToFollow.ResetClusters();
374      if (RefitAt(49.,&fTrackToFollow,t)) {
375         if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
376           Warning("PropagateBack",
377                   "failed to correct for the material in the dead zone !\n");
378           delete t;
379           continue;
380         }
381         fTrackToFollow.SetLabel(t->GetLabel());
382         //fTrackToFollow.CookdEdx();
383         CookLabel(&fTrackToFollow,0.); //For comparison only
384         fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
385         UseClusters(&fTrackToFollow);
386         ntrk++;
387      }
388      delete t;
389   }
390
391   Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
392
393   return 0;
394 }
395
396 Int_t AliITStrackerV2::RefitInward(AliESDEvent *event) {
397   //--------------------------------------------------------------------
398   // This functions refits ITS tracks using the 
399   // "inward propagated" TPC tracks
400   // The clusters must be loaded !
401   //--------------------------------------------------------------------
402   Int_t nentr=event->GetNumberOfTracks();
403   Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
404
405   Int_t ntrk=0;
406   for (Int_t i=0; i<nentr; i++) {
407     AliESDtrack *esd=event->GetTrack(i);
408
409     if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
410     if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
411     if (esd->GetStatus()&AliESDtrack::kTPCout)
412     if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
413
414     AliITStrackV2 *t = new AliITStrackV2(*esd);
415
416     if (CorrectForDeadZoneMaterial(t)!=0) {
417        Warning("RefitInward",
418                "failed to correct for the material in the dead zone !\n");
419        delete t;
420        continue;
421     }
422
423     ResetTrackToFollow(*t);
424     fTrackToFollow.ResetClusters();
425
426     //Refitting...
427     if (RefitAt(3.7, &fTrackToFollow, t, kTRUE)) {
428        fTrackToFollow.SetLabel(t->GetLabel());
429        fTrackToFollow.CookdEdx();
430        CookLabel(&fTrackToFollow,0.); //For comparison only
431
432        if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe 
433          fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
434          AliESDtrack  *esdTrack =fTrackToFollow.GetESDtrack();
435          Double_t r[3]={0.,0.,0.};
436          Double_t maxD=3.;
437          esdTrack->RelateToVertex(event->GetVertex(),GetBz(r),maxD);
438          ntrk++;
439        }
440     }
441     delete t;
442   }
443
444   Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
445
446   return 0;
447 }
448
449 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
450   //--------------------------------------------------------------------
451   //       Return pointer to a given cluster
452   //--------------------------------------------------------------------
453   Int_t l=(index & 0xf0000000) >> 28;
454   Int_t c=(index & 0x0fffffff) >> 00;
455   return fgLayers[l].GetCluster(c);
456 }
457
458
459 void AliITStrackerV2::FollowProlongation() {
460   //--------------------------------------------------------------------
461   //This function finds a track prolongation 
462   //--------------------------------------------------------------------
463   while (fI>fLastLayerToTrackTo) {
464     Int_t i=fI-1;
465
466     AliITSlayer &layer=fgLayers[i];
467     AliITStrackV2 &track=fTracks[i];
468
469     Double_t r=layer.GetR();
470
471     if (i==3 || i==1) {
472        Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
473        Double_t d=0.0034, x0=38.6;
474        if (i==1) {rs=9.; d=0.0097; x0=42;}
475        if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
476          //Warning("FollowProlongation","propagation failed !\n");
477          return;
478        }
479     }
480
481     //find intersection
482     Double_t phi,z;  
483     if (!fTrackToFollow.GetPhiZat(r,phi,z)) {
484       //Warning("FollowProlongation","failed to estimate track !\n");
485       return;
486     }
487
488     Int_t idet=layer.FindDetectorIndex(phi,z);
489     if (idet<0) {
490       //Warning("FollowProlongation","failed to find a detector !\n");
491       return;
492     }
493
494     //propagate to the intersection
495     const AliITSdetector &det=layer.GetDetector(idet);
496     phi=det.GetPhi();
497     if (!fTrackToFollow.Propagate(phi,det.GetR())) {
498       //Warning("FollowProlongation","propagation failed !\n");
499       return;
500     }
501     fTrackToFollow.SetDetectorIndex(idet);
502
503     //Select possible prolongations and store the current track estimation
504     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
505     Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
506     Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
507     Double_t road=layer.GetRoad();
508     if (dz*dy>road*road) {
509        Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
510        dz=road*scz; dy=road*scy;
511     } 
512
513     //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
514     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
515     if (dz > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
516       //Warning("FollowProlongation","too broad road in Z !\n");
517       return;
518     }
519
520     if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
521
522     //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
523     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
524     if (dy > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
525       //Warning("FollowProlongation","too broad road in Y !\n");
526       return;
527     }
528
529     fI--;
530
531     Double_t zmin=track.GetZ() - dz; 
532     Double_t zmax=track.GetZ() + dz;
533     Double_t ymin=track.GetY() + r*phi - dy;
534     Double_t ymax=track.GetY() + r*phi + dy;
535     if (layer.SelectClusters(zmin,zmax,ymin,ymax)==0) 
536        if (fLayersNotToSkip[fI]) return;  
537
538     if (!TakeNextProlongation()) 
539        if (fLayersNotToSkip[fI]) return;
540
541   } 
542
543   //deal with the best track
544   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
545   Int_t nclb=fBestTrack.GetNumberOfClusters();
546   if (ncl)
547   if (ncl >= nclb) {
548      Double_t chi2=fTrackToFollow.GetChi2();
549      if (chi2/ncl < AliITSReconstructor::GetRecoParam()->GetChi2PerCluster()) {        
550         if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
551            ResetBestTrack();
552         }
553      }
554   }
555
556 }
557
558 Int_t AliITStrackerV2::TakeNextProlongation() {
559   //--------------------------------------------------------------------
560   // This function takes another track prolongation 
561   //
562   //  dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
563   //--------------------------------------------------------------------
564   AliITSlayer &layer=fgLayers[fI];
565   ResetTrackToFollow(fTracks[fI]);
566
567   Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(fI));
568   Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(fI));
569   Double_t road=layer.GetRoad();
570   if (dz*dy>road*road) {
571      Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
572      dz=road*scz; dy=road*scy;
573   } 
574
575   const AliITSRecPoint *c=0; Int_t ci=-1;
576   const AliITSRecPoint *cc=0; Int_t cci=-1;
577   Double_t chi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
578   while ((c=layer.GetNextCluster(ci))!=0) {
579     Int_t idet=c->GetDetectorIndex();
580
581     if (fTrackToFollow.GetDetectorIndex()!=idet) {
582        const AliITSdetector &det=layer.GetDetector(idet);
583        ResetTrackToFollow(fTracks[fI]);
584        if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
585          //Warning("TakeNextProlongation","propagation failed !\n");
586          continue;
587        }
588        fTrackToFollow.SetDetectorIndex(idet);
589        if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
590     }
591
592     if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
593     if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
594
595     Double_t ch2=fTrackToFollow.GetPredictedChi2(c); 
596     if (ch2 > chi2) continue;
597     chi2=ch2;
598     cc=c; cci=ci;
599     break;
600   }
601
602   if (!cc) return 0;
603
604   {// Take into account the mis-alignment
605     Double_t x = fTrackToFollow.GetX() + cc->GetX();
606     if (!fTrackToFollow.PropagateTo(x,0.,0.)) return 0;
607   }
608   if (!fTrackToFollow.Update(cc,chi2,(fI<<28)+cci)) {
609      //Warning("TakeNextProlongation","filtering failed !\n");
610      return 0;
611   }
612
613   if (fTrackToFollow.GetNumberOfClusters()>1)
614     if (TMath::Abs(fTrackToFollow.GetD(GetX(),GetY()))>4) return 0;
615
616   fTrackToFollow.
617     SetSampledEdx(cc->GetQ(),fI-2); //b.b.
618
619   {
620   Double_t x0;
621  Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
622   fTrackToFollow.CorrectForMaterial(d,x0);
623   }
624
625   if (fConstraint[fPass]) {
626     Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
627     Double_t xyz[]={GetX(),GetY(),GetZ()};
628     Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
629     fTrackToFollow.Improve(d,xyz,ers);
630   }
631
632   return 1;
633 }
634
635
636 AliITStrackerV2::AliITSlayer::AliITSlayer():
637   fR(0.),
638   fPhiOffset(0.),
639   fNladders(0),
640   fZOffset(0.),
641   fNdetectors(0),
642   fDetectors(0),
643   fNsel(0),
644   fRoad(2*fR*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
645 {
646   //--------------------------------------------------------------------
647   //default AliITSlayer constructor
648   //--------------------------------------------------------------------
649   
650   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
651   for (Int_t i=0; i<AliITSRecoParam::fgkMaxClusterPerLayer; i++){
652     fClusters[i]=0;
653     fIndex[i]=0;
654   }
655 }
656
657 AliITStrackerV2::AliITSlayer::
658 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd): 
659   fR(r), 
660   fPhiOffset(p), 
661   fNladders(nl),
662   fZOffset(z),
663   fNdetectors(nd),
664   fDetectors(new AliITSdetector[nl*nd]),
665   fNsel(0),
666   fRoad(2*r*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
667 {
668   //--------------------------------------------------------------------
669   //main AliITSlayer constructor
670   //--------------------------------------------------------------------
671
672   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
673
674   for (Int_t i=0; i<AliITSRecoParam::fgkMaxClusterPerLayer; i++){
675     fClusters[i]=0;
676     fIndex[i]=0;
677   }
678 }
679
680 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
681   //--------------------------------------------------------------------
682   // AliITSlayer destructor
683   //--------------------------------------------------------------------
684   delete[] fDetectors;
685   ResetClusters();
686 }
687
688 void AliITStrackerV2::AliITSlayer::ResetClusters() {
689   //--------------------------------------------------------------------
690   // This function removes loaded clusters
691   //--------------------------------------------------------------------
692    for (Int_t s=0; s<kNsector; s++) {
693        Int_t &n=fN[s];
694        while (n) {
695           n--;
696           delete fClusters[s*kMaxClusterPerSector+n];
697        }
698    }
699 }
700
701 void AliITStrackerV2::AliITSlayer::ResetRoad() {
702   //--------------------------------------------------------------------
703   // This function calculates the road defined by the cluster density
704   //--------------------------------------------------------------------
705   Int_t n=0;
706   for (Int_t s=0; s<kNsector; s++) {
707     Int_t i=fN[s];
708     while (i--) 
709        if (TMath::Abs(fClusters[s*kMaxClusterPerSector+i]->GetZ())<fR) n++;
710   }
711   if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
712 }
713
714 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSRecPoint *c) {
715   //--------------------------------------------------------------------
716   // This function inserts a cluster to this layer in increasing
717   // order of the cluster's fZ
718   //--------------------------------------------------------------------
719   Float_t circ=TMath::TwoPi()*fR;
720   Int_t sec=Int_t(kNsector*c->GetPhiR()/circ);
721   if (sec>=kNsector) {
722      ::Error("InsertCluster","Wrong sector !\n");
723      return 1;
724   }
725   Int_t &n=fN[sec];
726   if (n>=kMaxClusterPerSector) {
727      ::Error("InsertCluster","Too many clusters !\n");
728      return 1;
729   }
730   if (n==0) fClusters[sec*kMaxClusterPerSector]=c;
731   else {
732      Int_t i=FindClusterIndex(c->GetZ(),sec);
733      Int_t k=n-i+sec*kMaxClusterPerSector;
734      memmove(fClusters+i+1 ,fClusters+i,k*sizeof(AliITSRecPoint*));
735      fClusters[i]=c;
736   }
737   n++;
738   return 0;
739 }
740
741 Int_t 
742 AliITStrackerV2::AliITSlayer::FindClusterIndex(Float_t z,Int_t s) const {
743   //--------------------------------------------------------------------
744   // For the sector "s", this function returns the index of the first 
745   // with its fZ >= "z". 
746   //--------------------------------------------------------------------
747   Int_t nc=fN[s];
748   if (nc==0) return kMaxClusterPerSector*s;
749
750   Int_t b=kMaxClusterPerSector*s;
751   if (z <= fClusters[b]->GetZ()) return b;
752
753   Int_t e=b+nc-1;
754   if (z > fClusters[e]->GetZ()) return e+1;
755
756   Int_t m=(b+e)/2;
757   for (; b<e; m=(b+e)/2) {
758     if (z > fClusters[m]->GetZ()) b=m+1;
759     else e=m; 
760   }
761   return m;
762 }
763
764 Int_t AliITStrackerV2::AliITSlayer::
765 SelectClusters(Float_t zmin,Float_t zmax,Float_t ymin, Float_t ymax) {
766   //--------------------------------------------------------------------
767   // This function selects clusters within the "window"
768   //--------------------------------------------------------------------
769     Float_t circ=fR*TMath::TwoPi();
770
771     if (ymin>circ) ymin-=circ; else if (ymin<0) ymin+=circ;
772     if (ymax>circ) ymax-=circ; else if (ymax<0) ymax+=circ;
773
774     Int_t i1=Int_t(kNsector*ymin/circ); if (i1==kNsector) i1--;
775     if (fN[i1]!=0) {
776        Float_t ym = (ymax<ymin) ? ymax+circ : ymax;
777        Int_t i=FindClusterIndex(zmin,i1), imax=i1*kMaxClusterPerSector+fN[i1];
778        for (; i<imax; i++) {
779            AliITSRecPoint *c=fClusters[i];
780            if (c->IsUsed()) continue;
781            if (c->GetZ()>zmax) break;
782            if (c->GetPhiR()<=ymin) continue;
783            if (c->GetPhiR()>ym) continue;
784            fIndex[fNsel++]=i;
785        }
786     }
787
788     Int_t i2=Int_t(kNsector*ymax/circ); if (i2==kNsector) i2--;
789     if (i2==i1) return fNsel;
790
791     if (fN[i2]!=0) {
792        Float_t ym = (ymin>ymax) ? ymin-circ : ymin;
793        Int_t i=FindClusterIndex(zmin,i2), imax=i2*kMaxClusterPerSector+fN[i2];
794        for (; i<imax; i++) {
795            AliITSRecPoint *c=fClusters[i];
796            if (c->IsUsed()) continue;
797            if (c->GetZ()>zmax) break;
798            if (c->GetPhiR()<=ym) continue;
799            if (c->GetPhiR()>ymax) continue;
800            fIndex[fNsel++]=i;
801        }
802     }
803
804     return fNsel;
805 }
806
807 const AliITSRecPoint *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
808   //--------------------------------------------------------------------
809   // This function returns clusters within the "window" 
810   //--------------------------------------------------------------------
811   AliITSRecPoint *c=0;
812   ci=-1;
813   if (fNsel) {
814      fNsel--;
815      ci=fIndex[fNsel]; 
816      c=fClusters[ci];
817   }
818   return c; 
819 }
820
821 Int_t AliITStrackerV2::AliITSlayer::GetNumberOfClusters() const {
822   Int_t n=0;
823   for (Int_t s=0; s<kNsector; s++) n+=fN[s];
824   return n; 
825 }
826
827 Int_t 
828 AliITStrackerV2::AliITSlayer::FindDetectorIndex(Double_t phi,Double_t z)const {
829   //--------------------------------------------------------------------
830   //This function finds the detector crossed by the track
831   //--------------------------------------------------------------------
832   Double_t dphi;
833   if (fZOffset<0)            // old geometry
834     dphi = -(phi-fPhiOffset);
835   else                       // new geometry
836     dphi = phi-fPhiOffset;
837
838   if      (dphi <  0) dphi += 2*TMath::Pi();
839   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
840   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
841   if (np>=fNladders) np-=fNladders;
842   if (np<0)          np+=fNladders;
843
844   Double_t dz=fZOffset-z;
845   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
846   if (nz>=fNdetectors) return -1;
847   if (nz<0)            return -1;
848
849   return np*fNdetectors + nz;
850 }
851
852 Double_t 
853 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
854 const {
855   //--------------------------------------------------------------------
856   //This function returns the layer thickness at this point (units X0)
857   //--------------------------------------------------------------------
858   Double_t d=0.0085;
859   x0=21.82;
860
861   if (43<fR&&fR<45) { //SSD2
862      Double_t dd=0.0034;
863      d=dd;
864      if (TMath::Abs(y-0.00)>3.40) d+=dd;
865      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
866      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
867      for (Int_t i=0; i<12; i++) {
868        if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
869           if (TMath::Abs(y-0.00)>3.40) d+=dd;
870           d+=0.0034; 
871           break;
872        }
873        if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
874           if (TMath::Abs(y-0.00)>3.40) d+=dd;
875           d+=0.0034; 
876           break;
877        }         
878        if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
879        if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
880      }
881   } else 
882   if (37<fR&&fR<41) { //SSD1
883      Double_t dd=0.0034;
884      d=dd;
885      if (TMath::Abs(y-0.00)>3.40) d+=dd;
886      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
887      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
888      for (Int_t i=0; i<11; i++) {
889        if (TMath::Abs(z-3.9*i)<0.15) {
890           if (TMath::Abs(y-0.00)>3.40) d+=dd;
891           d+=dd; 
892           break;
893        }
894        if (TMath::Abs(z+3.9*i)<0.15) {
895           if (TMath::Abs(y-0.00)>3.40) d+=dd;
896           d+=dd; 
897           break;
898        }         
899        if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
900        if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}         
901      }
902   } else
903   if (13<fR&&fR<26) { //SDD
904      Double_t dd=0.0033;
905      d=dd;
906      if (TMath::Abs(y-0.00)>3.30) d+=dd;
907
908      if (TMath::Abs(y-1.80)<0.55) {
909         d+=0.016;
910         for (Int_t j=0; j<20; j++) {
911           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
912           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
913         } 
914      }
915      if (TMath::Abs(y+1.80)<0.55) {
916         d+=0.016;
917         for (Int_t j=0; j<20; j++) {
918           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
919           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
920         } 
921      }
922
923      for (Int_t i=0; i<4; i++) {
924        if (TMath::Abs(z-7.3*i)<0.60) {
925           d+=dd;
926           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
927           break;
928        }
929        if (TMath::Abs(z+7.3*i)<0.60) {
930           d+=dd; 
931           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
932           break;
933        }
934      }
935   } else
936   if (6<fR&&fR<8) {   //SPD2
937      Double_t dd=0.0063; x0=21.5;
938      d=dd;
939      if (TMath::Abs(y-3.08)>0.5) d+=dd;
940      //if (TMath::Abs(y-3.08)>0.45) d+=dd;
941      if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
942   } else
943   if (3<fR&&fR<5) {   //SPD1
944      Double_t dd=0.0063; x0=21.5;
945      d=dd;
946      if (TMath::Abs(y+0.21)>0.6) d+=dd;
947      //if (TMath::Abs(y+0.21)>0.45) d+=dd;
948      if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
949   }
950
951   return d;
952 }
953
954 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
955 {
956   //--------------------------------------------------------------------
957   //Returns the thickness between the current layer and the vertex (units X0)
958   //--------------------------------------------------------------------
959   Double_t d=0.0028*3*3; //beam pipe
960   Double_t x0=0;
961
962   Double_t xn=fgLayers[fI].GetR();
963   for (Int_t i=0; i<fI; i++) {
964     Double_t xi=fgLayers[i].GetR();
965     d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
966   }
967
968   if (fI>1) {
969     Double_t xi=9.;
970     d+=0.0097*xi*xi;
971   }
972
973   if (fI>3) {
974     Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
975     d+=0.0034*xi*xi;
976   }
977
978   return d/(xn*xn);
979 }
980
981 Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
982                                 const AliITStrackV2 *c, Bool_t extra) {
983   //--------------------------------------------------------------------
984   // This function refits the track "t" at the position "x" using
985   // the clusters from "c"
986   // If "extra"==kTRUE, 
987   //    the clusters from overlapped modules get attached to "t" 
988   //--------------------------------------------------------------------
989   Int_t index[AliITSgeomTGeo::kNLayers];
990   Int_t k;
991   for (k=0; k<AliITSgeomTGeo::GetNLayers(); k++) index[k]=-1;
992   Int_t nc=c->GetNumberOfClusters();
993   for (k=0; k<nc; k++) { 
994     Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
995     index[nl]=idx; 
996   }
997
998   Int_t from, to, step;
999   if (xx > t->GetX()) {
1000       from=0; to=AliITSgeomTGeo::GetNLayers();
1001       step=+1;
1002   } else {
1003       from=AliITSgeomTGeo::GetNLayers()-1; to=-1;
1004       step=-1;
1005   }
1006
1007   for (Int_t i=from; i != to; i += step) {
1008      AliITSlayer &layer=fgLayers[i];
1009      Double_t r=layer.GetR();
1010  
1011      {
1012      Double_t hI=i-0.5*step; 
1013      if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {  
1014        Int_t iLay = i-step;
1015        Double_t rs = 0.;
1016        if(iLay<0 || iLay>= AliITSgeomTGeo::kNLayers){
1017          AliError(Form("Invalid layer %d ",iLay));
1018          return kFALSE;
1019        }
1020        else{
1021          rs=0.5*(fgLayers[i-step].GetR() + r);
1022        }
1023        Double_t d=0.0034, x0=38.6; 
1024        if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1025        if (!t->PropagateTo(rs,-step*d,x0)) {
1026          return kFALSE;
1027        }
1028      }
1029      }
1030
1031      // remember old position [SR, GSI 18.02.2003]
1032      Double_t oldX=0., oldY=0., oldZ=0.;
1033      if (t->IsStartedTimeIntegral() && step==1) {
1034         t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1035      }
1036      //
1037
1038      Double_t phi,z;
1039      if (!t->GetPhiZat(r,phi,z)) { 
1040        return kFALSE;
1041      }
1042
1043      Int_t idet=layer.FindDetectorIndex(phi,z);
1044      if (idet<0) { 
1045        return kFALSE;
1046      }
1047      const AliITSdetector &det=layer.GetDetector(idet);
1048      phi=det.GetPhi();
1049      if (!t->Propagate(phi,det.GetR())) {
1050        return kFALSE;
1051      }
1052      t->SetDetectorIndex(idet);
1053
1054      const AliITSRecPoint *cl=0;
1055      Double_t maxchi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
1056
1057      Int_t idx=index[i];
1058      if (idx>=0) {
1059         const AliITSRecPoint *ccc=(AliITSRecPoint *)GetCluster(idx); 
1060         if (idet != ccc->GetDetectorIndex()) {
1061            idet=ccc->GetDetectorIndex();
1062            const AliITSdetector &det2=layer.GetDetector(idet);
1063            if (!t->Propagate(det2.GetPhi(),det2.GetR())) {
1064              return kFALSE;
1065            }
1066            t->SetDetectorIndex(idet);
1067         }
1068         Double_t chi2=t->GetPredictedChi2(ccc);
1069         if (chi2<maxchi2) { 
1070           cl=ccc; 
1071           maxchi2=chi2; 
1072         } else {
1073           return kFALSE;
1074         }
1075      }
1076  
1077      if (cl) {
1078        // Take into account the mis-alignment
1079        Double_t x=t->GetX()+cl->GetX();
1080        if (!t->PropagateTo(x,0.,0.)) return kFALSE;
1081        if (!t->Update(cl,maxchi2,idx)) {
1082           return kFALSE;
1083        }
1084        t->SetSampledEdx(cl->GetQ(),i-2);
1085      }
1086
1087      {
1088      Double_t x0;
1089      Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1090      t->CorrectForMaterial(-step*d,x0);
1091      }
1092                  
1093      if (extra) { //search for extra clusters
1094         AliITStrackV2 tmp(*t);
1095         Double_t dz=4*TMath::Sqrt(tmp.GetSigmaZ2()+AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
1096         if (dz < 0.5*TMath::Abs(tmp.GetTgl())) dz=0.5*TMath::Abs(tmp.GetTgl());
1097         Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
1098         if (dy < 0.5*TMath::Abs(tmp.GetSnp())) dy=0.5*TMath::Abs(tmp.GetSnp());
1099         Double_t zmin=t->GetZ() - dz;
1100         Double_t zmax=t->GetZ() + dz;
1101         Double_t ymin=t->GetY() + phi*r - dy;
1102         Double_t ymax=t->GetY() + phi*r + dy;
1103         layer.SelectClusters(zmin,zmax,ymin,ymax);
1104
1105         const AliITSRecPoint *cx=0; Int_t ci=-1,cci=-1;
1106         maxchi2=1000.*AliITSReconstructor::GetRecoParam()->GetMaxChi2();
1107         Double_t tolerance=0.1;
1108         while ((cx=layer.GetNextCluster(ci))!=0) {
1109            if (idet == cx->GetDetectorIndex()) continue;
1110
1111            const AliITSdetector &detx=layer.GetDetector(cx->GetDetectorIndex());
1112
1113            if (!tmp.Propagate(detx.GetPhi(),detx.GetR())) continue;
1114            
1115            if (TMath::Abs(tmp.GetZ() - cx->GetZ()) > tolerance) continue;
1116            if (TMath::Abs(tmp.GetY() - cx->GetY()) > tolerance) continue;
1117
1118            Double_t chi2=tmp.GetPredictedChi2(cx);
1119            if (chi2<maxchi2) { maxchi2=chi2; cci=ci; }
1120         }
1121         if (cci>=0) t->SetExtraCluster(i,(i<<28)+cci);
1122      }
1123
1124      // track time update [SR, GSI 17.02.2003]
1125      if (t->IsStartedTimeIntegral() && step==1) {
1126         Double_t newX, newY, newZ;
1127         t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1128         Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) + 
1129                        (oldZ-newZ)*(oldZ-newZ);
1130         t->AddTimeStep(TMath::Sqrt(dL2));
1131      }
1132      //
1133
1134   }
1135
1136   if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1137   return kTRUE;
1138 }
1139
1140 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1141   //--------------------------------------------------------------------
1142   // This function marks clusters assigned to the track
1143   //--------------------------------------------------------------------
1144   AliTracker::UseClusters(t,from);
1145
1146   Int_t clusterIndex = t->GetClusterIndex(0);
1147   AliITSRecPoint *c= 0x0;
1148
1149   if (clusterIndex>-1)
1150     c = (AliITSRecPoint *)GetCluster(clusterIndex);
1151   if (c && c->GetSigmaZ2()>0.1) c->UnUse();
1152
1153   c = 0x0;
1154   clusterIndex = t->GetClusterIndex(1);
1155   if (clusterIndex>-1)
1156     c=(AliITSRecPoint *)GetCluster(clusterIndex);
1157   if (c && c->GetSigmaZ2()>0.1) c->UnUse();
1158
1159 }