]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerV2.cxx
Modified plots and made jet finder use SDigits
[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   Int_t detector=0;
121   for (Int_t i=0; i<kMaxLayer; i++) {
122     Int_t ndet=fgLayers[i].GetNdetectors();
123     Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
124     for (; j<jmax; j++) {           
125       if (!cTree->GetEvent(j)) continue;
126       Int_t ncl=clusters->GetEntriesFast();
127       while (ncl--) {
128         AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
129         detector = c->GetDetectorIndex();
130         fgLayers[i].InsertCluster(new AliITSclusterV2(*c));
131       }
132       clusters->Delete();
133       //add dead zone virtual "cluster"
134       
135       if (i<1){
136         for (Float_t ydead = -2.; ydead < 2. ; ydead+=0.05){     
137           Int_t lab[4] = {0,0,0,detector};
138           Int_t info[3] = {0,0,0};
139           Float_t hit[5]={ydead,0,1,0.01,0};
140           fgLayers[i].InsertCluster(new AliITSclusterV2(lab, hit, info));
141           //
142           hit[1]=-7.;
143           fgLayers[i].InsertCluster(new AliITSclusterV2(lab, hit, info));
144         }       
145       }
146       
147     }
148     //
149     fgLayers[i].ResetRoad(); //road defined by the cluster density
150   }
151
152   return 0;
153 }
154
155 void AliITStrackerV2::UnloadClusters() {
156   //--------------------------------------------------------------------
157   //This function unloads ITS clusters
158   //--------------------------------------------------------------------
159   for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
160 }
161
162 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
163   //--------------------------------------------------------------------
164   // Correction for the material between the TPC and the ITS
165   // (should it belong to the TPC code ?)
166   //--------------------------------------------------------------------
167   Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ? 
168   Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
169   Double_t yr=12.8, dr=0.03; // rods ?
170   Double_t zm=0.2, dm=0.40;  // membrane
171   //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
172   Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
173
174   if (t->GetX() > riw) {
175      if (!t->PropagateTo(riw,diw,x0iw)) return 1;
176      if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
177      if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
178      if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
179      //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
180      //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r); 
181      if (!t->PropagateTo(rs,ds)) return 1;
182   } else if (t->GetX() < rs) {
183      if (!t->PropagateTo(rs,-ds)) return 1;
184      //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
185      //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r); 
186      if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
187      if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
188   } else {
189   ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
190     return 1;
191   }
192   
193   return 0;
194 }
195
196 Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
197   //--------------------------------------------------------------------
198   // This functions reconstructs ITS tracks
199   // The clusters must be already loaded !
200   //--------------------------------------------------------------------
201   TObjArray itsTracks(15000);
202
203   {/* Read ESD tracks */
204     Int_t nentr=event->GetNumberOfTracks();
205     Info("Clusters2Tracks", "Number of ESD tracks: %d\n", nentr);
206     while (nentr--) {
207       AliESDtrack *esd=event->GetTrack(nentr);
208
209       if ((esd->GetStatus()&AliESDtrack::kTPCin)==0) continue;
210       if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
211       if (esd->GetStatus()&AliESDtrack::kITSin) continue;
212
213       AliITStrackV2 *t=0;
214       try {
215         t=new AliITStrackV2(*esd);
216       } catch (const Char_t *msg) {
217         Warning("Clusters2Tracks",msg);
218         delete t;
219         continue;
220       }
221       if (TMath::Abs(t->GetD())>5) {
222         delete t;
223         continue;
224       }
225
226       if (CorrectForDeadZoneMaterial(t)!=0) {
227          Warning("Clusters2Tracks",
228                  "failed to correct for the material in the dead zone !\n");
229          delete t;
230          continue;
231       }
232       t->fReconstructed = kFALSE;
233       itsTracks.AddLast(t);
234     }
235   } /* End Read ESD tracks */
236
237   itsTracks.Sort();
238   Int_t nentr=itsTracks.GetEntriesFast();
239   fTrackHypothesys.Expand(nentr);
240   Int_t ntrk=0;
241   for (fPass=0; fPass<2; fPass++) {
242      Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
243      for (Int_t i=0; i<nentr; i++) {
244        fCurrentEsdTrack = i;
245        AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
246        if (t==0) continue;              //this track has been already tracked
247        if (t->fReconstructed) continue;  //this track was  already  "succesfully" reconstructed
248        if ( (TMath::Abs(t->GetD(GetX(),GetY()))  >2.) && fConstraint[fPass]) continue;
249        if ( (TMath::Abs(t->GetZat(GetX())-GetZ())>2.) && fConstraint[fPass]) continue;
250
251        Int_t tpcLabel=t->GetLabel(); //save the TPC track label
252
253        ResetTrackToFollow(*t);
254        ResetBestTrack();
255
256        for (FollowProlongation(); fI<kMaxLayer; fI++) {
257          while (TakeNextProlongation()) {
258             FollowProlongation();
259          }
260        }
261
262        SortTrackHypothesys(fCurrentEsdTrack,0.9);  //MI change
263        CompressTrackHypothesys(fCurrentEsdTrack,0.90,2);  //MI change
264
265        /*
266        if (fBestTrack.GetNumberOfClusters() == 0) continue;
267        
268        if (fConstraint[fPass]) {
269          ResetTrackToFollow(*t);
270          if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
271          ResetBestTrack();
272          }
273        
274        fBestTrack.SetLabel(tpcLabel);
275        fBestTrack.CookdEdx();
276        CookLabel(&fBestTrack,0.); //For comparison only
277        fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
278        //UseClusters(&fBestTrack);
279        delete itsTracks.RemoveAt(i);
280        ntrk++;
281        */
282        AliITStrackV2 * besttrack = GetBestHypothesys(fCurrentEsdTrack,t,5);
283        if (!besttrack) continue;
284        besttrack->SetLabel(tpcLabel);
285        besttrack->CookdEdx();
286        besttrack->fFakeRatio=1.;
287        CookLabel(besttrack,0.); //For comparison only
288        //       besttrack->UpdateESDtrack(AliESDtrack::kITSin);
289        
290        if (besttrack->GetChi2()/besttrack->GetNumberOfClusters()>3.5){
291          if ( (TMath::Abs(besttrack->GetD(GetX(),GetY()))>0.4) && fConstraint[fPass]) {
292            CompressTrackHypothesys(fCurrentEsdTrack,0.0,0);
293            continue;
294          }
295          if ( (TMath::Abs(besttrack->GetZat(GetX()) -GetZ() )>0.4) && fConstraint[fPass]){
296            CompressTrackHypothesys(fCurrentEsdTrack,0.0,0);
297            continue;
298          }
299        }
300        
301        //delete itsTracks.RemoveAt(i);
302        t->fReconstructed = kTRUE;
303        ntrk++;              
304        
305      }
306   }
307   
308   for (Int_t i=0; i<nentr; i++) {
309     fCurrentEsdTrack = i;
310     AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
311     if (t==0) continue;         
312     Int_t tpcLabel=t->GetLabel(); //save the TPC track label
313     AliITStrackV2 * besttrack = GetBestHypothesysMIP(fCurrentEsdTrack,t);
314     if (!besttrack) continue;
315
316     besttrack->SetLabel(tpcLabel);
317     besttrack->CookdEdx();
318     besttrack->fFakeRatio=1.;
319     CookLabel(besttrack,0.); //For comparison only
320     besttrack->UpdateESDtrack(AliESDtrack::kITSin);
321   }
322   //
323
324   itsTracks.Delete();
325   //
326   Int_t entries = fTrackHypothesys.GetEntriesFast();
327   for (Int_t ientry=0;ientry<entries;ientry++){
328     TObjArray * array =(TObjArray*)fTrackHypothesys.UncheckedAt(ientry);
329     if (array) array->Delete();
330     delete fTrackHypothesys.RemoveAt(ientry); 
331   }
332
333   fTrackHypothesys.Delete();
334   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
335
336   return 0;
337 }
338
339 Int_t AliITStrackerV2::Clusters2Tracks(TTree *tpcTree, TTree *itsTree) {
340   //--------------------------------------------------------------------
341   // This functions reconstructs ITS tracks
342   // The clusters must be already loaded !
343   //--------------------------------------------------------------------
344   Int_t nentr=0; TObjArray itsTracks(15000);
345
346    Warning("Clusters2Tracks(TTree *, TTree *)",
347       "Will be removed soon !   Use Clusters2Tracks(AliESD *) instead.");
348
349   {/* Read TPC tracks */ 
350     AliTPCtrack *itrack=new AliTPCtrack; 
351     TBranch *branch=tpcTree->GetBranch("tracks");
352     if (!branch) {
353        Error("Clusters2Tracks","Can't get the branch !");
354        return 1;
355     }
356     tpcTree->SetBranchAddress("tracks",&itrack);
357     nentr=(Int_t)tpcTree->GetEntries();
358
359     Info("Clusters2Tracks","Number of TPC tracks: %d\n",nentr);
360
361     for (Int_t i=0; i<nentr; i++) {
362        tpcTree->GetEvent(i);
363        AliITStrackV2 *t=0;
364        try {
365            t=new AliITStrackV2(*itrack);
366        } catch (const Char_t *msg) {
367            Warning("Clusters2Tracks",msg);
368            delete t;
369            continue;
370        }
371        if (TMath::Abs(t->GetD())>4) {
372          delete t;
373          continue;
374        }
375
376        if (CorrectForDeadZoneMaterial(t)!=0) {
377          Warning("Clusters2Tracks",
378                  "failed to correct for the material in the dead zone !\n");
379          delete t;
380          continue;
381        }
382
383        itsTracks.AddLast(t);
384     }
385     delete itrack;
386   }
387   itsTracks.Sort();
388   nentr=itsTracks.GetEntriesFast();
389
390
391   AliITStrackV2 *otrack=&fBestTrack;
392   TBranch *branch=itsTree->GetBranch("tracks");
393   if (!branch) itsTree->Branch("tracks","AliITStrackV2",&otrack,32000,3);
394   else branch->SetAddress(&otrack);
395
396   for (fPass=0; fPass<2; fPass++) {
397      Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
398      for (Int_t i=0; i<nentr; i++) {
399        AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
400        if (t==0) continue;           //this track has been already tracked
401        Int_t tpcLabel=t->GetLabel(); //save the TPC track label
402
403        ResetTrackToFollow(*t);
404        ResetBestTrack();
405
406        for (FollowProlongation(); fI<kMaxLayer; fI++) {
407           while (TakeNextProlongation()) FollowProlongation();
408        }
409
410        if (fBestTrack.GetNumberOfClusters() == 0) continue;
411
412        if (fConstraint[fPass]) {
413           ResetTrackToFollow(*t);
414           if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
415           ResetBestTrack();
416        }
417
418        fBestTrack.SetLabel(tpcLabel);
419        fBestTrack.CookdEdx();
420        CookLabel(&fBestTrack,0.); //For comparison only
421        itsTree->Fill();
422        //UseClusters(&fBestTrack);
423        delete itsTracks.RemoveAt(i);
424      }
425   }
426
427   nentr=(Int_t)itsTree->GetEntries();
428   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",nentr);
429
430   itsTracks.Delete();
431
432   return 0;
433 }
434
435 Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
436   //--------------------------------------------------------------------
437   // This functions propagates reconstructed ITS tracks back
438   // The clusters must be loaded !
439   //--------------------------------------------------------------------
440   Int_t nentr=event->GetNumberOfTracks();
441   Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
442
443   Int_t ntrk=0;
444   for (Int_t i=0; i<nentr; i++) {
445      AliESDtrack *esd=event->GetTrack(i);
446
447      if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
448      if (esd->GetStatus()&AliESDtrack::kITSout) continue;
449
450      AliITStrackV2 *t=0;
451      try {
452         t=new AliITStrackV2(*esd);
453      } catch (const Char_t *msg) {
454         Warning("PropagateBack",msg);
455         delete t;
456         continue;
457      }
458
459      ResetTrackToFollow(*t);
460
461      // propagete to vertex [SR, GSI 17.02.2003]
462      // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
463      if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
464        if (fTrackToFollow.PropagateToVertex()) {
465           fTrackToFollow.StartTimeIntegral();
466        }
467        fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
468      }
469
470      fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
471      if (RefitAt(49.,&fTrackToFollow,t)) {
472         if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
473           Warning("PropagateBack",
474                   "failed to correct for the material in the dead zone !\n");
475           delete t;
476           continue;
477         }
478         fTrackToFollow.SetLabel(t->GetLabel());
479         //fTrackToFollow.CookdEdx();
480         CookLabel(&fTrackToFollow,0.); //For comparison only
481         fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
482         //UseClusters(&fTrackToFollow);
483         ntrk++;
484      }
485      delete t;
486   }
487
488   Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
489
490   return 0;
491 }
492
493 Int_t AliITStrackerV2::RefitInward(AliESD *event) {
494   //--------------------------------------------------------------------
495   // This functions refits ITS tracks using the 
496   // "inward propagated" TPC tracks
497   // The clusters must be loaded !
498   //--------------------------------------------------------------------
499   Int_t nentr=event->GetNumberOfTracks();
500   Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
501
502   Int_t ntrk=0;
503   for (Int_t i=0; i<nentr; i++) {
504     AliESDtrack *esd=event->GetTrack(i);
505
506     if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
507     if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
508     if (esd->GetStatus()&AliESDtrack::kTPCout)
509       if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
510
511     AliITStrackV2 *t=0;
512     try {
513         t=new AliITStrackV2(*esd);
514     } catch (const Char_t *msg) {
515         Warning("RefitInward",msg);
516         delete t;
517         continue;
518     }
519
520     if (CorrectForDeadZoneMaterial(t)!=0) {
521        Warning("RefitInward",
522                "failed to correct for the material in the dead zone !\n");
523        delete t;
524        continue;
525     }
526
527     ResetTrackToFollow(*t);
528     fTrackToFollow.ResetClusters();
529
530     if ((esd->GetStatus()&AliESDtrack::kTPCin)==0)
531       fTrackToFollow.ResetCovariance();
532
533     //Refitting...
534     if (RefitAt(3.7, &fTrackToFollow, t)) {
535        fTrackToFollow.SetLabel(t->GetLabel());
536        fTrackToFollow.CookdEdx();
537        CookLabel(&fTrackToFollow,0.0); //For comparison only
538
539        if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe    
540          Double_t a=fTrackToFollow.GetAlpha();
541          Double_t cs=TMath::Cos(a),sn=TMath::Sin(a);
542          Double_t xv= GetX()*cs + GetY()*sn;
543          Double_t yv=-GetX()*sn + GetY()*cs;
544          
545          Double_t c=fTrackToFollow.GetC(), snp=fTrackToFollow.GetSnp();
546          Double_t x=fTrackToFollow.GetX(), y=fTrackToFollow.GetY();
547          Double_t tgfv=-(c*(x-xv)-snp)/(c*(y-yv) + TMath::Sqrt(1.-snp*snp));
548          Double_t fv=TMath::ATan(tgfv);
549
550          cs=TMath::Cos(fv); sn=TMath::Sin(fv);
551          x = xv*cs + yv*sn;
552          yv=-xv*sn + yv*cs; xv=x;
553
554          if (fTrackToFollow.Propagate(fv+a,xv)) {
555             fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
556             //UseClusters(&fTrackToFollow);
557             {
558             AliITSclusterV2 c; c.SetY(yv); c.SetZ(GetZ());
559             c.SetSigmaY2(GetSigmaY()*GetSigmaY());
560             c.SetSigmaZ2(GetSigmaZ()*GetSigmaZ());
561             Double_t chi2=fTrackToFollow.GetPredictedChi2(&c);
562             if (chi2<kMaxChi2)
563                if (fTrackToFollow.Update(&c,-chi2,0))
564                    fTrackToFollow.SetConstrainedESDtrack(chi2);            
565             }
566             ntrk++;
567          }
568        }
569     }
570     delete t;
571   }
572
573   Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
574
575   return 0;
576 }
577
578 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
579   //--------------------------------------------------------------------
580   //       Return pointer to a given cluster
581   //--------------------------------------------------------------------
582   Int_t l=(index & 0xf0000000) >> 28;
583   Int_t c=(index & 0x0fffffff) >> 00;
584   return fgLayers[l].GetCluster(c);
585 }
586
587
588 void AliITStrackerV2::FollowProlongation() {
589   //--------------------------------------------------------------------
590   //This function finds a track prolongation 
591   //--------------------------------------------------------------------
592   while (fI>fLastLayerToTrackTo) {
593     Int_t i=fI-1;
594
595     AliITSlayer &layer=fgLayers[i];
596     AliITStrackV2 &track=fTracks[i];
597
598     Double_t r=layer.GetR();
599
600     if (i==3 || i==1) {
601        Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
602        Double_t d=0.0034, x0=38.6;
603        if (i==1) {rs=9.; d=0.0097; x0=42;}
604        if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
605          //Warning("FollowProlongation","propagation failed !\n");
606          return;
607        }
608     }
609
610     //find intersection
611     Double_t x,y,z;  
612     if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
613       //Warning("FollowProlongation","failed to estimate track !\n");
614       return;
615     }
616     Double_t phi=TMath::ATan2(y,x);
617
618     Int_t idet=layer.FindDetectorIndex(phi,z);
619     if (idet<0) {
620       //Warning("FollowProlongation","failed to find a detector !\n");
621       return;
622     }
623
624     //propagate to the intersection
625     const AliITSdetector &det=layer.GetDetector(idet);
626     phi=det.GetPhi();
627     if (!fTrackToFollow.Propagate(phi,det.GetR())) {
628       //Warning("FollowProlongation","propagation failed !\n");
629       return;
630     }
631     fTrackToFollow.SetDetectorIndex(idet);
632
633     //Select possible prolongations and store the current track estimation
634     track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
635     Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
636     Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
637     Double_t road=layer.GetRoad();
638     if (dz*dy>road*road) {
639        Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
640        dz=road*scz; dy=road*scy;
641     } 
642
643     //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
644     if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
645     if (dz > kMaxRoad) {
646       //Warning("FollowProlongation","too broad road in Z !\n");
647       return;
648     }
649
650     //    if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
651
652     //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
653     if (TMath::Abs(track.GetSnp()>kMaxSnp)) {
654       fI--;
655       continue;   // MI
656     }
657     if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
658     if (dy > kMaxRoad) {
659       //Warning("FollowProlongation","too broad road in Y !\n");
660       return;
661     }
662
663     Double_t zmin=track.GetZ() - dz; 
664     Double_t zmax=track.GetZ() + dz;
665     Double_t ymin=track.GetY() + r*phi - dy;
666     Double_t ymax=track.GetY() + r*phi + dy;
667     layer.SelectClusters(zmin,zmax,ymin,ymax); 
668     fI--;
669
670     //take another prolongation
671     if (!TakeNextProlongation()){ 
672       //skipped++;
673       fTrackToFollow.fNSkipped++;
674       if (fLayersNotToSkip[fI]||fTrackToFollow.fNSkipped>1) return;
675     }
676     if (fTrackToFollow.fNUsed>1) return;
677     if (fTrackToFollow.fNUsed+fTrackToFollow.fNSkipped>1) return;    
678     if ( (fI<3) && ( fTrackToFollow.GetChi2()/fTrackToFollow.GetNumberOfClusters()>kChi2PerCluster)) return;
679   } 
680
681   //deal with the best track
682   Int_t ncl=fTrackToFollow.GetNumberOfClusters();
683   Int_t nclb=fBestTrack.GetNumberOfClusters();
684   if (ncl){
685     if (ncl<4) return;
686     if ( (ncl<6) && (fTrackToFollow.GetChi2()/float(ncl))>3) return;
687     if (fTrackToFollow.GetChi2()/ncl>5.5) return;
688     fTrackToFollow.CookdEdx();
689     if (fTrackToFollow.fESDtrack->GetTPCsignal()>80.)
690       if ((fTrackToFollow.GetdEdx()/fTrackToFollow.fESDtrack->GetTPCsignal())<0.35){
691         // mismatch in dedx
692         return;
693       }
694     //
695     fTrackToFollow.SetLabel(fTrackToFollow.fESDtrack->GetLabel());
696     CookLabel(&fTrackToFollow,0.); //
697     //
698     if ( (nclb>3) && ((fTrackToFollow.GetChi2()/ncl)<(3*fBestTrack.GetChi2()/(nclb))))
699       AddTrackHypothesys(new AliITStrackV2(fTrackToFollow), fCurrentEsdTrack);
700     else 
701       if (ncl>3) AddTrackHypothesys(new AliITStrackV2(fTrackToFollow), fCurrentEsdTrack);
702     if (ncl >= nclb) {
703       Double_t chi2=fTrackToFollow.GetChi2();
704       if (chi2/ncl < kChi2PerCluster) {        
705         if (ncl > nclb ) {
706           ResetBestTrack();
707         }
708         if ( (ncl == nclb) && chi2 < fBestTrack.GetChi2()) {
709           ResetBestTrack();
710         }       
711       }
712     }
713   }
714 }
715
716 Int_t AliITStrackerV2::TakeNextProlongation() {
717   //--------------------------------------------------------------------
718   // This function takes another track prolongation 
719   //
720   //  dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
721   //--------------------------------------------------------------------
722   AliITSlayer &layer=fgLayers[fI];
723   ResetTrackToFollow(fTracks[fI]);
724
725   Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
726   Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
727   Double_t road=layer.GetRoad();
728   if (dz*dy>road*road) {
729      Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
730      dz=road*scz; dy=road*scy;
731   } 
732
733   const AliITSclusterV2 *c=0; Int_t ci=-1;
734   Double_t chi2=12345.;
735   while ((c=layer.GetNextCluster(ci))!=0) {
736     Int_t idet=c->GetDetectorIndex();
737     if (fTrackToFollow.GetDetectorIndex()!=idet) {
738        const AliITSdetector &det=layer.GetDetector(idet);
739        ResetTrackToFollow(fTracks[fI]);
740        if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
741          //Warning("TakeNextProlongation","propagation failed !\n");
742          continue;
743        }
744        fTrackToFollow.SetDetectorIndex(idet);
745        if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
746     }
747
748     if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
749     if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
750
751     chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
752   }
753
754   if (chi2>=kMaxChi2) return 0;
755   if (!c) return 0;
756   if (c->IsUsed()&&c->GetNy()<5) {  //shared factor    
757     chi2+=1;
758     chi2*=2*(1./(TMath::Max(c->GetNy(),1)));
759   }
760   if (c->GetQ()==0){  //dead zone virtual cluster
761     chi2*=4.;
762     chi2+=20;
763     fTrackToFollow.fNUsed++;
764     return 1;
765   }
766   //if ((fI<2)&&chi2>kMaxChi2In) return 0;
767
768   if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
769      //Warning("TakeNextProlongation","filtering failed !\n");
770      return 0;
771   }
772   if (c->IsUsed()&&c->GetNy()<5) fTrackToFollow.fNUsed++;
773   if (fTrackToFollow.GetNumberOfClusters()>1)
774   if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
775
776   fTrackToFollow.
777     SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
778
779   {
780   Double_t x0;
781  Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
782   fTrackToFollow.CorrectForMaterial(d,x0);
783   }
784
785   if (fConstraint[fPass]) {
786     Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
787     Double_t xyz[]={GetX(),GetY(),GetZ()};
788     Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
789     Double_t deltad = TMath::Abs(fTrackToFollow.GetD(GetX(),GetY()));
790     Double_t deltaz = TMath::Abs(fTrackToFollow.GetZat(GetX())-GetZ());
791
792     if ( (fI==4) &&  (deltad>2.0 || deltaz>1.5))  return 0; // don't "improve" secondaries     
793     if ( (fI==3) &&  (deltad>1.5 || deltaz>0.9))  return 0; // don't "improve" secondaries 
794     if ( (fI==2) &&  (deltad>0.9 || deltaz>0.6))  return 1; // don't "improve" secondaries 
795     if ( (fI==1) &&  (deltad>0.3 || deltaz>0.3))  return 1; // don't "improve" secondaries 
796     if ( (fI==0) &&  (deltad>0.1 || deltaz>0.1))  return 1; // don't "improve" secondaries 
797
798     fTrackToFollow.Improve(d,xyz,ers);
799   }
800
801   return 1;
802 }
803
804
805 AliITStrackerV2::AliITSlayer::AliITSlayer() {
806   //--------------------------------------------------------------------
807   //default AliITSlayer constructor
808   //--------------------------------------------------------------------
809   fN=0;
810   fDetectors=0;
811 }
812
813 AliITStrackerV2::AliITSlayer::
814 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
815   //--------------------------------------------------------------------
816   //main AliITSlayer constructor
817   //--------------------------------------------------------------------
818   fR=r; fPhiOffset=p; fZOffset=z;
819   fNladders=nl; fNdetectors=nd;
820   fDetectors=new AliITSdetector[fNladders*fNdetectors];
821
822   fN=0;
823   fI=0;
824
825   fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
826 }
827
828 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
829   //--------------------------------------------------------------------
830   // AliITSlayer destructor
831   //--------------------------------------------------------------------
832   delete[] fDetectors;
833   for (Int_t i=0; i<fN; i++) delete fClusters[i];
834   for (Int_t i=0; i<kMaxClusterPerLayer;i++) fClusterWeight[i]=0;
835 }
836
837 void AliITStrackerV2::AliITSlayer::ResetClusters() {
838   //--------------------------------------------------------------------
839   // This function removes loaded clusters
840   //--------------------------------------------------------------------
841   for (Int_t i=0; i<fN; i++) delete fClusters[i];
842   for (Int_t i=0; i<kMaxClusterPerLayer;i++) fClusterWeight[i]=0;
843   fN=0;
844   fI=0;
845 }
846
847 void AliITStrackerV2::AliITSlayer::ResetWeights() {
848   //--------------------------------------------------------------------
849   // This function reset weights of the clusters
850   //--------------------------------------------------------------------
851   for (Int_t i=0; i<kMaxClusterPerLayer;i++) fClusterWeight[i]=0;
852 }
853
854 void AliITStrackerV2::AliITSlayer::ResetRoad() {
855   //--------------------------------------------------------------------
856   // This function calculates the road defined by the cluster density
857   //--------------------------------------------------------------------
858   Int_t n=0;
859   for (Int_t i=0; i<fN; i++) {
860      if (TMath::Abs(fClusters[i]->GetZ())<fR) n++;
861   }
862   if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
863 }
864
865 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
866   //--------------------------------------------------------------------
867   //This function adds a cluster to this layer
868   //--------------------------------------------------------------------
869   if (fN==kMaxClusterPerLayer) {
870     ::Error("InsertCluster","Too many clusters !\n");
871     return 1;
872   }
873
874   if (fN==0) {fClusters[fN++]=c; return 0;}
875   Int_t i=FindClusterIndex(c->GetZ());
876   memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
877   fClusters[i]=c; fN++;
878
879   return 0;
880 }
881
882 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
883   //--------------------------------------------------------------------
884   // This function returns the index of the nearest cluster 
885   //--------------------------------------------------------------------
886   if (fN==0) return 0;
887   if (z <= fClusters[0]->GetZ()) return 0;
888   if (z > fClusters[fN-1]->GetZ()) return fN;
889   Int_t b=0, e=fN-1, m=(b+e)/2;
890   for (; b<e; m=(b+e)/2) {
891     if (z > fClusters[m]->GetZ()) b=m+1;
892     else e=m; 
893   }
894   return m;
895 }
896
897 void AliITStrackerV2::AliITSlayer::
898 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
899   //--------------------------------------------------------------------
900   // This function sets the "window"
901   //--------------------------------------------------------------------
902   fI=FindClusterIndex(zmin); fZmax=zmax;
903   Double_t circle=2*TMath::Pi()*fR;
904   if (ymax>circle) { ymax-=circle; ymin-=circle; }
905   fYmin=ymin; fYmax=ymax;
906 }
907
908 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
909   //--------------------------------------------------------------------
910   // This function returns clusters within the "window" 
911   //--------------------------------------------------------------------
912   const AliITSclusterV2 *cluster=0;
913   for (Int_t i=fI; i<fN; i++) {
914     const AliITSclusterV2 *c=fClusters[i];
915     if (c->GetZ() > fZmax) break;
916     //    if (c->IsUsed()) continue;
917     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
918     Double_t y=fR*det.GetPhi() + c->GetY();
919
920     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
921     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
922
923     if (y<fYmin) continue;
924     if (y>fYmax) continue;
925     cluster=c; ci=i;
926     fI=i+1;
927     break; 
928   }
929
930   return cluster;
931 }
932
933 Int_t AliITStrackerV2::AliITSlayer::
934 FindDetectorIndex(Double_t phi, Double_t z) const {
935   //--------------------------------------------------------------------
936   //This function finds the detector crossed by the track
937   //--------------------------------------------------------------------
938   Double_t dphi=-(phi-fPhiOffset);
939   if      (dphi <  0) dphi += 2*TMath::Pi();
940   else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
941   Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
942   if (np>=fNladders) np-=fNladders;
943   if (np<0)          np+=fNladders;
944
945   Double_t dz=fZOffset-z;
946   Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
947   if (nz>=fNdetectors) return -1;
948   if (nz<0)            return -1;
949
950   return np*fNdetectors + nz;
951 }
952
953 Double_t 
954 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
955 const {
956   //--------------------------------------------------------------------
957   //This function returns the layer thickness at this point (units X0)
958   //--------------------------------------------------------------------
959   Double_t d=0.0085;
960   x0=21.82;
961
962   if (43<fR&&fR<45) { //SSD2
963      Double_t dd=0.0034;
964      d=dd;
965      if (TMath::Abs(y-0.00)>3.40) d+=dd;
966      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
967      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
968      for (Int_t i=0; i<12; i++) {
969        if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
970           if (TMath::Abs(y-0.00)>3.40) d+=dd;
971           d+=0.0034; 
972           break;
973        }
974        if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
975           if (TMath::Abs(y-0.00)>3.40) d+=dd;
976           d+=0.0034; 
977           break;
978        }         
979        if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
980        if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
981      }
982   } else 
983   if (37<fR&&fR<41) { //SSD1
984      Double_t dd=0.0034;
985      d=dd;
986      if (TMath::Abs(y-0.00)>3.40) d+=dd;
987      if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
988      if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
989      for (Int_t i=0; i<11; i++) {
990        if (TMath::Abs(z-3.9*i)<0.15) {
991           if (TMath::Abs(y-0.00)>3.40) d+=dd;
992           d+=dd; 
993           break;
994        }
995        if (TMath::Abs(z+3.9*i)<0.15) {
996           if (TMath::Abs(y-0.00)>3.40) d+=dd;
997           d+=dd; 
998           break;
999        }         
1000        if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
1001        if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}         
1002      }
1003   } else
1004   if (13<fR&&fR<26) { //SDD
1005      Double_t dd=0.0033;
1006      d=dd;
1007      if (TMath::Abs(y-0.00)>3.30) d+=dd;
1008
1009      if (TMath::Abs(y-1.80)<0.55) {
1010         d+=0.016;
1011         for (Int_t j=0; j<20; j++) {
1012           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1013           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1014         } 
1015      }
1016      if (TMath::Abs(y+1.80)<0.55) {
1017         d+=0.016;
1018         for (Int_t j=0; j<20; j++) {
1019           if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1020           if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1021         } 
1022      }
1023
1024      for (Int_t i=0; i<4; i++) {
1025        if (TMath::Abs(z-7.3*i)<0.60) {
1026           d+=dd;
1027           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
1028           break;
1029        }
1030        if (TMath::Abs(z+7.3*i)<0.60) {
1031           d+=dd; 
1032           if (TMath::Abs(y-0.00)>3.30) d+=dd; 
1033           break;
1034        }
1035      }
1036   } else
1037   if (6<fR&&fR<8) {   //SPD2
1038      Double_t dd=0.0063; x0=21.5;
1039      d=dd;
1040      if (TMath::Abs(y-3.08)>0.5) d+=dd;
1041      //if (TMath::Abs(y-3.08)>0.45) d+=dd;
1042      if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
1043   } else
1044   if (3<fR&&fR<5) {   //SPD1
1045      Double_t dd=0.0063; x0=21.5;
1046      d=dd;
1047      if (TMath::Abs(y+0.21)>0.6) d+=dd;
1048      //if (TMath::Abs(y+0.21)>0.45) d+=dd;
1049      if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
1050   }
1051
1052   return d;
1053 }
1054
1055 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
1056 {
1057   //--------------------------------------------------------------------
1058   //Returns the thickness between the current layer and the vertex (units X0)
1059   //--------------------------------------------------------------------
1060   Double_t d=0.0028*3*3; //beam pipe
1061   Double_t x0=0;
1062
1063   Double_t xn=fgLayers[fI].GetR();
1064   for (Int_t i=0; i<fI; i++) {
1065     Double_t xi=fgLayers[i].GetR();
1066     d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
1067   }
1068
1069   if (fI>1) {
1070     Double_t xi=9.;
1071     d+=0.0097*xi*xi;
1072   }
1073
1074   if (fI>3) {
1075     Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
1076     d+=0.0034*xi*xi;
1077   }
1078
1079   return d/(xn*xn);
1080 }
1081
1082 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
1083   //--------------------------------------------------------------------
1084   // This function returns number of clusters within the "window" 
1085   //--------------------------------------------------------------------
1086   Int_t ncl=0;
1087   for (Int_t i=fI; i<fN; i++) {
1088     const AliITSclusterV2 *c=fClusters[i];
1089     if (c->GetZ() > fZmax) break;
1090     if (c->IsUsed()) continue;
1091     const AliITSdetector &det=GetDetector(c->GetDetectorIndex());    
1092     Double_t y=fR*det.GetPhi() + c->GetY();
1093
1094     if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
1095     if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
1096
1097     if (y<fYmin) continue;
1098     if (y>fYmax) continue;
1099     ncl++;
1100   }
1101   return ncl;
1102 }
1103
1104 Bool_t 
1105 AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
1106   //--------------------------------------------------------------------
1107   // This function refits the track "t" at the position "x" using
1108   // the clusters from "c"
1109   //--------------------------------------------------------------------
1110   Int_t index[kMaxLayer];
1111   Int_t k;
1112   for (k=0; k<kMaxLayer; k++) index[k]=-1;
1113   Int_t nc=c->GetNumberOfClusters();
1114   for (k=0; k<nc; k++) { 
1115     Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
1116     index[nl]=idx; 
1117   }
1118
1119   Int_t from, to, step;
1120   if (xx > t->GetX()) {
1121       from=0; to=kMaxLayer;
1122       step=+1;
1123   } else {
1124       from=kMaxLayer-1; to=-1;
1125       step=-1;
1126   }
1127
1128   for (Int_t i=from; i != to; i += step) {
1129      AliITSlayer &layer=fgLayers[i];
1130      Double_t r=layer.GetR();
1131  
1132      {
1133      Double_t hI=i-0.5*step; 
1134      if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {             
1135         Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
1136         Double_t d=0.0034, x0=38.6; 
1137         if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1138         if (!t->PropagateTo(rs,-step*d,x0)) {
1139           return kFALSE;
1140         }
1141      }
1142      }
1143
1144      // remember old position [SR, GSI 18.02.2003]
1145      Double_t oldX=0., oldY=0., oldZ=0.;
1146      if (t->IsStartedTimeIntegral() && step==1) {
1147         t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1148      }
1149      //
1150
1151      Double_t x,y,z;
1152      if (!t->GetGlobalXYZat(r,x,y,z)) { 
1153        return kFALSE;
1154      }
1155      Double_t phi=TMath::ATan2(y,x);
1156      Int_t idet=layer.FindDetectorIndex(phi,z);
1157      if (idet<0) { 
1158        return kFALSE;
1159      }
1160      const AliITSdetector &det=layer.GetDetector(idet);
1161      phi=det.GetPhi();
1162      if (!t->Propagate(phi,det.GetR())) {
1163        return kFALSE;
1164      }
1165      t->SetDetectorIndex(idet);
1166
1167      const AliITSclusterV2 *cl=0;
1168      Double_t maxchi2=kMaxChi2;
1169
1170      Int_t idx=index[i];
1171      if (idx>0) {
1172         const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx); 
1173         if (idet != c->GetDetectorIndex()) {
1174            idet=c->GetDetectorIndex();
1175            const AliITSdetector &det=layer.GetDetector(idet);
1176            if (!t->Propagate(det.GetPhi(),det.GetR())) {
1177              return kFALSE;
1178            }
1179            t->SetDetectorIndex(idet);
1180         }
1181         Double_t chi2=t->GetPredictedChi2(c);
1182         if (chi2<maxchi2) { 
1183           cl=c; 
1184           maxchi2=chi2; 
1185         } else {
1186           return kFALSE;
1187         }
1188      }
1189      /*
1190      if (cl==0)
1191      if (t->GetNumberOfClusters()>2) {
1192         Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1193         Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1194         Double_t zmin=t->GetZ() - dz;
1195         Double_t zmax=t->GetZ() + dz;
1196         Double_t ymin=t->GetY() + phi*r - dy;
1197         Double_t ymax=t->GetY() + phi*r + dy;
1198         layer.SelectClusters(zmin,zmax,ymin,ymax);
1199
1200         const AliITSclusterV2 *c=0; Int_t ci=-1;
1201         while ((c=layer.GetNextCluster(ci))!=0) {
1202            if (idet != c->GetDetectorIndex()) continue;
1203            Double_t chi2=t->GetPredictedChi2(c);
1204            if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1205         }
1206      }
1207      */
1208      if (cl) {
1209        if (!t->Update(cl,maxchi2,idx)) {
1210           return kFALSE;
1211        }
1212        t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1213      }
1214
1215      {
1216      Double_t x0;
1217      Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1218      t->CorrectForMaterial(-step*d,x0);
1219      }
1220                  
1221      // track time update [SR, GSI 17.02.2003]
1222      if (t->IsStartedTimeIntegral() && step==1) {
1223         Double_t newX, newY, newZ;
1224         t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1225         Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) + 
1226                        (oldZ-newZ)*(oldZ-newZ);
1227         t->AddTimeStep(TMath::Sqrt(dL2));
1228      }
1229      //
1230
1231   }
1232
1233   if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1234   return kTRUE;
1235 }
1236
1237
1238 Float_t  *AliITStrackerV2::GetWeight(Int_t index) {
1239   //--------------------------------------------------------------------
1240   //       Return pointer to a given cluster
1241   //--------------------------------------------------------------------
1242   Int_t l=(index & 0xf0000000) >> 28;
1243   Int_t c=(index & 0x0fffffff) >> 00;
1244   return fgLayers[l].GetWeight(c);
1245 }
1246
1247
1248
1249 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1250   //--------------------------------------------------------------------
1251   // This function marks clusters assigned to the track
1252   //--------------------------------------------------------------------
1253   AliTracker::UseClusters(t,from);
1254
1255   AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1256   //if (c->GetQ()>2) c->Use();
1257   if (c->GetSigmaZ2()>0.1) c->Use();
1258   c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1259   //if (c->GetQ()>2) c->Use();
1260   if (c->GetSigmaZ2()>0.1) c->Use();
1261
1262 }
1263
1264
1265 void AliITStrackerV2::AddTrackHypothesys(AliITStrackV2 * track, Int_t esdindex)
1266 {
1267   //------------------------------------------------------------------
1268   // add track to the list of hypothesys
1269   //------------------------------------------------------------------
1270
1271   if (esdindex>=fTrackHypothesys.GetEntriesFast()) fTrackHypothesys.Expand(esdindex*2+10);
1272   //
1273   TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1274   if (!array) {
1275     array = new TObjArray(10);
1276     fTrackHypothesys.AddAt(array,esdindex);
1277   }
1278   array->AddLast(track);
1279 }
1280
1281 void AliITStrackerV2::SortTrackHypothesys(Int_t esdindex, Float_t likelihoodlevel)
1282 {
1283   //-------------------------------------------------------------------
1284   // compress array of track hypothesys
1285   // keep only maxsize best hypothesys
1286   //-------------------------------------------------------------------
1287   if (esdindex>fTrackHypothesys.GetEntriesFast()) return;
1288   if (! (fTrackHypothesys.At(esdindex)) ) return;
1289   TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1290   Int_t entries = array->GetEntriesFast();
1291
1292   Float_t * chi2        = new Float_t[entries];
1293   Float_t * probability = new Float_t[entries];
1294   Float_t sumprobabilityall=0;
1295   Int_t * index         = new Int_t[entries];
1296   //
1297   //
1298   for (Int_t itrack=0;itrack<array->GetEntriesFast();itrack++){
1299     AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1300     //
1301     if (track && track->GetNumberOfClusters()>(track->fNUsed+track->fNSkipped)){
1302       //
1303       chi2[itrack] = track->GetChi2()/(track->GetNumberOfClusters()-track->fNUsed-track->fNSkipped);      
1304       if (track->fESDtrack)
1305         if (track->fESDtrack->GetTPCsignal()>80){
1306           track->CookdEdx();
1307           if ((track->GetdEdx()/track->fESDtrack->GetTPCsignal())<0.4){
1308             Float_t factor = 2.+10.*(0.6-track->GetdEdx()/track->fESDtrack->GetTPCsignal());
1309             chi2[itrack]*= factor;  //mismatch in dEdx
1310           }
1311         }      
1312     }
1313     else
1314       chi2[itrack] = 10000000.;
1315     probability[itrack] = 1./(0.3+chi2[itrack]);
1316     sumprobabilityall+=probability[itrack];
1317   }
1318
1319   TMath::Sort(entries,chi2,index,kFALSE);
1320   TObjArray * newarray = new TObjArray();
1321   Float_t     sumprobability = 0.;
1322   for (Int_t i=0;i<entries;i++){
1323     AliITStrackV2 * track = (AliITStrackV2*)array->At(index[i]);
1324     if (!track) break;
1325     if (chi2[index[i]]<30){
1326       newarray->AddLast(array->RemoveAt(index[i]));      
1327       track->fChi2MIP[0] = chi2[index[i]];     // base chi 2
1328       sumprobability+= probability[index[i]]/sumprobabilityall;
1329       if (sumprobability> likelihoodlevel) break;
1330     }
1331     else{
1332       delete array->RemoveAt(index[i]);
1333     }
1334   }
1335
1336   array->Delete();
1337   delete fTrackHypothesys.RemoveAt(esdindex);
1338   fTrackHypothesys.AddAt(newarray,esdindex);
1339
1340   delete [] chi2;
1341   delete [] index;
1342
1343 }
1344
1345
1346 void AliITStrackerV2::CompressTrackHypothesys(Int_t esdindex, Float_t likelihoodlevel, Int_t maxsize)
1347 {
1348   //
1349   //
1350   if (esdindex>fTrackHypothesys.GetEntriesFast()) return;
1351   if (! (fTrackHypothesys.At(esdindex)) ) return;
1352   TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1353   Int_t entries = array->GetEntriesFast();
1354   //
1355   if (likelihoodlevel>0.000001){
1356     Float_t *probability = new Float_t[entries];
1357     for (Int_t i=0;i<entries;i++) probability[i]=0.;
1358     Float_t  sumprobabilityall=0;
1359     Int_t   *index         = new Int_t[entries];
1360     
1361     for (Int_t itrack=0;itrack<entries;itrack++){
1362       AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1363       probability[itrack]=0;
1364       if (!track) continue;
1365       probability[itrack] = 1./(0.3+track->fChi2MIP[0]);
1366       sumprobabilityall  += probability[itrack];
1367       //    
1368     }
1369     if (sumprobabilityall<=0.000000000000001){
1370       return;
1371     }
1372     for (Int_t itrack=0;itrack<entries;itrack++){
1373       probability[itrack]/=sumprobabilityall;
1374     }
1375     
1376     TMath::Sort(entries, probability, index, kTRUE);
1377     //
1378     TObjArray * newarray = new TObjArray();
1379     Float_t     sumprobability = 0.;
1380     for (Int_t i=0;i<entries;i++){
1381       AliITStrackV2 * track = (AliITStrackV2*)array->At(index[i]);
1382       if (!track) continue;    
1383       newarray->AddLast(array->RemoveAt(index[i]));      
1384       sumprobability+= probability[index[i]];
1385       if (sumprobability> likelihoodlevel) break;
1386       if (i>maxsize) break;
1387     }
1388     
1389     array->Delete();
1390     delete fTrackHypothesys.RemoveAt(esdindex);
1391     fTrackHypothesys.AddAt(newarray,esdindex);
1392     //    
1393     delete []index;
1394     delete []probability;
1395   }
1396   else{
1397     array->Delete();
1398     delete fTrackHypothesys.RemoveAt(esdindex);
1399   }
1400 }
1401
1402
1403 AliITStrackV2 * AliITStrackerV2::GetBestHypothesys(Int_t esdindex, AliITStrackV2 * original, Int_t checkmax)
1404 {
1405   //-------------------------------------------------------------
1406   // try to find best hypothesy
1407   // currently - minimal chi2 of track+backpropagated track+matching to the tpc track
1408   //-------------------------------------------------------------
1409   if (fTrackHypothesys.GetEntriesFast()<=esdindex) return 0;
1410   TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1411   if (!array) return 0;
1412   Int_t entries = array->GetEntriesFast();
1413   if (!entries) return 0;  
1414   Float_t minchi2 = 100000;
1415   Int_t   maxn    = 3;
1416   AliITStrackV2 * besttrack=0;
1417   Int_t accepted =0;
1418   Int_t maxindex=0;
1419   //
1420   Float_t sumz2=0;
1421   Float_t sumy2=0;
1422   Float_t sumall=0;
1423   for (Int_t itrack=0;itrack<entries; itrack++){
1424     AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1425     if (!track) continue;
1426     sumall++;
1427     sumz2+=track->GetSigmaZ2();
1428     sumy2+=track->GetSigmaY2();
1429   }
1430   sumz2/=sumall;
1431   sumy2/=sumall;
1432
1433   Float_t dedxmismatch=1;
1434   for (Int_t i=0;i<entries;i++){    
1435     maxindex = i;
1436     AliITStrackV2 * track = (AliITStrackV2*)array->At(i);    
1437     if (!track) continue;
1438     track->fChi2MIP[1] = 1000000;
1439     track->fChi2MIP[2] = 1000000;
1440
1441     if ( (track->GetNumberOfClusters()-track->fNSkipped-track->fNUsed)<2) continue;
1442     //
1443     if (track->fESDtrack)
1444       if (track->fESDtrack->GetTPCsignal()>80){
1445         track->CookdEdx();      
1446         if ((track->GetdEdx()/track->fESDtrack->GetTPCsignal())<0.4){
1447           //mismatch in dEdx 
1448           dedxmismatch= 2.+ 10.*(0.6-track->GetdEdx()/track->fESDtrack->GetTPCsignal());
1449         }
1450       }
1451
1452     //    track->SetLabel(original->GetLabel());
1453     //CookLabel(track,0.0);
1454     //if (track->GetFakeRatio()>0.01) continue;
1455     //
1456     //
1457     // backtrack
1458     AliITStrackV2 * backtrack = new AliITStrackV2(*track);
1459     backtrack->ResetCovariance();
1460     backtrack->ResetClusters();
1461     Double_t x = original->GetX();
1462     if (!RefitAt(x,backtrack,track)){
1463       delete backtrack;
1464       delete array->RemoveAt(i);
1465       continue;
1466     }
1467     if (  (backtrack->GetChi2() / float(backtrack->GetNumberOfClusters()-track->fNSkipped-track->fNUsed-0.5))>6) 
1468       {
1469         delete backtrack; 
1470         delete array->RemoveAt(i);
1471         continue;
1472       }
1473     Double_t deltac   = backtrack->GetC()-original->GetC();
1474     Double_t deltatgl = backtrack->GetTgl()-original->GetTgl();
1475     //
1476     Double_t poolc2      = (deltac*deltac)/(original->fC44+backtrack->fC44);
1477     Double_t pooltgl2    = (deltatgl*deltatgl)/(original->fC33+backtrack->fC33);
1478     if ((poolc2+pooltgl2)>32){  //4 sigma
1479       delete backtrack; 
1480       delete array->RemoveAt(i);
1481       continue;      
1482     }
1483     //Double_t bpoolc      = (deltac*deltac)/(original->fC44);
1484     //Double_t bpooltgl    = (deltatgl*deltatgl)/(original->fC33);
1485
1486     //
1487     //forward track - without constraint
1488     AliITStrackV2 * forwardtrack = new AliITStrackV2(*original);
1489     //    forwardtrack->ResetCovariance();
1490     forwardtrack->ResetClusters();
1491     x = track->GetX();
1492     if (!RefitAt(x,forwardtrack,track)){
1493       delete forwardtrack;
1494       delete backtrack;
1495       delete array->RemoveAt(i);
1496       continue;
1497     }
1498     if ( (forwardtrack->GetChi2()/float(forwardtrack->GetNumberOfClusters()-track->fNSkipped-track->fNUsed))>6) 
1499       {
1500         delete forwardtrack; 
1501         delete array->RemoveAt(i);
1502         continue;
1503       }
1504     //
1505     accepted++;
1506     if (accepted>checkmax){
1507       delete backtrack;
1508       delete forwardtrack;
1509       break;
1510     }
1511     Double_t chi2 = (backtrack->GetChi2()/(backtrack->GetNumberOfClusters()-1-track->fNSkipped-track->fNUsed)+
1512                      forwardtrack->GetChi2()/(forwardtrack->GetNumberOfClusters()-track->fNSkipped-track->fNUsed));
1513       //                     bpoolc+bpooltgl;
1514     //    chi2 *= (forwardtrack->GetSigmaZ2()/sumz2+forwardtrack->GetSigmaY2()/sumy2);
1515     chi2 *= dedxmismatch;
1516     //
1517     //
1518     track->fChi2MIP[1] = backtrack->GetChi2()/(backtrack->GetNumberOfClusters()-1-track->fNSkipped-track->fNUsed);
1519     track->fChi2MIP[2] = forwardtrack->GetChi2()/(forwardtrack->GetNumberOfClusters()-track->fNSkipped-track->fNUsed);
1520     track->fChi2MIP[3] = poolc2+pooltgl2;
1521     //
1522     
1523     if (track->GetNumberOfClusters()>maxn){
1524       besttrack =  new AliITStrackV2(*forwardtrack);
1525       maxn      =  track->GetNumberOfClusters();
1526       minchi2   =  chi2;
1527       delete backtrack;
1528       delete forwardtrack;
1529       continue;
1530     }   
1531     //
1532     if (chi2 < minchi2){
1533       besttrack = new AliITStrackV2(*forwardtrack);
1534       minchi2   = chi2;      
1535     }    
1536     delete backtrack;
1537     delete forwardtrack;
1538   }
1539   //
1540   //
1541   if (!besttrack || besttrack->GetNumberOfClusters()<4) {
1542     return 0;
1543   }
1544
1545   //
1546   besttrack->SetLabel(original->GetLabel());
1547   CookLabel(besttrack,0.0);
1548   //
1549   // calculate "weight of the cluster"
1550   //
1551
1552   {
1553     //sign usage information for clusters
1554     Int_t clusterindex[6][100];
1555     Double_t clusterweight[6][100];
1556     for (Int_t ilayer=0;ilayer<6;ilayer++)
1557       for (Int_t icluster=0;icluster<100;icluster++){
1558         clusterindex[ilayer][icluster]   = -1;
1559         clusterweight[ilayer][icluster]  = 0;
1560       }
1561     //printf("%d\t%d\n",esdindex, entries);
1562     //
1563     Float_t sumchi2=0;
1564     for (Int_t itrack=0;itrack<entries; itrack++){
1565       AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1566       if (!track) continue;
1567       if (track->fChi2MIP[1]>1000) continue;  //not accepted
1568       sumchi2 +=1./(0.3+track->fChi2MIP[1]+track->fChi2MIP[2]);
1569     }
1570     for (Int_t itrack=0;itrack<entries;itrack++){
1571       AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1572       if (!track) continue;
1573       if (track->fChi2MIP[1]>1000) continue;  //not accepted  
1574       for (Int_t icluster=0;icluster<track->GetNumberOfClusters();icluster++){     
1575         Int_t tindex = track->GetClusterIndex(icluster);
1576         Int_t ilayer = (tindex & 0xf0000000) >> 28;
1577         if (tindex<0) continue;
1578         Int_t cindex =0;
1579         //
1580         for (cindex=0;cindex<100;cindex++){
1581           if (clusterindex[ilayer][cindex]<0) break;
1582           if (clusterindex[ilayer][cindex]==tindex) break;
1583         }
1584         if (cindex>100) break;
1585         if (clusterindex[ilayer][cindex]!=tindex) clusterindex[ilayer][cindex] = tindex;
1586         clusterweight[ilayer][cindex]+= (1./(0.3+track->fChi2MIP[1]+track->fChi2MIP[2]))* (1./sumchi2); 
1587         Float_t *weight = GetWeight(tindex);
1588         
1589         if (weight){
1590           *weight+= (1./(0.3+track->fChi2MIP[1]+track->fChi2MIP[2]))* (1./sumchi2);
1591         }
1592       }
1593     }
1594     
1595     if (besttrack->GetChi2()/besttrack->GetNumberOfClusters()>3.5) return besttrack; //don't sign clusters
1596     Int_t current=0;
1597     Double_t deltad = besttrack->GetD(GetX(),GetY());
1598     Double_t deltaz = besttrack->GetZat(GetX()) - GetZ();
1599     Double_t deltaprim = TMath::Sqrt(deltad*deltad+deltaz*deltaz);
1600
1601     for (Int_t icluster=0;icluster<besttrack->GetNumberOfClusters();icluster++){
1602       Int_t index = besttrack->GetClusterIndex(icluster);
1603       Int_t ilayer =  (index & 0xf0000000) >> 28;
1604       AliITSclusterV2 *c = (AliITSclusterV2*)GetCluster(index);
1605       if (!c) continue;
1606       // 
1607       for (Int_t icluster=0;icluster<100;icluster++){
1608         // sign non "doubt" clusters as used
1609         if (clusterindex[ilayer][icluster]!=index) continue;
1610         //      Float_t * weight = GetWeight(index);
1611         //if (weight) if (*weight>1){
1612         //  if (c->IsUsed()) continue;
1613         //  c->Use();
1614         //}      
1615         if ( (ilayer*0.2+0.2)<deltaprim) continue; // secondaries
1616         if (c->GetNy()>4) continue; // don sign cluster
1617         if ( (ilayer>1&&clusterweight[ilayer][icluster]>0.7) || (ilayer<2&&clusterweight[ilayer][icluster]>0.8) ){
1618           current++;
1619           if (c->IsUsed()) continue;
1620           c->Use();
1621         }
1622       }
1623     }
1624   }
1625  
1626   //
1627   return besttrack;
1628
1629
1630
1631 AliITStrackV2 * AliITStrackerV2::GetBestHypothesysMIP(Int_t esdindex, AliITStrackV2 * original)
1632 {
1633   //-------------------------------------------------------------
1634   // try to find best hypothesy
1635   // currently - minimal chi2 of track+backpropagated track+matching to the tpc track
1636   //-------------------------------------------------------------
1637   if (fTrackHypothesys.GetEntriesFast()<=esdindex) return 0;
1638   TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1639   if (!array) return 0;
1640   Int_t entries = array->GetEntriesFast();
1641   if (!entries) return 0;  
1642   AliITStrackV2 * besttrack=0;
1643   //
1644   //sign usage information for clusters
1645   Int_t clusterindex[6][100];
1646   Double_t clusterweight[6][100];
1647   for (Int_t ilayer=0;ilayer<6;ilayer++)
1648     for (Int_t icluster=0;icluster<100;icluster++){
1649       clusterindex[ilayer][icluster]   = -1;
1650       clusterweight[ilayer][icluster]  = 0;
1651     }
1652   //printf("%d\t%d\n",esdindex, entries);
1653   //
1654   Float_t sumchi2=0;
1655   for (Int_t itrack=0;itrack<entries; itrack++){
1656     AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1657     if (!track) continue;
1658     if (track->fChi2MIP[1]>1000) continue;  //not accepted - before
1659     sumchi2 +=1./(0.3+track->fChi2MIP[1]+track->fChi2MIP[2]);
1660   }
1661   //
1662   // get cluster weight 
1663   for (Int_t itrack=0;itrack<entries;itrack++){
1664     AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1665     if (!track) continue;
1666     if (track->fChi2MIP[1]>1000) continue;  // track not accepted in previous iterration  
1667     for (Int_t icluster=0;icluster<track->GetNumberOfClusters();icluster++){     
1668       Int_t tindex = track->GetClusterIndex(icluster);
1669       Int_t ilayer = (tindex & 0xf0000000) >> 28;
1670       if (tindex<0) continue;
1671       Int_t cindex =0;
1672       //
1673       for (cindex=0;cindex<100;cindex++){
1674         if (clusterindex[ilayer][cindex]<0) break;
1675         if (clusterindex[ilayer][cindex]==tindex) break;
1676       }
1677       if (cindex>100) break;
1678       if (clusterindex[ilayer][cindex]!=tindex) clusterindex[ilayer][cindex] = tindex;
1679       clusterweight[ilayer][cindex]+= (1./(0.3+track->fChi2MIP[1]+track->fChi2MIP[2]))* (1./sumchi2); 
1680     }
1681   }
1682   //  
1683   // get cluster relative sharing - factor
1684   //
1685   // 
1686   for (Int_t ilayer=0;ilayer<6;ilayer++)
1687     for (Int_t cindex=0;cindex<100;cindex++){
1688       if (clusterindex[ilayer][cindex]<0) continue;
1689       Int_t tindex = clusterindex[ilayer][cindex];
1690       Float_t *weight = GetWeight(tindex);
1691       if (!weight){
1692         printf("Problem 1\n");  // not existing track
1693         continue;
1694       }
1695       if (*weight<(clusterweight[ilayer][cindex]-0.00001)){
1696         printf("Problem 2\n");  // not normalized probability
1697         continue;
1698       }   
1699       AliITSclusterV2 *c = (AliITSclusterV2*)GetCluster(tindex);
1700       if (c->GetNy()<5){
1701         clusterweight[ilayer][cindex]/= *weight; 
1702       }
1703       else{
1704         Float_t weight2 = TMath::Max(*weight-0.5*clusterweight[ilayer][cindex],0.0000001);
1705         clusterweight[ilayer][cindex]/= weight2;
1706         if (    clusterweight[ilayer][cindex]>1)   clusterweight[ilayer][cindex] =1.;
1707       }
1708     }
1709   //
1710   //take to the account sharing factor
1711   //
1712   Float_t chi2 =10000000;
1713   Float_t sharefactor=0;
1714   Float_t minchi2 = 100000000;
1715   Float_t secchi2 = 100000000;
1716   Float_t norm=0;
1717   for (Int_t itrack=0;itrack<entries; itrack++){
1718     AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1719     if (!track) continue;
1720     if (track->fChi2MIP[1]>1000) continue;  //not accepted - before
1721     chi2 = track->fChi2MIP[1];
1722     Float_t newchi2=0;
1723     sharefactor=0;
1724     norm =0;
1725     //
1726     for (Int_t icluster=0;icluster<track->GetNumberOfClusters();icluster++){     
1727       Int_t tindex = track->GetClusterIndex(icluster);
1728       Int_t ilayer = (tindex & 0xf0000000) >> 28;
1729       if (tindex<0) continue;
1730       Int_t cindex =0;
1731       Float_t cchi2 = (track->fDy[ilayer]*track->fDy[ilayer])/(track->fSigmaY[ilayer]*track->fSigmaY[ilayer]) +
1732         (track->fDz[ilayer]*track->fDz[ilayer])/(track->fSigmaZ[ilayer]*track->fSigmaZ[ilayer]) ;
1733       //
1734       for (cindex=0;cindex<100;cindex++){
1735         if (clusterindex[ilayer][cindex]<0) break;
1736         if (clusterindex[ilayer][cindex]==tindex) break;
1737       }
1738       if (cindex>100) continue;
1739       if (clusterweight[ilayer][cindex]>0.00001){
1740         sharefactor+= clusterweight[ilayer][cindex];
1741         cchi2/=clusterweight[ilayer][cindex];
1742         norm++;
1743       }
1744       newchi2+=cchi2;
1745     }
1746     newchi2/=(norm-track->fNSkipped-track->fNUsed);
1747     track->fChi2MIP[4] = newchi2;
1748     //if (norm>0) sharefactor/=norm;
1749     //if (sharefactor<0.5) return 0;
1750     //    chi2*=1./(0.5+sharefactor);
1751     if (newchi2<minchi2){
1752       besttrack = track;
1753       minchi2   = newchi2;
1754     }
1755     else{
1756       if (newchi2<secchi2){
1757         secchi2 = newchi2;
1758       }
1759     }
1760     //
1761   }
1762
1763   //
1764   //
1765   if (!besttrack || besttrack->GetNumberOfClusters()<4) {
1766     return 0;
1767   }
1768   //
1769   if ((minchi2/secchi2)<0.7){
1770     //
1771     //increase the weight for clusters if the probability of other hypothesys is small
1772     Double_t deltad = besttrack->GetD(GetX(),GetY());
1773     Double_t deltaz = besttrack->GetZat(GetX()) - GetZ();
1774     Double_t deltaprim = TMath::Sqrt(deltad*deltad+deltaz*deltaz);
1775     //
1776     for (Int_t icluster=0;icluster<besttrack->GetNumberOfClusters();icluster++){
1777       Int_t index = besttrack->GetClusterIndex(icluster);
1778       Int_t ilayer =  (index & 0xf0000000) >> 28;
1779       AliITSclusterV2 *c = (AliITSclusterV2*)GetCluster(index);
1780       if (!c) continue; 
1781       if (c->GetNy()>3) continue;
1782       if ( (ilayer*0.2+0.2)<deltaprim) continue; // secondaries     
1783       Float_t * weight = GetWeight(index);
1784       *weight*=2.;
1785       *weight+=1.;
1786     }   
1787   }
1788   
1789
1790   besttrack->SetLabel(original->GetLabel());
1791   CookLabel(besttrack,0.0);
1792   
1793   //
1794   return besttrack;
1795
1796