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