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