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