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