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