]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerSA.cxx
Updated flags
[u/mrichter/AliRoot.git] / ITS / AliITStrackerSA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, 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 /* $Id$ */
17
18 ////////////////////////////////////////////////////
19 //  Stand alone tracker class                     //
20 //  Origin:  Elisabetta Crescio                   //
21 //  e-mail:  crescio@to.infn.it                   //
22 //  tracks are saved as AliITStrackV2 objects     //
23 ////////////////////////////////////////////////////
24
25 #include <stdlib.h>
26
27 #include <TArrayI.h>
28 #include <TBranch.h>
29 #include <TObjArray.h>
30 #include <TTree.h>
31
32 #include "AliESD.h"
33 #include "AliESDVertex.h"
34 #include "AliESDtrack.h"
35 #include "AliITSRiemannFit.h"
36 #include "AliITSVertexer.h"
37 #include "AliITSclusterTable.h"
38 #include "AliITSclusterV2.h"
39 #include "AliITSgeom.h"
40 #include "AliITStrackSA.h"
41 #include "AliITStrackerSA.h"
42 #include "AliRun.h"
43
44 ClassImp(AliITStrackerSA)
45
46 //____________________________________________________________________________
47 AliITStrackerSA::AliITStrackerSA():AliITStrackerMI(){
48   // Default constructor
49   Init();
50  
51 }
52 //____________________________________________________________________________
53 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom):AliITStrackerMI(geom) 
54 {
55   // Standard constructor (Vertex is known and passed to this obj.)
56   Init();
57   fVert = 0;
58   fGeom = geom;
59  
60 }
61
62 //____________________________________________________________________________
63 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliESDVertex *vert):AliITStrackerMI(geom) 
64 {
65   // Standard constructor (Vertex is known and passed to this obj.)
66   Init();
67   fVert = vert;
68   fGeom = geom;
69  
70 }
71
72 //____________________________________________________________________________
73 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerMI(geom) 
74 {
75   // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
76   Init();
77   fVertexer = vertexer;
78   fGeom = geom;
79  
80 }
81
82 //____________________________________________________________________________
83 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA& tracker):AliITStrackerMI(){
84   // Copy constructor
85   fPhiEstimate = tracker.fPhiEstimate;
86   for(Int_t i=0;i<2;i++){
87     fPoint1[i]=tracker.fPoint1[i];
88     fPoint2[i]=tracker.fPoint2[i];
89     fPoint3[i]=tracker.fPoint3[i];
90     fPointc[i]=tracker.fPointc[i];
91   }
92   fLambdac = tracker.fLambdac;
93   fPhic = tracker.fPhic;
94   fCoef1 = tracker.fCoef1;
95   fCoef2 = tracker.fCoef2;
96   fCoef3 = tracker.fCoef3;
97   fNloop = tracker.fNloop;
98   fPhiWin = tracker.fPhiWin;
99   fLambdaWin = tracker.fLambdaWin;
100   if(tracker.fVertexer && tracker.fVert){
101     fVert = new AliESDVertex(*tracker.fVert);
102   }
103   else {
104     fVert = tracker.fVert;
105   }
106   fVertexer = tracker.fVertexer;
107   fGeom = tracker.fGeom;
108   fTable = tracker.fTable;
109   fListOfTracks = tracker.fListOfTracks;
110 }
111 //______________________________________________________________________
112 AliITStrackerSA& AliITStrackerSA::operator=(const AliITStrackerSA& /*source*/){
113     // Assignment operator. This is a function which is not allowed to be
114     // done.
115     Error("operator=","Assignment operator not allowed\n");
116     return *this; 
117 }
118
119 //____________________________________________________________________________
120 AliITStrackerSA::~AliITStrackerSA(){
121   // destructor
122   // if fVertexer is not null, the AliESDVertex obj. is owned by this class
123   // and is deleted here
124   if(fVertexer){
125     if(fVert)delete fVert;
126   }
127   fVert = 0;
128   fVertexer = 0;
129   fGeom = 0;
130   if(fPhiWin)delete []fPhiWin;
131   if(fLambdaWin)delete []fLambdaWin;
132   fTable =0;
133   fListOfTracks->Delete();
134   }
135
136 //____________________________________________________________________________
137 void AliITStrackerSA::Init(){
138   //  Reset all data members
139     fPhiEstimate=0;
140     for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
141     fLambdac=0;
142     fPhic=0;
143     fCoef1=0;
144     fCoef2=0;
145     fCoef3=0;
146     fPointc[0]=0;
147     fPointc[1]=0;
148     fVert = 0;
149     fVertexer = 0;
150     fGeom = 0;
151     SetWindowSizes();
152     fTable = 0;
153     fITSclusters = 0;
154     SetSixPoints();
155     fListOfTracks=new TObjArray(0,0);
156  }
157 //_______________________________________________________________________
158 void AliITStrackerSA::ResetForFinding(){
159   //  Reset data members used in all loops during track finding
160     fPhiEstimate=0;
161     for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
162     fLambdac=0;
163     fPhic=0;
164     fCoef1=0;
165     fCoef2=0;
166     fCoef3=0;
167     fPointc[0]=0;
168     fPointc[1]=0;
169     fListOfTracks->Delete();
170 }
171 //____________________________________________________________________________
172 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
173
174   /**************************************************************************
175    * This function finds primary tracks.
176    *                                                                        * 
177    *                                                                        *
178    * Example: to execute function with only the ITS (no combined tracking   *
179    *          with TPC+ITS) and requiring 5/6 points to define a good track *
180    *          call SetSixPoinbts(kFALSE) in advance and then                *
181    *          use: FindTracks(treein,treeout,evnumber)                      *
182    *          to execute combined tracking, before using FindTracks, use    *
183    *          UseFoundTracksV2                                              *
184    *************************************************************************/
185
186   if(!fITSclusters){
187     Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
188     exit(1);
189   }
190   //Get primary vertex
191   if(fVertexer){
192     if(fVert)delete fVert;
193     fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
194   }
195   else {
196     gAlice->GetEvent(evnumber);
197     if(!fVert){
198       Fatal("FindTracks","Vertex is missing\n");
199       return;
200     }
201   }
202   Double_t primaryVertex[3];
203   Double_t errorsprimvert[3];
204   fVert->GetXYZ(primaryVertex);
205   fVert->GetSigmaXYZ(errorsprimvert);
206   if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
207     Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
208     errorsprimvert[0]=0.0001;
209     errorsprimvert[1]=0.0001;
210   }
211   fVert->PrintStatus();
212
213
214   //Fill array with cluster indices for each module
215   if(!fTable){
216     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
217     fTable->FillArray(fITSclusters);
218     fTable->FillArrayCoorAngles(); 
219   }
220
221    
222 //Fill tree for found tracks
223   AliITStrackV2* outrack=0;
224   TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
225   if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
226   else branch->SetAddress(&outrack);
227
228   
229   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
230   for(Int_t i=0;i<fGeom->GetNlayers();i++){
231     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
232   }
233   // firstmod [i] number of the first module in the ITS layer i.
234     
235   AliITSlayer &layer=fgLayers[0];   // first layer
236   Int_t ntrack=0;
237   Int_t dim=layer.GetNumberOfClusters();
238   //loop on the different windows
239   for(Int_t nloop=0;nloop<fNloop;nloop++){
240     for(Int_t ncl=0;ncl<dim;ncl++){ 
241       //loop starting from layer 0
242       ResetForFinding();
243       Int_t pflag=0;
244       AliITSclusterV2* cl = layer.GetCluster(ncl);
245       if(cl->IsUsed()==1) continue;
246       if(cl->TestBit(kSAflag)==kTRUE) continue;
247       
248       fPhic = fTable->GetPhiCluster(0,ncl);
249       fLambdac = fTable->GetLambdaCluster(0,ncl);
250       fPhiEstimate = fPhic;
251       AliITStrackSA* trs = new AliITStrackSA();      
252       fPoint1[0]=primaryVertex[0];
253       fPoint1[1]=primaryVertex[1];
254       fPoint2[0]=fTable->GetXCluster(0,ncl);
255       fPoint2[1]=fTable->GetYCluster(0,ncl);
256       
257       Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
258       for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
259       nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
260       nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
261
262       if(nn[1]>0){
263         pflag=1;
264         fPoint3[0] = fPointc[0];
265         fPoint3[1] = fPointc[1];
266       }
267       nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
268       if(nn[1]==0 && nn[2]==0) pflag=0;
269       if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
270       if(nn[2]!=0 && nn[1]==0){
271         pflag=1;
272         fPoint3[0]=fPointc[0];
273         fPoint3[1]=fPointc[1];
274       }
275
276       nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
277       pflag=1;
278       if(nn[3]!=0) UpdatePoints();
279       nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag); 
280       pflag=1;
281       if(nn[4]!=0) UpdatePoints();
282       nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag); 
283           
284
285       Int_t layOK=0;
286       Int_t numberofpoints;
287       if(fSixPoints) numberofpoints=6;  //check of the candidate track
288       else numberofpoints=5;           //if track is good (with the required number        
289       for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){    //of points) it is written on file
290         if(nn[nnp]!=0) layOK+=1;
291       }
292       if(layOK>=numberofpoints){
293         AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
294         if(tr2==0){
295           Int_t nct = trs->GetNumberOfClustersSA();
296           while(nct--){
297             Int_t index = trs->GetClusterIndexSA(nct);
298             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
299             if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
300           }
301           delete [] nn;
302           delete trs;
303           continue;
304         }
305         outrack=tr2;  
306         out->Fill();
307         ntrack++;
308         Int_t nct = tr2->GetNumberOfClusters();
309
310         while(nct--){
311           Int_t index = tr2->GetClusterIndex(nct);     
312           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);      
313           kl->SetBit(kSAflag);
314
315         }
316       } 
317       else{
318         Int_t nct = trs->GetNumberOfClustersSA();
319         while(nct--){
320           Int_t index = trs->GetClusterIndexSA(nct);
321           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
322           if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
323         }
324       }
325       delete [] nn;
326       delete trs;
327
328     }//end loop on clusters of layer1
329
330   }//end loop2
331
332   //if 5/6 points are required, second loop starting 
333   //from second layer, to find tracks with point of 
334   //layer 1 missing
335    
336   if(!fSixPoints){
337     //   counter for clusters on each layer  
338     Int_t * nn = new Int_t[fGeom->GetNlayers()-1];      
339     for(Int_t nloop=0;nloop<fNloop;nloop++){
340       AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
341       Int_t ncl2=layer2.GetNumberOfClusters();
342       while(ncl2--){ //loop starting from layer 2
343         ResetForFinding();
344         Int_t pflag=0;
345         AliITSclusterV2* cl = layer2.GetCluster(ncl2);
346         if(cl->IsUsed()==1) continue;
347         if(cl->TestBit(kSAflag)==kTRUE) continue;
348         fPhic = fTable->GetPhiCluster(1,ncl2);
349         fLambdac = fTable->GetLambdaCluster(1,ncl2);
350         fPhiEstimate = fPhic;
351         AliITStrackSA* trs = new AliITStrackSA(); 
352      
353         fPoint1[0]=primaryVertex[0];
354         fPoint1[1]=primaryVertex[1];
355         fPoint2[0]=fTable->GetXCluster(1,ncl2);;
356         fPoint2[1]=fTable->GetYCluster(1,ncl2);;
357    
358         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
359         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
360           nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
361                    trs,primaryVertex[2],pflag);
362           if(nn[kk]==0)break;
363           if(kk>0){
364             UpdatePoints();
365             pflag = 1;
366           }
367         }
368         Int_t fl=0;
369         for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
370           if(nn[nnp]!=0) fl+=1;
371         }
372         if(fl>=5){  // 5/6       
373           AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
374           if(tr2==0){
375             Int_t nct = trs->GetNumberOfClustersSA();
376             while(nct--){
377               Int_t index = trs->GetClusterIndexSA(nct);
378               AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
379               if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
380             }
381             delete trs;
382             continue;
383           }
384           outrack=tr2;
385           out->Fill();
386           Int_t nct = tr2->GetNumberOfClusters();
387           while(nct--){
388             Int_t index = tr2->GetClusterIndex(nct);     
389             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
390             if(kl==0) continue;
391             kl->SetBit(kSAflag);
392           }
393         }       
394         else{
395           Int_t nct = trs->GetNumberOfClustersSA();
396           while(nct--){
397             Int_t index = trs->GetClusterIndexSA(nct);
398             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
399             if(kl==0) continue;
400             if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
401           }
402         }
403         delete trs;
404       }//end loop on clusters of layer2
405     }
406     delete [] nn;
407   }  // if(!fSixPoints....  
408
409   delete [] firstmod;
410   delete fTable; fTable=0;
411 }
412
413
414 //______________________________________________________________________
415 Int_t AliITStrackerSA::FindTracks(AliESD* event){
416
417   // Track finder using the ESD object
418
419
420   if(!fITSclusters){
421     Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
422     return -1;
423   }
424  
425   //Get primary vertex
426   Double_t errorsprimvert[3];
427   Double_t primaryVertex[3];
428   event->GetVertex()->GetXYZ(primaryVertex);
429   event->GetVertex()->GetSigmaXYZ(errorsprimvert);
430
431   if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
432     //    Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
433     errorsprimvert[0]=0.005;
434     errorsprimvert[1]=0.005;
435   }
436   //Fill array with cluster indices for each module
437   if(!fTable){
438     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
439     fTable->FillArray(fITSclusters);
440     fTable->FillArrayCoorAngles();
441   }
442   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
443   for(Int_t i=0;i<fGeom->GetNlayers();i++){
444     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
445   }
446   // firstmod [i] number of the first module in the ITS layer i.
447   
448     
449   AliITSlayer &layer=fgLayers[0];   
450   Int_t ntrack=0;
451    Int_t dim=layer.GetNumberOfClusters();
452   //loop on the different windows
453   for(Int_t nloop=0;nloop<fNloop;nloop++){
454     for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
455
456       ResetForFinding();
457       Int_t pflag=0;
458       AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
459       if(cl==0) continue;
460       if(cl->IsUsed()==1) continue;
461       if(cl->TestBit(kSAflag)==kTRUE) continue;
462       if (cl->GetQ()<=0) continue;
463       
464       fPhic = fTable->GetPhiCluster(0,ncl);
465       fLambdac = fTable->GetLambdaCluster(0,ncl);
466
467       if (TMath::Abs(fLambdac)>0.26*TMath::Pi()) continue;
468
469       fPhiEstimate = fPhic;
470       AliITStrackSA* trs = new AliITStrackSA();      
471       fPoint1[0]=primaryVertex[0];
472       fPoint1[1]=primaryVertex[1];
473       fPoint2[0]=fTable->GetXCluster(0,ncl);
474       fPoint2[1]=fTable->GetYCluster(0,ncl);
475       Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
476       for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
477       nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
478           
479       nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
480       if(nn[1]>0){
481         pflag=1;
482         fPoint3[0] = fPointc[0];
483         fPoint3[1] = fPointc[1];
484       }
485       nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
486       if(nn[1]==0 && nn[2]==0) pflag=0;
487       if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
488       if(nn[2]!=0 && nn[1]==0){
489         pflag=1;
490         fPoint3[0]=fPointc[0];
491         fPoint3[1]=fPointc[1];
492       }
493
494       nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag);
495       pflag=1;
496       if(nn[3]!=0) UpdatePoints();
497       nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag); 
498       pflag=1;
499       if(nn[4]!=0) UpdatePoints();
500       nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag); 
501           
502
503       Int_t layOK=0;
504       Int_t numberofpoints;
505       if(fSixPoints) numberofpoints=6;  //check of the candidate track
506       else numberofpoints=5;           //if track is good (with the required number        
507       for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){    //of points) it is written on file
508         if(nn[nnp]!=0) layOK+=1;
509       }
510       if(layOK>=numberofpoints){
511         AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
512         if(tr2==0){
513           Int_t nct = trs->GetNumberOfClustersSA();
514           while(nct--){
515             Int_t index = trs->GetClusterIndexSA(nct);
516             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
517             if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
518             
519           }
520           continue;
521         }
522         
523         AliESDtrack outtrack;
524         outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
525         event->AddTrack(&outtrack);
526         ntrack++;
527         Int_t nct = tr2->GetNumberOfClusters();
528         while(nct--){
529           Int_t index = tr2->GetClusterIndex(nct);     
530           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);      
531           kl->SetBit(kSAflag);
532           
533         } 
534       }
535       else{
536         Int_t nct = trs->GetNumberOfClustersSA();
537         while(nct--){
538           Int_t index = trs->GetClusterIndexSA(nct);
539           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
540           if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
541           
542         }
543       }
544       delete trs;
545       delete[] nn;
546       
547     }//end loop on clusters of layer1
548     
549   }//end loop2
550
551
552
553   //if 5/6 points are required, second loop starting 
554   //from second layer, to find tracks with point of 
555   //layer 1 missing
556    
557   if(!fSixPoints){
558     //   counter for clusters on each layer  
559     Int_t * nn = new Int_t[fGeom->GetNlayers()-1];      
560     for(Int_t nloop=0;nloop<fNloop;nloop++){
561       AliITSlayer &layer2=fgLayers[1]; 
562       Int_t ncl2=layer2.GetNumberOfClusters();
563       while(ncl2--){ //loop starting from layer 2
564         ResetForFinding();
565         Int_t pflag=0;
566         AliITSclusterV2* cl = layer2.GetCluster(ncl2);
567         if(cl->IsUsed()==1) continue;
568         if(cl->TestBit(kSAflag)==kTRUE) continue;
569         fPhic = fTable->GetPhiCluster(1,ncl2);
570         fLambdac = fTable->GetLambdaCluster(1,ncl2);
571         fPhiEstimate = fPhic;
572         AliITStrackSA* trs = new AliITStrackSA(); 
573         fPoint1[0]=primaryVertex[0];
574         fPoint1[1]=primaryVertex[1];
575         fPoint2[0]=fTable->GetXCluster(1,ncl2);
576         fPoint2[1]=fTable->GetYCluster(1,ncl2);
577    
578         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
579         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
580           nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
581                    trs,primaryVertex[2],pflag);
582           if(nn[kk]==0)break;
583           if(kk>0){
584             UpdatePoints();
585             pflag = 1;
586           }
587         }
588         Int_t fl=0;
589         for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
590           if(nn[nnp]!=0) fl+=1;
591         }
592         if(fl>=5){  // 5/6       
593           AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
594           if(tr2==0){
595             Int_t nct = trs->GetNumberOfClustersSA();
596             while(nct--){
597               Int_t index = trs->GetClusterIndexSA(nct);
598               AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
599               if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
600               
601             }
602
603             continue;
604           }
605
606           AliESDtrack outtrack;
607           outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
608           event->AddTrack(&outtrack);
609           ntrack++;
610           Int_t nct = tr2->GetNumberOfClusters();
611           while(nct--){
612             Int_t index = tr2->GetClusterIndex(nct);     
613             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
614             if(kl==0) continue;
615             kl->SetBit(kSAflag);
616           }
617         }       
618         else{
619           Int_t nct = trs->GetNumberOfClustersSA();
620           while(nct--){
621             Int_t index = trs->GetClusterIndexSA(nct);
622             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
623             if(kl==0) continue;
624             if(kl->TestBit(kSAflag)==kTRUE) kl->ResetBit(kSAflag);
625             
626           }
627         }
628         delete trs;
629       }//end loop on clusters of layer2
630     }
631     delete [] nn;
632   }  //end opt="5/6"  
633
634   delete [] firstmod;
635   delete fTable;fTable=0;   
636   Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
637   return 0;
638
639 }
640
641
642 //________________________________________________________________________
643
644 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
645   //fit of the found track
646
647   
648   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
649   for(Int_t i=0;i<fGeom->GetNlayers();i++){
650     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
651   }  
652
653   Int_t nclusters = tr->GetNumberOfClustersSA();
654   TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
655   for(Int_t i=0;i<fGeom->GetNlayers();i++){
656     listlayer[i] = new TObjArray(0,0);
657   }
658
659   TArrayI clind0(20);
660   TArrayI clind1(20);
661   TArrayI clind2(20);
662   TArrayI clind3(20);
663   TArrayI clind4(20);
664   TArrayI clind5(20);
665
666   Int_t * nnn = new Int_t[fGeom->GetNlayers()];
667   for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
668   
669   for(Int_t ncl=0;ncl<nclusters;ncl++){
670     Int_t index = tr->GetClusterIndexSA(ncl);   
671     AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
672
673     if(cl->TestBit(kSAflag)==kTRUE) cl->ResetBit(kSAflag);
674     Int_t lay = (index & 0xf0000000) >> 28;
675     if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
676     if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
677     if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
678     if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
679     if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
680     if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}    
681   }
682   delete [] nnn;
683
684  
685   Int_t * end = new Int_t[fGeom->GetNlayers()];
686   for(Int_t i=0;i<fGeom->GetNlayers();i++){
687     if(listlayer[i]->GetEntries()==0) end[i]=1;
688     else end[i]=listlayer[i]->GetEntries();
689   }
690
691   for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
692     AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1); 
693     TVector3** recp = new TVector3*[3];
694     TVector3** errs = new TVector3*[3];
695     recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
696     errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
697     Double_t x1,y1,z1,sx1,sy1,sz1;
698     Double_t x2,y2,z2,sx2,sy2,sz2;
699     AliITSclusterV2* p1=0;
700     AliITSclusterV2* p2=0;
701     Int_t index1=clind0[l1];
702     Int_t index2=0;
703     for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
704       AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2); 
705       index2=clind1[l2];
706       for(Int_t l3=0;l3<end[2];l3++){  //loop on layer 3
707         AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
708
709         if(cl0==0 && cl1!=0) {
710           p2 = cl2;index1=clind2[l3];
711           p1=cl1;
712                 
713         }
714         if(cl0!=0 && cl1==0){
715           p1=cl0;
716           p2=cl2;index2=clind2[l3];
717         }
718         if(cl0!=0 && cl1!=0){
719           p1=cl0;
720           p2=cl1;
721         }
722         Int_t lay1=(index1 & 0xf0000000) >> 28;
723         Int_t cln1=(index1 & 0x0fffffff) >> 00;
724         Int_t lay2=(index2 & 0xf0000000) >> 28;
725         Int_t cln2=(index2 & 0x0fffffff) >> 00;
726         x1 = fTable->GetXCluster(lay1,cln1);
727         x2 = fTable->GetXCluster(lay2,cln2);
728         y1 = fTable->GetYCluster(lay1,cln1);
729         y2 = fTable->GetYCluster(lay2,cln2);
730         z1 = fTable->GetZCluster(lay1,cln1);
731         z2 = fTable->GetZCluster(lay2,cln2);
732         sx1 = fTable->GetXClusterError(lay1,cln1);
733         sx2 = fTable->GetXClusterError(lay2,cln2);
734         sy1 = fTable->GetYClusterError(lay1,cln1);
735         sy2 = fTable->GetYClusterError(lay2,cln2);
736         sz1 = fTable->GetZClusterError(lay1,cln1);
737         sz2 = fTable->GetZClusterError(lay2,cln2);
738         Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
739         Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
740         recp[1] = new TVector3(x1,y1,z1);
741         errs[1] = new TVector3(sx1,sy1,sz1);
742         recp[2] = new TVector3(x2,y2,z2);
743         errs[2] = new TVector3(sx2,sy2,sz2);
744         
745         //fit on the Riemann sphere
746         Float_t seed1,seed2,seed3;
747         AliITSRiemannFit fit;
748         Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
749         if(rf==0) {
750           for(Int_t i=1;i<3;i++){
751             delete recp[i];
752             delete errs[i];
753           }     
754           continue;
755         }  
756         Double_t phi=seed1;
757         Double_t tgl=seed2;
758         
759         if(phi1>0){ 
760           if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
761             phi=seed1+1.5*TMath::Pi();
762             tgl=seed2; 
763           }
764           if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
765             phi=seed1+0.5*TMath::Pi();
766             tgl=(-1)*seed2; 
767           }
768           if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
769             phi=seed1-0.5*TMath::Pi();
770             tgl=seed2; 
771           }
772         }      
773         if(phi1<0){
774           if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
775             phi=seed1+0.5*TMath::Pi();
776             tgl=(-1)*seed2; 
777           }
778           if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
779             phi=seed1-0.5*TMath::Pi();
780             tgl=seed2; 
781           }
782           if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
783             phi=seed1-1.5*TMath::Pi();
784             tgl=(-1)*seed2; 
785             }
786           }
787         
788         Int_t layer,ladder,detector;
789         fGeom->GetModuleId(module1,layer,ladder,detector);
790         Float_t yclu1 = p1->GetY();
791         Float_t zclu1 = p1->GetZ();
792         Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
793               
794         for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4   
795           AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
796           for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
797             AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
798             for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6  
799               AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
800               AliITStrackSA* trac = new AliITStrackSA(fGeom,layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
801                               
802               if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
803               if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
804               if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
805               if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
806               if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
807               if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
808             
809               //fit with Kalman filter using AliITStrackerMI::RefitAt()
810           
811               AliITStrackMI* ot = new AliITStrackSA(*trac);
812               
813               ot->ResetCovariance();
814               ot->ResetClusters();
815               
816               if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
817
818                 AliITStrackMI *otrack2 = new AliITStrackMI(*ot);
819                 otrack2->ResetCovariance(); 
820                 otrack2->ResetClusters();
821                 //fit from layer 6 to layer 1
822                 if(RefitAt(3.7,otrack2,ot)) {
823                   fListOfTracks->AddLast(otrack2);
824                 } else {
825                   delete otrack2;
826                 }
827                               
828               }       
829           
830               delete ot;
831               delete trac;
832             }//end loop layer 6
833           }//end loop layer 5
834         }//end loop layer 4
835         
836         for(Int_t i=1;i<3;i++){
837           delete recp[i];
838           delete errs[i];
839         }
840       }//end loop layer 3
841     }//end loop layer 2 
842     delete recp[0];
843     delete errs[0];
844     delete[] recp;
845     delete[] errs;    
846   }//end loop layer 1
847
848   delete [] end;
849
850   Int_t dim=fListOfTracks->GetEntries();
851   if(dim==0){
852     for(Int_t i=0;i<fGeom->GetNlayers();i++){
853       delete listlayer[i];
854     }
855     delete [] listlayer;
856     delete [] firstmod;
857     return 0;
858   }
859
860   AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
861
862   if(otrack==0) {
863     for(Int_t i=0;i<fGeom->GetNlayers();i++){
864       delete listlayer[i];
865     }
866     delete [] listlayer; 
867     delete [] firstmod;
868     return 0;
869   }
870   Int_t * indexc = new Int_t[fGeom->GetNlayers()];
871   for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
872   for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
873     indexc[nind] = otrack->GetClusterIndex(nind);
874   }      
875   AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
876   AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);     
877   AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);     
878   AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
879   AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
880   Int_t labl[3]={-1,-1,-1};
881   if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
882     AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
883     labl[0]=cl5->GetLabel(0);
884     labl[1]=cl5->GetLabel(1);
885     labl[2]=cl5->GetLabel(2);
886   }
887   delete [] indexc;
888   if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
889     labl[0]=-1;
890     labl[1]=-1;
891     labl[2]=-1;
892   }
893   Int_t numberofpoints;
894   if(fSixPoints) numberofpoints=6;
895   else numberofpoints=5;
896   CookLabel(otrack,0.); //MI change - to see fake ratio
897   Int_t label =  Label(cl0->GetLabel(0),cl1->GetLabel(0), 
898                        cl2->GetLabel(0),cl3->GetLabel(0),
899                        cl4->GetLabel(0),labl[0],
900                        cl0->GetLabel(1),cl1->GetLabel(1),
901                        cl2->GetLabel(1),cl3->GetLabel(1),
902                        cl4->GetLabel(1),labl[1],
903                        cl0->GetLabel(2),cl1->GetLabel(2),
904                        cl2->GetLabel(2),cl3->GetLabel(2),
905                        cl4->GetLabel(2),labl[2],numberofpoints);
906   
907   otrack->SetLabel(label);  
908   for(Int_t i=0;i<fGeom->GetNlayers();i++){
909     delete listlayer[i];
910   }
911   delete [] listlayer; 
912   delete [] firstmod;
913   return otrack;
914
915 }
916
917 //_______________________________________________________________________
918 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
919   // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
920   //(or AliITStrackV2 tracks found with function FindTracks of this class)
921   
922
923   //Get primary vertex
924   if(fVertexer){
925     if(fVert)delete fVert;
926     fVert = fVertexer->FindVertexForCurrentEvent(evnum);
927   }
928   else {
929     gAlice->GetEvent(evnum);
930     if(!fVert){
931       Fatal("FindTracks","Vertex is missing\n");
932       return;
933     }
934   }
935   Double_t primaryVertex[3];
936   fVert->GetXYZ(primaryVertex);
937
938   if(!fTable){
939     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
940     fTable->FillArray(fITSclusters);
941     fTable->FillArrayCoorAngles();
942   }
943
944   TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
945   if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
946   AliITStrackV2* ttrrt = new AliITStrackV2;
947   bra->SetAddress(&ttrrt);
948
949   for(Int_t nj=0;nj<treev2->GetEntries();nj++){
950     treev2->GetEvent(nj);
951     Int_t ncl = ttrrt->GetNumberOfClusters();
952     for(Int_t k=0;k<ncl;k++){
953       Int_t index = ttrrt->GetClusterIndex(k);
954       AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
955       if(clui->IsUsed()==0) clui->Use();
956       
957     }
958   }
959   delete ttrrt;
960   
961 }
962
963 //_______________________________________________________________________
964 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
965   // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
966
967   //Get primary vertex
968
969   Double_t primaryVertex[3];
970   event->GetVertex()->GetXYZ(primaryVertex);
971
972   if(!fTable){
973     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
974     fTable->FillArray(fITSclusters);
975     fTable->FillArrayCoorAngles(); 
976   }
977
978   Int_t ntracks = event->GetNumberOfTracks();
979   while (ntracks--) {
980     AliESDtrack *esd=event->GetTrack(ntracks);
981     if ((esd->GetStatus()&
982          AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue; 
983     UInt_t idx[6];
984     Int_t ncl = esd->GetITSclusters(idx);
985     for(Int_t  clu=0; clu<ncl; clu++){
986       AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
987       if(cl->IsUsed()==0) cl->Use();
988     }
989   }
990   
991   Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
992   
993
994 }
995
996 /*
997 //_______________________________________________________
998 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
999   //function used to to find the clusters associated to the track
1000   Int_t nc=0;
1001   AliITSlayer &lay = fgLayers[layer];
1002   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
1003   for(Int_t i=0;i<fGeom->GetNlayers();i++){
1004     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
1005   }
1006   if(pflag==1){
1007       
1008     Float_t cx1,cx2,cy1,cy2;
1009     FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1010     Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
1011     if(fun==0) {
1012       delete[] firstmod;
1013       return 0;
1014     }
1015
1016     Double_t fi1 =TMath::ATan2(cy1,cx1);
1017     Double_t fi2 =TMath::ATan2(cy2,cx2);
1018     fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1019   }
1020
1021   Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1022   Double_t zed1  =  TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1023   Double_t zed2  =  TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1024   
1025   Double_t fi = fPhiEstimate;
1026   Int_t nmod  = lay.FindDetectorIndex(fi,zed);
1027   if (nmod < 0) {
1028     delete[] firstmod;
1029     return 0;
1030   }
1031   nmod += firstmod[layer];
1032  
1033   Int_t nm[8]={0,0,0,0,0,0,0,0};
1034   nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1035   nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1036   nm[2] = lay.FindDetectorIndex(fi,zed1);
1037   nm[3] = lay.FindDetectorIndex(fi,zed2);
1038   nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1039   nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1040   nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1041   nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1042
1043
1044   Int_t nn=0;
1045   TArrayI* array =(TArrayI*)fTable->GetListOfClusters(nmod);
1046   TArrayI* listc = new TArrayI(array->GetSize());
1047   for(Int_t i=0;i<array->GetSize();i++){
1048     Int_t in=(Int_t)array->At(i);
1049     listc->AddAt(in,nn);
1050     nn++;
1051   }
1052     
1053   Int_t k=0;
1054   Int_t val;
1055   while(k<8){
1056     for(Int_t h=k+1;h<8;h++){
1057       if(nm[k]>nm[h]){
1058         val=nm[k];
1059         nm[k]=nm[h];
1060         nm[h]=val;
1061       }
1062      
1063    }
1064     k++;
1065   }
1066  
1067   Int_t value=-5;
1068  
1069   for(Int_t ii=0;ii<8;ii++){
1070     if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1071       TArrayI* ar =(TArrayI*)fTable->GetListOfClusters(nm[ii]+firstmod[layer]);
1072       listc->Set(listc->GetSize()+ar->GetSize());
1073       for(Int_t j=0;j<ar->GetSize();j++){
1074         Int_t in=(Int_t)ar->At(j);
1075         listc->AddAt(in,nn);
1076         nn++;
1077         value=nm[ii];
1078       }
1079     }
1080   }
1081  
1082   
1083   for(Int_t i=0;i<listc->GetSize();i++){
1084     Int_t index = (Int_t)listc->At(i);
1085     AliITSclusterV2* cllay = lay.GetCluster(index);
1086     if(cllay==0) continue;
1087     if(cllay->IsUsed()==1) continue;
1088     if(cllay->TestBit(kSAflag)==kTRUE) continue;
1089     Double_t phi   = fTable->GetPhiCluster(layer,index);
1090     Double_t lambda= fTable->GetLambdaCluster(layer,index);
1091  
1092     if(TMath::Abs(fLambdac-lambda)<lambdawindow && 
1093        TMath::Abs(fPhiEstimate-phi)<phiwindow){
1094       nc+=1;
1095       fLambdac = lambda;
1096       if(trs->GetNumberOfClustersSA()==15){
1097         delete[] firstmod;
1098         delete listc;
1099         return 0;
1100       }
1101       trs->AddClusterSA(layer,index);
1102       cllay->SetBit(kSAflag);
1103       fPhiEstimate=phi;
1104       fPointc[0]=fTable->GetXCluster(layer,index);
1105       fPointc[1]=fTable->GetYCluster(layer,index);
1106     }
1107
1108   }
1109   delete listc;
1110   delete [] firstmod;
1111   return nc;
1112
1113 }
1114 */
1115
1116
1117 //_______________________________________________________
1118 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag){
1119   //function used to to find the clusters associated to the track
1120   Int_t nc=0;
1121   AliITSlayer &lay = fgLayers[layer];
1122   Double_t r=lay.GetR(),tgl=TMath::Tan(fLambdac);
1123
1124   if(pflag==1){      
1125     Float_t cx1,cx2,cy1,cy2;
1126     FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1127     if (FindIntersection(fCoef1,fCoef2,fCoef3,-r*r,cx1,cy1,cx2,cy2)==0)
1128        return 0;
1129     Double_t fi1=TMath::ATan2(cy1,cx1);
1130     Double_t fi2=TMath::ATan2(cy2,cx2);
1131     fPhiEstimate=ChoosePoint(fi1,fi2,fPhic);
1132   }
1133
1134   Double_t dz=r*lambdawindow*TMath::Sqrt(1+tgl*tgl) + 0.3*TMath::Abs(tgl);
1135   Double_t zmax=r*tgl + zvertex + dz;
1136   Double_t zmin=r*tgl + zvertex - dz;
1137
1138   Int_t indmin = lay.FindClusterIndex(zmin);
1139   Int_t indmax = lay.FindClusterIndex(zmax);
1140   for (Int_t index=indmin; index<indmax; index++) {
1141      AliITSclusterV2 *c=lay.GetCluster(index);
1142      if (!c) continue;
1143      if (c->IsUsed()) continue;
1144      if (c->GetQ()<=0) continue;
1145      if (c->TestBit(kSAflag)==kTRUE) continue;
1146
1147      Double_t phi   =fTable->GetPhiCluster(layer,index);
1148      if (TMath::Abs(phi-fPhiEstimate)>phiwindow) continue;
1149
1150      Double_t lambda=fTable->GetLambdaCluster(layer,index);
1151      if (TMath::Abs(lambda-fLambdac)>lambdawindow) continue;
1152
1153      if(trs->GetNumberOfClustersSA()==15) return 0;
1154
1155      trs->AddClusterSA(layer,index);
1156      nc++;
1157      fLambdac=lambda;
1158      fPhiEstimate=phi;
1159      fPointc[0]=fTable->GetXCluster(layer,index);
1160      fPointc[1]=fTable->GetYCluster(layer,index);
1161
1162      c->SetBit(kSAflag);
1163   }
1164   return nc;
1165 }
1166
1167
1168 //________________________________________________________________
1169 void AliITStrackerSA::UpdatePoints(){
1170   //update of points for the estimation of the curvature  
1171
1172   //fPoint1[0]=fPoint2[0]; 
1173   //fPoint1[1]=fPoint2[1];
1174   fPoint2[0]=fPoint3[0];
1175   fPoint2[1]=fPoint3[1];
1176   fPoint3[0]=fPointc[0];
1177   fPoint3[1]=fPointc[1];
1178
1179   
1180 }
1181
1182 //___________________________________________________________________
1183 Int_t AliITStrackerSA::FindEquation(Float_t x1, Float_t y1, Float_t x2, Float_t y2, Float_t x3, Float_t y3,Float_t& a, Float_t& b, Float_t& c){
1184
1185    //given (x,y) of three recpoints (in global coordinates) 
1186    //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1187
1188    Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1189    if(den==0) return 0;
1190    a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1191    b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1192    c = -x1*x1-y1*y1-a*x1-b*y1;
1193    return 1;
1194  }
1195 //__________________________________________________________________________
1196  Int_t AliITStrackerSA::FindIntersection(Float_t a1, Float_t b1, Float_t c1, Float_t c2,Float_t& x1,Float_t& y1, Float_t& x2, Float_t& y2){
1197  
1198  //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1199  //c2 is -rlayer*rlayer
1200
1201   if(a1==0) return 0;
1202  Double_t m = c2-c1; 
1203  Double_t aA = (b1*b1)/(a1*a1)+1;
1204  Double_t bB = (-2*m*b1/(a1*a1));
1205  Double_t cC = c2+(m*m)/(a1*a1);
1206  Double_t dD = bB*bB-4*aA*cC;
1207  if(dD<0) return 0;
1208  
1209  y1 = (-bB+TMath::Sqrt(dD))/(2*aA); 
1210  y2 = (-bB-TMath::Sqrt(dD))/(2*aA); 
1211  x1 = (c2-c1-b1*y1)/a1;
1212  x2 = (c2-c1-b1*y2)/a1;
1213
1214  return 1; 
1215 }
1216 //____________________________________________________________________
1217 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t 
1218 x2,Double_t y2,Double_t x3,Double_t y3){
1219
1220   //calculates the curvature of track  
1221   Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1222   if(den==0) return 0;
1223   Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1224   Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1225   Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1226   Double_t xc=-a/2.;
1227
1228   if((a*a+b*b-4*c)<0) return 0;
1229   Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1230   if(rad==0) return 0;
1231   
1232   if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1233   if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1234   //  if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1235   // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1236   
1237   return 1/rad;
1238  
1239 }
1240 //____________________________________________________________________
1241 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1242
1243   //Returns the point closest to pp
1244
1245   Double_t diff1 = p1-pp;
1246   Double_t diff2 = p2-pp;
1247   
1248   if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1249   else fPhiEstimate=p2;  
1250   return fPhiEstimate;
1251   
1252 }
1253
1254
1255 //_________________________________________________________________
1256 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1257   // returns track with lowes chi square  
1258   if(dim==1){
1259     AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1260     return trk;
1261   }
1262   if(dim==0) return 0;
1263   Double_t * chi2 = new Double_t[dim];
1264   Int_t * index = new Int_t[dim];
1265   for(Int_t i=0;i<dim;i++){
1266     AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1267     chi2[i]=trk->GetChi2();
1268     index[i]=i;
1269   }
1270
1271   Int_t w=0;Double_t value;
1272   Int_t lp;
1273   while(w<dim){
1274     for(Int_t j=w+1;j<dim;j++){
1275       if(chi2[w]<chi2[j]){
1276         value=chi2[w];
1277         chi2[w]=chi2[j];
1278         chi2[j]=value;
1279         lp=index[w];
1280         index[w]=index[j];
1281         index[j]=lp;
1282       }
1283     }
1284     w++;
1285   }
1286
1287   AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1288   delete [] chi2;
1289   delete [] index;
1290   return trk;
1291   
1292 }
1293
1294 //__________________________________________________________
1295 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1296
1297   //function used to determine the track label
1298   
1299   Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1300   Int_t aa[6]={1,1,1,1,1,1};
1301   Int_t ff=0; 
1302   Int_t ll=0;
1303   Int_t k=0;Int_t w=0;Int_t num=6;
1304   if(lb[5]==-1) num=5;
1305   
1306   while(k<num){
1307   
1308     for(Int_t i=k+1;i<num;i++){
1309     
1310       if(lb[k]==lb[i] && aa[k]!=0){
1311       
1312         aa[k]+=1;
1313         aa[i]=0;
1314       }
1315     }
1316   k++;
1317   }
1318
1319   while(w<num){
1320   
1321     for(Int_t j=0;j<6;j++){
1322       if(aa[w]<aa[j]){
1323       ff=aa[w];
1324       aa[w]=aa[j];
1325       aa[j]=ff;
1326       ll=lb[w];
1327       lb[w]=lb[j];
1328       lb[j]=ll;
1329      }
1330     }
1331   w++;
1332   }
1333   if(num==6)  return lb[5];
1334   else return lb[4];
1335 }
1336
1337 //_____________________________________________________________________________
1338 Int_t AliITStrackerSA::Label(Int_t gl1, Int_t gl2, Int_t gl3, Int_t gl4, Int_t gl5, Int_t gl6,Int_t gl7, Int_t gl8, Int_t gl9, Int_t gl10,Int_t gl11,
1339 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1340
1341  
1342   //function used to assign label to the found track. If track is fake, the label is negative
1343
1344   Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1345   Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1346   Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1347   Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1348   Int_t lflag=0;Int_t num=6;
1349   if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1350
1351   for(Int_t i=0;i<num;i++){
1352     if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1353   }
1354
1355   if(lflag>=numberofpoints) return ll;
1356   else return -ll;
1357
1358   
1359 }
1360
1361 //_____________________________________________________________________________
1362 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1363   // Set sizes of the phi and lambda windows used for track finding
1364   fNloop = n;
1365   if(phi){ // user defined values
1366     fPhiWin = new Double_t[fNloop];
1367     fLambdaWin = new Double_t[fNloop];
1368     for(Int_t k=0;k<fNloop;k++){
1369       fPhiWin[k]=phi[k];
1370       fLambdaWin[k]=lam[k];
1371     }
1372   }
1373   else {  // default values
1374             
1375     Double_t phid[33]   = {0.002,0.003,0.004,0.0045,0.0047,
1376                            0.005,0.0053,0.0055,
1377                            0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1378                            0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1379                            0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,0.0140,0.0145};
1380     Double_t lambdad[33] = {0.003,0.004,0.005,0.005,0.005,
1381                             0.005,0.005,0.006,
1382                             0.006,0.006,0.006,0.007,0.007,0.007,0.007,
1383                             0.007,0.007,0.007,0.007,0.007,0.007,0.007,
1384                             0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008};
1385     
1386     if(fNloop!=33){
1387       fNloop = 33;
1388     }
1389     
1390     
1391     fPhiWin = new Double_t[fNloop];
1392     fLambdaWin = new Double_t[fNloop];
1393    
1394     for(Int_t k=0;k<fNloop;k++){
1395       fPhiWin[k]=phid[k];
1396       fLambdaWin[k]=lambdad[k];
1397     }
1398   
1399   }
1400
1401 }
1402