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