]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
Bug fix. This eliminates a compilation warning. (R. Shahoyan)
[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=0;
286       try {
287         t=new AliITStrackV2(*esd);
288       } catch (const Char_t *msg) {
289         Warning("Clusters2Tracks",msg);
290         delete t;
291         continue;
292       }
293       if (TMath::Abs(t->GetD(GetX(),GetY()))>4) {
294         delete t;
295         continue;
296       }
297
298       if (CorrectForDeadZoneMaterial(t)!=0) {
299          Warning("Clusters2Tracks",
300                  "failed to correct for the material in the dead zone !\n");
301          delete t;
302          continue;
303       }
304       itsTracks.AddLast(t);
305     }
306   } /* End Read ESD tracks */
307
308   itsTracks.Sort();
309   Int_t nentr=itsTracks.GetEntriesFast();
310
311   Int_t ntrk=0;
312   for (fPass=0; fPass<2; fPass++) {
313      Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
314      for (Int_t i=0; i<nentr; i++) {
315        AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
316        if (t==0) continue;           //this track has been already tracked
317        Int_t tpcLabel=t->GetLabel(); //save the TPC track label
318
319        ResetTrackToFollow(*t);
320        ResetBestTrack();
321
322        for (FollowProlongation(); fI<AliITSgeomTGeo::GetNLayers(); fI++) {
323           while (TakeNextProlongation()) FollowProlongation();
324        }
325
326        if (fBestTrack.GetNumberOfClusters() == 0) continue;
327
328        if (fConstraint[fPass]) {
329           ResetTrackToFollow(*t);
330           if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
331           ResetBestTrack();
332        }
333
334        fBestTrack.SetLabel(tpcLabel);
335        fBestTrack.CookdEdx();
336        CookLabel(&fBestTrack,0.); //For comparison only
337        fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
338        UseClusters(&fBestTrack);
339        delete itsTracks.RemoveAt(i);
340        ntrk++;
341      }
342   }
343
344   itsTracks.Delete();
345
346   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
347
348   return 0;
349 }
350
351 Int_t AliITStrackerV2::PropagateBack(AliESDEvent *event) {
352   //--------------------------------------------------------------------
353   // This functions propagates reconstructed ITS tracks back
354   // The clusters must be loaded !
355   //--------------------------------------------------------------------
356   Int_t nentr=event->GetNumberOfTracks();
357   Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
358
359   Int_t ntrk=0;
360   for (Int_t i=0; i<nentr; i++) {
361      AliESDtrack *esd=event->GetTrack(i);
362
363      if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
364      if (esd->GetStatus()&AliESDtrack::kITSout) continue;
365
366      AliITStrackV2 *t=0;
367      try {
368         t=new AliITStrackV2(*esd);
369      } catch (const Char_t *msg) {
370         Warning("PropagateBack",msg);
371         delete t;
372         continue;
373      }
374
375      ResetTrackToFollow(*t);
376
377      // propagete to vertex [SR, GSI 17.02.2003]
378      // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
379      if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
380        if (fTrackToFollow.PropagateToVertex(event->GetVertex())) {
381           fTrackToFollow.StartTimeIntegral();
382        }
383        fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
384      }
385
386      fTrackToFollow.ResetCovariance(10.); fTrackToFollow.ResetClusters();
387      if (RefitAt(49.,&fTrackToFollow,t)) {
388         if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
389           Warning("PropagateBack",
390                   "failed to correct for the material in the dead zone !\n");
391           delete t;
392           continue;
393         }
394         fTrackToFollow.SetLabel(t->GetLabel());
395         //fTrackToFollow.CookdEdx();
396         CookLabel(&fTrackToFollow,0.); //For comparison only
397         fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
398         UseClusters(&fTrackToFollow);
399         ntrk++;
400      }
401      delete t;
402   }
403
404   Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
405
406   return 0;
407 }
408
409 Int_t AliITStrackerV2::RefitInward(AliESDEvent *event) {
410   //--------------------------------------------------------------------
411   // This functions refits ITS tracks using the 
412   // "inward propagated" TPC tracks
413   // The clusters must be loaded !
414   //--------------------------------------------------------------------
415   Int_t nentr=event->GetNumberOfTracks();
416   Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
417
418   Int_t ntrk=0;
419   for (Int_t i=0; i<nentr; i++) {
420     AliESDtrack *esd=event->GetTrack(i);
421
422     if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
423     if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
424     if (esd->GetStatus()&AliESDtrack::kTPCout)
425     if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
426
427     AliITStrackV2 *t=0;
428     try {
429         t=new AliITStrackV2(*esd);
430     } catch (const Char_t *msg) {
431         Warning("RefitInward",msg);
432         delete t;
433         continue;
434     }
435
436     if (CorrectForDeadZoneMaterial(t)!=0) {
437        Warning("RefitInward",
438                "failed to correct for the material in the dead zone !\n");
439        delete t;
440        continue;
441     }
442
443     ResetTrackToFollow(*t);
444     fTrackToFollow.ResetClusters();
445
446     //Refitting...
447     if (RefitAt(3.7, &fTrackToFollow, t, kTRUE)) {
448        fTrackToFollow.SetLabel(t->GetLabel());
449        fTrackToFollow.CookdEdx();
450        CookLabel(&fTrackToFollow,0.); //For comparison only
451
452        if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe 
453          AliESDtrack  *esdTrack =fTrackToFollow.GetESDtrack();
454          esdTrack->UpdateTrackParams(&fTrackToFollow,AliESDtrack::kITSrefit);
455          Double_t r[3]={0.,0.,0.};
456          Double_t maxD=3.;
457          esdTrack->RelateToVertex(event->GetVertex(),GetBz(r),maxD);
458          ntrk++;
459        }
460     }
461     delete t;
462   }
463
464   Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
465
466   return 0;
467 }
468
469 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
470   //--------------------------------------------------------------------
471   //       Return pointer to a given cluster
472   //--------------------------------------------------------------------
473   Int_t l=(index & 0xf0000000) >> 28;
474   Int_t c=(index & 0x0fffffff) >> 00;
475   return fgLayers[l].GetCluster(c);
476 }
477
478
479 void AliITStrackerV2::FollowProlongation() {
480   //--------------------------------------------------------------------
481   //This function finds a track prolongation 
482   //--------------------------------------------------------------------
483   while (fI>fLastLayerToTrackTo) {
484     Int_t i=fI-1;
485
486     AliITSlayer &layer=fgLayers[i];
487     AliITStrackV2 &track=fTracks[i];
488
489     Double_t r=layer.GetR();
490
491     if (i==3 || i==1) {
492        Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
493        Double_t d=0.0034, x0=38.6;
494        if (i==1) {rs=9.; d=0.0097; x0=42;}
495        if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
496          //Warning("FollowProlongation","propagation failed !\n");
497          return;
498        }
499     }
500
501     //find intersection
502     Double_t phi,z;  
503     if (!fTrackToFollow.GetPhiZat(r,phi,z)) {
504       //Warning("FollowProlongation","failed to estimate track !\n");
505       return;
506     }
507
508     Int_t idet=layer.FindDetectorIndex(phi,z);
509     if (idet<0) {
510       //Warning("FollowProlongation","failed to find a detector !\n");
511       return;
512     }
513
514     //propagate to the intersection
515     const AliITSdetector &det=layer.GetDetector(idet);
516     phi=det.GetPhi();
517     if (!fTrackToFollow.Propagate(phi,det.GetR())) {
518       //Warning("FollowProlongation","propagation failed !\n");
519       return;
520     }
521     fTrackToFollow.SetDetectorIndex(idet);
522
523     //Select possible prolongations and store the current track estimation
524     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
525     Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
526     Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
527     Double_t road=layer.GetRoad();
528     if (dz*dy>road*road) {
529        Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
530        dz=road*scz; dy=road*scy;
531     } 
532
533     //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
534     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
535     if (dz > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
536       //Warning("FollowProlongation","too broad road in Z !\n");
537       return;
538     }
539
540     if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
541
542     //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
543     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
544     if (dy > AliITSReconstructor::GetRecoParam()->GetMaxRoad()) {
545       //Warning("FollowProlongation","too broad road in Y !\n");
546       return;
547     }
548
549     fI--;
550
551     Double_t zmin=track.GetZ() - dz; 
552     Double_t zmax=track.GetZ() + dz;
553     Double_t ymin=track.GetY() + r*phi - dy;
554     Double_t ymax=track.GetY() + r*phi + dy;
555     if (layer.SelectClusters(zmin,zmax,ymin,ymax)==0) 
556        if (fLayersNotToSkip[fI]) return;  
557
558     if (!TakeNextProlongation()) 
559        if (fLayersNotToSkip[fI]) return;
560
561   } 
562
563   //deal with the best track
564   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
565   Int_t nclb=fBestTrack.GetNumberOfClusters();
566   if (ncl)
567   if (ncl >= nclb) {
568      Double_t chi2=fTrackToFollow.GetChi2();
569      if (chi2/ncl < AliITSReconstructor::GetRecoParam()->GetChi2PerCluster()) {        
570         if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
571            ResetBestTrack();
572         }
573      }
574   }
575
576 }
577
578 Int_t AliITStrackerV2::TakeNextProlongation() {
579   //--------------------------------------------------------------------
580   // This function takes another track prolongation 
581   //
582   //  dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
583   //--------------------------------------------------------------------
584   AliITSlayer &layer=fgLayers[fI];
585   ResetTrackToFollow(fTracks[fI]);
586
587   Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + AliITSReconstructor::GetRecoParam()->GetSigmaZ2(fI));
588   Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + AliITSReconstructor::GetRecoParam()->GetSigmaY2(fI));
589   Double_t road=layer.GetRoad();
590   if (dz*dy>road*road) {
591      Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
592      dz=road*scz; dy=road*scy;
593   } 
594
595   const AliITSRecPoint *c=0; Int_t ci=-1;
596   const AliITSRecPoint *cc=0; Int_t cci=-1;
597   Double_t chi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
598   while ((c=layer.GetNextCluster(ci))!=0) {
599     Int_t idet=c->GetDetectorIndex();
600
601     if (fTrackToFollow.GetDetectorIndex()!=idet) {
602        const AliITSdetector &det=layer.GetDetector(idet);
603        ResetTrackToFollow(fTracks[fI]);
604        if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
605          //Warning("TakeNextProlongation","propagation failed !\n");
606          continue;
607        }
608        fTrackToFollow.SetDetectorIndex(idet);
609        if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
610     }
611
612     if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
613     if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
614
615     Double_t ch2=fTrackToFollow.GetPredictedChi2(c); 
616     if (ch2 > chi2) continue;
617     chi2=ch2;
618     cc=c; cci=ci;
619     break;
620   }
621
622   if (!cc) return 0;
623
624   {// Take into account the mis-alignment
625     Double_t x = fTrackToFollow.GetX() + cc->GetX();
626     if (!fTrackToFollow.PropagateTo(x,0.,0.)) return 0;
627   }
628   if (!fTrackToFollow.Update(cc,chi2,(fI<<28)+cci)) {
629      //Warning("TakeNextProlongation","filtering failed !\n");
630      return 0;
631   }
632
633   if (fTrackToFollow.GetNumberOfClusters()>1)
634     if (TMath::Abs(fTrackToFollow.GetD(GetX(),GetY()))>4) return 0;
635
636   fTrackToFollow.
637     SetSampledEdx(cc->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
638
639   {
640   Double_t x0;
641  Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
642   fTrackToFollow.CorrectForMaterial(d,x0);
643   }
644
645   if (fConstraint[fPass]) {
646     Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
647     Double_t xyz[]={GetX(),GetY(),GetZ()};
648     Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
649     fTrackToFollow.Improve(d,xyz,ers);
650   }
651
652   return 1;
653 }
654
655
656 AliITStrackerV2::AliITSlayer::AliITSlayer():
657   fR(0.),
658   fPhiOffset(0.),
659   fNladders(0),
660   fZOffset(0.),
661   fNdetectors(0),
662   fDetectors(0),
663   fNsel(0),
664   fRoad(2*fR*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
665 {
666   //--------------------------------------------------------------------
667   //default AliITSlayer constructor
668   //--------------------------------------------------------------------
669   
670   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
671
672 }
673
674 AliITStrackerV2::AliITSlayer::
675 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd): 
676   fR(r), 
677   fPhiOffset(p), 
678   fNladders(nl),
679   fZOffset(z),
680   fNdetectors(nd),
681   fDetectors(new AliITSdetector[nl*nd]),
682   fNsel(0),
683   fRoad(2*r*TMath::Sqrt(3.14/1.)) //assuming that there's only one cluster
684 {
685   //--------------------------------------------------------------------
686   //main AliITSlayer constructor
687   //--------------------------------------------------------------------
688
689   for (Int_t i=0; i<kNsector; i++) fN[i]=0;
690
691   for (Int_t i=0; i<AliITSRecoParam::fgkMaxClusterPerLayer; i++) fClusters[i]=0;
692
693 }
694
695 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
696   //--------------------------------------------------------------------
697   // AliITSlayer destructor
698   //--------------------------------------------------------------------
699   delete[] fDetectors;
700   ResetClusters();
701 }
702
703 void AliITStrackerV2::AliITSlayer::ResetClusters() {
704   //--------------------------------------------------------------------
705   // This function removes loaded clusters
706   //--------------------------------------------------------------------
707    for (Int_t s=0; s<kNsector; s++) {
708        Int_t &n=fN[s];
709        while (n) {
710           n--;
711           delete fClusters[s*kMaxClusterPerSector+n];
712        }
713    }
714 }
715
716 void AliITStrackerV2::AliITSlayer::ResetRoad() {
717   //--------------------------------------------------------------------
718   // This function calculates the road defined by the cluster density
719   //--------------------------------------------------------------------
720   Int_t n=0;
721   for (Int_t s=0; s<kNsector; s++) {
722     Int_t i=fN[s];
723     while (i--) 
724        if (TMath::Abs(fClusters[s*kMaxClusterPerSector+i]->GetZ())<fR) n++;
725   }
726   if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
727 }
728
729 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSRecPoint *c) {
730   //--------------------------------------------------------------------
731   // This function inserts a cluster to this layer in increasing
732   // order of the cluster's fZ
733   //--------------------------------------------------------------------
734   Float_t circ=TMath::TwoPi()*fR;
735   Int_t sec=Int_t(kNsector*c->GetPhiR()/circ);
736   if (sec>=kNsector) {
737      ::Error("InsertCluster","Wrong sector !\n");
738      return 1;
739   }
740   Int_t &n=fN[sec];
741   if (n>=kMaxClusterPerSector) {
742      ::Error("InsertCluster","Too many clusters !\n");
743      return 1;
744   }
745   if (n==0) fClusters[sec*kMaxClusterPerSector]=c;
746   else {
747      Int_t i=FindClusterIndex(c->GetZ(),sec);
748      Int_t k=n-i+sec*kMaxClusterPerSector;
749      memmove(fClusters+i+1 ,fClusters+i,k*sizeof(AliITSRecPoint*));
750      fClusters[i]=c;
751   }
752   n++;
753   return 0;
754 }
755
756 Int_t 
757 AliITStrackerV2::AliITSlayer::FindClusterIndex(Float_t z,Int_t s) const {
758   //--------------------------------------------------------------------
759   // For the sector "s", this function returns the index of the first 
760   // with its fZ >= "z". 
761   //--------------------------------------------------------------------
762   Int_t nc=fN[s];
763   if (nc==0) return kMaxClusterPerSector*s;
764
765   Int_t b=kMaxClusterPerSector*s;
766   if (z <= fClusters[b]->GetZ()) return b;
767
768   Int_t e=b+nc-1;
769   if (z > fClusters[e]->GetZ()) return e+1;
770
771   Int_t m=(b+e)/2;
772   for (; b<e; m=(b+e)/2) {
773     if (z > fClusters[m]->GetZ()) b=m+1;
774     else e=m; 
775   }
776   return m;
777 }
778
779 Int_t AliITStrackerV2::AliITSlayer::
780 SelectClusters(Float_t zmin,Float_t zmax,Float_t ymin, Float_t ymax) {
781   //--------------------------------------------------------------------
782   // This function selects clusters within the "window"
783   //--------------------------------------------------------------------
784     Float_t circ=fR*TMath::TwoPi();
785
786     if (ymin>circ) ymin-=circ; else if (ymin<0) ymin+=circ;
787     if (ymax>circ) ymax-=circ; else if (ymax<0) ymax+=circ;
788
789     Int_t i1=Int_t(kNsector*ymin/circ); if (i1==kNsector) i1--;
790     if (fN[i1]!=0) {
791        Float_t ym = (ymax<ymin) ? ymax+circ : ymax;
792        Int_t i=FindClusterIndex(zmin,i1), imax=i1*kMaxClusterPerSector+fN[i1];
793        for (; i<imax; i++) {
794            AliITSRecPoint *c=fClusters[i];
795            if (c->IsUsed()) continue;
796            if (c->GetZ()>zmax) break;
797            if (c->GetPhiR()<=ymin) continue;
798            if (c->GetPhiR()>ym) continue;
799            fIndex[fNsel++]=i;
800        }
801     }
802
803     Int_t i2=Int_t(kNsector*ymax/circ); if (i2==kNsector) i2--;
804     if (i2==i1) return fNsel;
805
806     if (fN[i2]!=0) {
807        Float_t ym = (ymin>ymax) ? ymin-circ : ymin;
808        Int_t i=FindClusterIndex(zmin,i2), imax=i2*kMaxClusterPerSector+fN[i2];
809        for (; i<imax; i++) {
810            AliITSRecPoint *c=fClusters[i];
811            if (c->IsUsed()) continue;
812            if (c->GetZ()>zmax) break;
813            if (c->GetPhiR()<=ym) continue;
814            if (c->GetPhiR()>ymax) continue;
815            fIndex[fNsel++]=i;
816        }
817     }
818
819     return fNsel;
820 }
821
822 const AliITSRecPoint *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
823   //--------------------------------------------------------------------
824   // This function returns clusters within the "window" 
825   //--------------------------------------------------------------------
826   AliITSRecPoint *c=0;
827   ci=-1;
828   if (fNsel) {
829      fNsel--;
830      ci=fIndex[fNsel]; 
831      c=fClusters[ci];
832   }
833   return c; 
834 }
835
836 Int_t AliITStrackerV2::AliITSlayer::GetNumberOfClusters() const {
837   Int_t n=0;
838   for (Int_t s=0; s<kNsector; s++) n+=fN[s];
839   return n; 
840 }
841
842 Int_t 
843 AliITStrackerV2::AliITSlayer::FindDetectorIndex(Double_t phi,Double_t z)const {
844   //--------------------------------------------------------------------
845   //This function finds the detector crossed by the track
846   //--------------------------------------------------------------------
847   Double_t dphi;
848   if (fZOffset<0)            // old geometry
849     dphi = -(phi-fPhiOffset);
850   else                       // new geometry
851     dphi = phi-fPhiOffset;
852
853   if      (dphi <  0) dphi += 2*TMath::Pi();
854   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
855   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
856   if (np>=fNladders) np-=fNladders;
857   if (np<0)          np+=fNladders;
858
859   Double_t dz=fZOffset-z;
860   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
861   if (nz>=fNdetectors) return -1;
862   if (nz<0)            return -1;
863
864   return np*fNdetectors + nz;
865 }
866
867 Double_t 
868 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
869 const {
870   //--------------------------------------------------------------------
871   //This function returns the layer thickness at this point (units X0)
872   //--------------------------------------------------------------------
873   Double_t d=0.0085;
874   x0=21.82;
875
876   if (43<fR&&fR<45) { //SSD2
877      Double_t dd=0.0034;
878      d=dd;
879      if (TMath::Abs(y-0.00)>3.40) d+=dd;
880      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
881      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
882      for (Int_t i=0; i<12; i++) {
883        if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
884           if (TMath::Abs(y-0.00)>3.40) d+=dd;
885           d+=0.0034; 
886           break;
887        }
888        if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
889           if (TMath::Abs(y-0.00)>3.40) d+=dd;
890           d+=0.0034; 
891           break;
892        }         
893        if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
894        if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
895      }
896   } else 
897   if (37<fR&&fR<41) { //SSD1
898      Double_t dd=0.0034;
899      d=dd;
900      if (TMath::Abs(y-0.00)>3.40) d+=dd;
901      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
902      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
903      for (Int_t i=0; i<11; i++) {
904        if (TMath::Abs(z-3.9*i)<0.15) {
905           if (TMath::Abs(y-0.00)>3.40) d+=dd;
906           d+=dd; 
907           break;
908        }
909        if (TMath::Abs(z+3.9*i)<0.15) {
910           if (TMath::Abs(y-0.00)>3.40) d+=dd;
911           d+=dd; 
912           break;
913        }         
914        if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
915        if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}         
916      }
917   } else
918   if (13<fR&&fR<26) { //SDD
919      Double_t dd=0.0033;
920      d=dd;
921      if (TMath::Abs(y-0.00)>3.30) d+=dd;
922
923      if (TMath::Abs(y-1.80)<0.55) {
924         d+=0.016;
925         for (Int_t j=0; j<20; j++) {
926           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
927           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
928         } 
929      }
930      if (TMath::Abs(y+1.80)<0.55) {
931         d+=0.016;
932         for (Int_t j=0; j<20; j++) {
933           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
934           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
935         } 
936      }
937
938      for (Int_t i=0; i<4; i++) {
939        if (TMath::Abs(z-7.3*i)<0.60) {
940           d+=dd;
941           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
942           break;
943        }
944        if (TMath::Abs(z+7.3*i)<0.60) {
945           d+=dd; 
946           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
947           break;
948        }
949      }
950   } else
951   if (6<fR&&fR<8) {   //SPD2
952      Double_t dd=0.0063; x0=21.5;
953      d=dd;
954      if (TMath::Abs(y-3.08)>0.5) d+=dd;
955      //if (TMath::Abs(y-3.08)>0.45) d+=dd;
956      if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
957   } else
958   if (3<fR&&fR<5) {   //SPD1
959      Double_t dd=0.0063; x0=21.5;
960      d=dd;
961      if (TMath::Abs(y+0.21)>0.6) d+=dd;
962      //if (TMath::Abs(y+0.21)>0.45) d+=dd;
963      if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
964   }
965
966   return d;
967 }
968
969 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
970 {
971   //--------------------------------------------------------------------
972   //Returns the thickness between the current layer and the vertex (units X0)
973   //--------------------------------------------------------------------
974   Double_t d=0.0028*3*3; //beam pipe
975   Double_t x0=0;
976
977   Double_t xn=fgLayers[fI].GetR();
978   for (Int_t i=0; i<fI; i++) {
979     Double_t xi=fgLayers[i].GetR();
980     d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
981   }
982
983   if (fI>1) {
984     Double_t xi=9.;
985     d+=0.0097*xi*xi;
986   }
987
988   if (fI>3) {
989     Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
990     d+=0.0034*xi*xi;
991   }
992
993   return d/(xn*xn);
994 }
995
996 Bool_t AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,
997                                 const AliITStrackV2 *c, Bool_t extra) {
998   //--------------------------------------------------------------------
999   // This function refits the track "t" at the position "x" using
1000   // the clusters from "c"
1001   // If "extra"==kTRUE, 
1002   //    the clusters from overlapped modules get attached to "t" 
1003   //--------------------------------------------------------------------
1004   Int_t index[AliITSgeomTGeo::kNLayers];
1005   Int_t k;
1006   for (k=0; k<AliITSgeomTGeo::GetNLayers(); k++) index[k]=-1;
1007   Int_t nc=c->GetNumberOfClusters();
1008   for (k=0; k<nc; k++) { 
1009     Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
1010     index[nl]=idx; 
1011   }
1012
1013   Int_t from, to, step;
1014   if (xx > t->GetX()) {
1015       from=0; to=AliITSgeomTGeo::GetNLayers();
1016       step=+1;
1017   } else {
1018       from=AliITSgeomTGeo::GetNLayers()-1; to=-1;
1019       step=-1;
1020   }
1021
1022   for (Int_t i=from; i != to; i += step) {
1023      AliITSlayer &layer=fgLayers[i];
1024      Double_t r=layer.GetR();
1025  
1026      {
1027      Double_t hI=i-0.5*step; 
1028      if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {             
1029         Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
1030         Double_t d=0.0034, x0=38.6; 
1031         if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1032         if (!t->PropagateTo(rs,-step*d,x0)) {
1033           return kFALSE;
1034         }
1035      }
1036      }
1037
1038      // remember old position [SR, GSI 18.02.2003]
1039      Double_t oldX=0., oldY=0., oldZ=0.;
1040      if (t->IsStartedTimeIntegral() && step==1) {
1041         t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1042      }
1043      //
1044
1045      Double_t phi,z;
1046      if (!t->GetPhiZat(r,phi,z)) { 
1047        return kFALSE;
1048      }
1049
1050      Int_t idet=layer.FindDetectorIndex(phi,z);
1051      if (idet<0) { 
1052        return kFALSE;
1053      }
1054      const AliITSdetector &det=layer.GetDetector(idet);
1055      phi=det.GetPhi();
1056      if (!t->Propagate(phi,det.GetR())) {
1057        return kFALSE;
1058      }
1059      t->SetDetectorIndex(idet);
1060
1061      const AliITSRecPoint *cl=0;
1062      Double_t maxchi2=AliITSReconstructor::GetRecoParam()->GetMaxChi2();
1063
1064      Int_t idx=index[i];
1065      if (idx>=0) {
1066         const AliITSRecPoint *ccc=(AliITSRecPoint *)GetCluster(idx); 
1067         if (idet != ccc->GetDetectorIndex()) {
1068            idet=ccc->GetDetectorIndex();
1069            const AliITSdetector &det2=layer.GetDetector(idet);
1070            if (!t->Propagate(det2.GetPhi(),det2.GetR())) {
1071              return kFALSE;
1072            }
1073            t->SetDetectorIndex(idet);
1074         }
1075         Double_t chi2=t->GetPredictedChi2(ccc);
1076         if (chi2<maxchi2) { 
1077           cl=ccc; 
1078           maxchi2=chi2; 
1079         } else {
1080           return kFALSE;
1081         }
1082      }
1083  
1084      if (cl) {
1085        // Take into account the mis-alignment
1086        Double_t x=t->GetX()+cl->GetX();
1087        if (!t->PropagateTo(x,0.,0.)) return kFALSE;
1088        if (!t->Update(cl,maxchi2,idx)) {
1089           return kFALSE;
1090        }
1091        t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1092      }
1093
1094      {
1095      Double_t x0;
1096      Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1097      t->CorrectForMaterial(-step*d,x0);
1098      }
1099                  
1100      if (extra) { //search for extra clusters
1101         AliITStrackV2 tmp(*t);
1102         Double_t dz=4*TMath::Sqrt(tmp.GetSigmaZ2()+AliITSReconstructor::GetRecoParam()->GetSigmaZ2(i));
1103         if (dz < 0.5*TMath::Abs(tmp.GetTgl())) dz=0.5*TMath::Abs(tmp.GetTgl());
1104         Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+AliITSReconstructor::GetRecoParam()->GetSigmaY2(i));
1105         if (dy < 0.5*TMath::Abs(tmp.GetSnp())) dy=0.5*TMath::Abs(tmp.GetSnp());
1106         Double_t zmin=t->GetZ() - dz;
1107         Double_t zmax=t->GetZ() + dz;
1108         Double_t ymin=t->GetY() + phi*r - dy;
1109         Double_t ymax=t->GetY() + phi*r + dy;
1110         layer.SelectClusters(zmin,zmax,ymin,ymax);
1111
1112         const AliITSRecPoint *cx=0; Int_t ci=-1,cci=-1;
1113         maxchi2=1000.*AliITSReconstructor::GetRecoParam()->GetMaxChi2();
1114         Double_t tolerance=0.1;
1115         while ((cx=layer.GetNextCluster(ci))!=0) {
1116            if (idet == cx->GetDetectorIndex()) continue;
1117
1118            const AliITSdetector &detx=layer.GetDetector(cx->GetDetectorIndex());
1119
1120            if (!tmp.Propagate(detx.GetPhi(),detx.GetR())) continue;
1121            
1122            if (TMath::Abs(tmp.GetZ() - cx->GetZ()) > tolerance) continue;
1123            if (TMath::Abs(tmp.GetY() - cx->GetY()) > tolerance) continue;
1124
1125            Double_t chi2=tmp.GetPredictedChi2(cx);
1126            if (chi2<maxchi2) { maxchi2=chi2; cci=ci; }
1127         }
1128         if (cci>=0) t->SetExtraCluster(i,(i<<28)+cci);
1129      }
1130
1131      // track time update [SR, GSI 17.02.2003]
1132      if (t->IsStartedTimeIntegral() && step==1) {
1133         Double_t newX, newY, newZ;
1134         t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1135         Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) + 
1136                        (oldZ-newZ)*(oldZ-newZ);
1137         t->AddTimeStep(TMath::Sqrt(dL2));
1138      }
1139      //
1140
1141   }
1142
1143   if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1144   return kTRUE;
1145 }
1146
1147 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1148   //--------------------------------------------------------------------
1149   // This function marks clusters assigned to the track
1150   //--------------------------------------------------------------------
1151   AliTracker::UseClusters(t,from);
1152
1153   Int_t clusterIndex = t->GetClusterIndex(0);
1154   AliITSRecPoint *c= 0x0;
1155
1156   if (clusterIndex>-1)
1157     c = (AliITSRecPoint *)GetCluster(clusterIndex);
1158   if (c && c->GetSigmaZ2()>0.1) c->UnUse();
1159
1160   c = 0x0;
1161   clusterIndex = t->GetClusterIndex(1);
1162   if (clusterIndex>-1)
1163     c=(AliITSRecPoint *)GetCluster(clusterIndex);
1164   if (c && c->GetSigmaZ2()>0.1) c->UnUse();
1165
1166 }