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