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