]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackerSA.cxx
SA tracking adapted to ESD
[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 "AliITSVertex.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 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertex *vert):AliITStrackerV2(geom) 
50 {
51   // Standard constructor (Vertex is known and passed to this obj.)
52   Init();
53   fVert = vert;
54   fGeom = geom;
55 }
56
57 //______________________________________________________________________
58 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) : 
59                     AliITStrackerV2(trkr) {
60   // Copy constructor
61   // Copies are not allowed. The method is protected to avoid misuse.
62   Error("AliITStrackerSA","Copy constructor not allowed\n");
63 }
64
65 //______________________________________________________________________
66 AliITStrackerSA& AliITStrackerSA::operator=(const 
67                     AliITStrackerSA& /* trkr */){
68   // Assignment operator
69   // Assignment is not allowed. The method is protected to avoid misuse.
70   Error("= operator","Assignment operator not allowed\n");
71   return *this;
72 }
73
74 //____________________________________________________________________________
75 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom) 
76 {
77   // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
78   Init();
79   fVertexer = vertexer;
80   fGeom = geom;
81 }
82
83 //____________________________________________________________________________
84 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){
85   // Copy constructor
86   fPhiEstimate = tracker.fPhiEstimate;
87   for(Int_t i=0;i<2;i++){
88     fPoint1[i]=tracker.fPoint1[i];
89     fPoint2[i]=tracker.fPoint2[i];
90     fPoint3[i]=tracker.fPoint3[i];
91     fPointc[i]=tracker.fPointc[i];
92   }
93   fLambdac = tracker.fLambdac;
94   fPhic = tracker.fPhic;
95   fCoef1 = tracker.fCoef1;
96   fCoef2 = tracker.fCoef2;
97   fCoef3 = tracker.fCoef3;
98   fNloop = tracker.fNloop;
99   fPhiWin = tracker.fPhiWin;
100   fLambdaWin = tracker.fLambdaWin;
101   if(tracker.fVertexer && tracker.fVert){
102     fVert = new AliITSVertex(*tracker.fVert);
103   }
104   else {
105     fVert = tracker.fVert;
106   }
107   fVertexer = tracker.fVertexer;
108   fGeom = tracker.fGeom;
109   fFlagLoad = tracker.fFlagLoad;
110   fTable = tracker.fTable;
111 }
112
113 //____________________________________________________________________________
114 AliITStrackerSA::~AliITStrackerSA(){
115   // destructor
116   // if fVertexer is not null, the AliITSVertex obj. is owned by this class
117   // and is deleted here
118   if(fVertexer){
119     if(fVert)delete fVert;
120   }
121   fVert = 0;
122   fVertexer = 0;
123   fGeom = 0;
124   fFlagLoad = 0;
125   if(fPhiWin)delete []fPhiWin;
126   if(fLambdaWin)delete []fLambdaWin;
127   fTable =0;
128 }
129
130 //____________________________________________________________________________
131 void AliITStrackerSA::Init(){
132   //  Reset all data members
133     fPhiEstimate=0;
134     for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
135     fLambdac=0;
136     fPhic=0;
137     fCoef1=0;
138     fCoef2=0;
139     fCoef3=0;
140     fPointc[0]=0;
141     fPointc[1]=0;
142     fVert = 0;
143     fVertexer = 0;
144     fGeom = 0;
145     fFlagLoad = 0;
146     SetWindowSizes();
147     fTable = 0;
148 }
149 //_______________________________________________________________________
150 void AliITStrackerSA::ResetForFinding(){
151   //  Reset data members used in all loops during track finding
152     fPhiEstimate=0;
153     for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
154     fLambdac=0;
155     fPhic=0;
156     fCoef1=0;
157     fCoef2=0;
158     fCoef3=0;
159     fPointc[0]=0;
160     fPointc[1]=0;
161 }
162 //____________________________________________________________________________
163 void AliITStrackerSA::FindTracks(TTree *clusterTree, TTree *out,Int_t evnumber,char *opt){
164
165   /**************************************************************************
166    * This function finds primary tracks.
167    * Options: to use only ITS execute FindTracks                            *
168    *          to define good tracks with 6 points out of 6 in the ITS write *
169    *          *opt="6/6". The default is *opt="6/6"                         *
170    *                                                                        * 
171    *                                                                        *
172    * Example: to execute function with only the ITS (no combined tracking   *
173    *          with TPC+ITS) and requiring 5/6 points to define a good track *
174    *          use: FindTracks(treein,treeout,evnumber,"5/6")                *
175    *          to execute combined tracking, before using FindTracks use     *
176    *          UseFoundTracksV2                                              *
177    *************************************************************************/
178
179   //options
180   TString choice(opt);
181   Bool_t sixpoints= choice.Contains("6/6");
182
183   //Get primary vertex
184   if(fVertexer){
185     if(fVert)delete fVert;
186     fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
187   }
188   else {
189     gAlice->GetEvent(evnumber);
190     if(!fVert){
191       Fatal("FindTracks","Vertex is missing\n");
192       return;
193     }
194   }
195   Double_t primaryVertex[3];
196   Double_t errorsprimvert[3];
197   fVert->GetXYZ(primaryVertex);
198   fVert->GetSigmaXYZ(errorsprimvert);
199   if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
200     Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
201     errorsprimvert[0]=0.0001;
202     errorsprimvert[1]=0.0001;
203   }
204   fVert->PrintStatus();
205
206
207   //Fill array with cluster indices for each module
208   if(!fTable){
209     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
210     fTable->FillArray(clusterTree,evnumber);
211     fTable->FillArrayCoorAngles(clusterTree,evnumber); 
212   }
213
214   SetEventNumber(evnumber);
215   if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree);
216    
217 //Fill tree for found tracks
218   AliITStrackV2* outrack=0;
219   TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
220   if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
221   else branch->SetAddress(&outrack);
222
223   
224   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
225   for(Int_t i=0;i<fGeom->GetNlayers();i++){
226     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
227   }
228   // firstmod [i] number of the first module in the ITS layer i.
229     
230   AliITSlayer &layer=fgLayers[0];   // first layer
231   Int_t ntrack=0;
232   Int_t dim=layer.GetNumberOfClusters();
233   //loop on the different windows
234   for(Int_t nloop=0;nloop<fNloop;nloop++){
235     for(Int_t ncl=0;ncl<dim;ncl++){ 
236       //loop starting from layer 0
237       ResetForFinding();
238       Int_t pflag=0;
239       AliITSclusterV2* cl = layer.GetCluster(ncl);
240       if(cl->IsUsed()==1) continue;     
241       fPhic = fTable->GetPhiCluster(0,ncl);
242       fLambdac = fTable->GetLambdaCluster(0,ncl);
243       fPhiEstimate = fPhic;
244       AliITStrackSA* trs = new AliITStrackSA();      
245       fPoint1[0]=primaryVertex[0];
246       fPoint1[1]=primaryVertex[1];
247       fPoint2[0]=fTable->GetXCluster(0,ncl);;
248       fPoint2[1]=fTable->GetYCluster(0,ncl);;
249       
250       Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
251       for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
252       nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
253       nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
254
255       if(nn[1]>0){
256         pflag=1;
257         fPoint3[0] = fPointc[0];
258         fPoint3[1] = fPointc[1];
259       }
260       nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
261       if(nn[1]==0 && nn[2]==0) pflag=0;
262       if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
263       if(nn[2]!=0 && nn[1]==0){
264         pflag=1;
265         fPoint3[0]=fPointc[0];
266         fPoint3[1]=fPointc[1];
267       }
268
269       nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
270       pflag=1;
271       if(nn[3]!=0) UpdatePoints();
272       nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable); 
273       pflag=1;
274       if(nn[4]!=0) UpdatePoints();
275       nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable); 
276           
277
278       Int_t layOK=0;
279       Int_t numberofpoints;
280       if(sixpoints) numberofpoints=6;  //check of the candidate track
281       else numberofpoints=5;           //if track is good (with the required number        
282       for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){    //of points) it is written on file
283         if(nn[nnp]!=0) layOK+=1;
284       }
285       if(layOK>=numberofpoints){
286         AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
287         if(tr2==0){
288           Int_t nct = trs->GetNumberOfClustersSA();
289           while(nct--){
290             Int_t index = trs->GetClusterIndexSA(nct);
291             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
292             if(kl->IsUsed()==1) kl->Use();
293           }
294           continue;
295         }
296         outrack=tr2;  
297         out->Fill();
298         ntrack++;
299         Int_t nct = tr2->GetNumberOfClusters();
300
301         while(nct--){
302           Int_t index = tr2->GetClusterIndex(nct);     
303           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);      
304           if(kl->IsUsed()==0) kl->Use();
305
306         }
307       } 
308       else{
309         Int_t nct = trs->GetNumberOfClustersSA();
310         while(nct--){
311           Int_t index = trs->GetClusterIndexSA(nct);
312           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
313           if(kl->IsUsed()==1) kl->Use();
314         }
315       }
316       delete trs;
317
318     }//end loop on clusters of layer1
319
320   }//end loop2
321
322   //if 5/6 points are required, second loop starting 
323   //from second layer, to find tracks with point of 
324   //layer 1 missing
325    
326   if(!sixpoints){
327     //   counter for clusters on each layer  
328     Int_t * nn = new Int_t[fGeom->GetNlayers()-1];      
329     for(Int_t nloop=0;nloop<fNloop;nloop++){
330       AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
331       Int_t ncl2=layer2.GetNumberOfClusters();
332       while(ncl2--){ //loop starting from layer 2
333         ResetForFinding();
334         Int_t pflag=0;
335         AliITSclusterV2* cl = layer2.GetCluster(ncl2);
336         if(cl->IsUsed()==1) continue;
337         fPhic = fTable->GetPhiCluster(1,ncl2);
338         fLambdac = fTable->GetLambdaCluster(1,ncl2);
339         fPhiEstimate = fPhic;
340         AliITStrackSA* trs = new AliITStrackSA(); 
341      
342         fPoint1[0]=primaryVertex[0];
343         fPoint1[1]=primaryVertex[1];
344         fPoint2[0]=fTable->GetXCluster(1,ncl2);;
345         fPoint2[1]=fTable->GetYCluster(1,ncl2);;
346    
347         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
348         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
349           nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
350                    trs,primaryVertex[2],pflag,fTable);
351           if(nn[kk]==0)break;
352           if(kk>0){
353             UpdatePoints();
354             pflag = 1;
355           }
356         }
357         Int_t fl=0;
358         for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
359           if(nn[nnp]!=0) fl+=1;
360         }
361         if(fl>=5){  // 5/6       
362           AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
363           if(tr2==0){
364             Int_t nct = trs->GetNumberOfClustersSA();
365             while(nct--){
366               Int_t index = trs->GetClusterIndexSA(nct);
367               AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
368               if(kl->IsUsed()==1) kl->Use();
369             }
370             continue;
371           }
372           outrack=tr2;
373           out->Fill();
374           Int_t nct = tr2->GetNumberOfClusters();
375           while(nct--){
376             Int_t index = tr2->GetClusterIndex(nct);     
377             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
378             if(kl==0) continue;
379             if(kl->IsUsed()==0) kl->Use();
380           }
381         }       
382         else{
383           Int_t nct = trs->GetNumberOfClustersSA();
384           while(nct--){
385             Int_t index = trs->GetClusterIndexSA(nct);
386             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
387             if(kl==0) continue;
388             if(kl->IsUsed()==1) kl->Use();
389           }
390         }
391         delete trs;
392       }//end loop on clusters of layer2
393     }
394     delete [] nn;
395   }  //end opt="5/6"  
396
397   delete [] firstmod;
398
399   UnloadClusters();
400 }
401
402 //______________________________________________________________________
403 void AliITStrackerSA::FindTracks(TTree *clusterTree, AliESD* event, Int_t evnumber,char *opt){
404   /**************************************************************************
405    * Options: to use only ITS execute FindTracks                            *
406    *          to define good tracks with 6 points out of 6 in the ITS write *
407    *          *opt="6/6". The default is *opt="6/6"                         *
408    *                                                                        * 
409    *                                                                        *
410    * Example: to execute function with only the ITS (no combined tracking   *
411    *          with TPC+ITS) and requiring 5/6 points to define a good track *
412    *          use: FindTracks(treein,treeout,evnumber,"5/6")                *
413    *          to execute combined tracking, before using FindTracks use     *
414    *          UseFoundTracksV2                                              *
415    *************************************************************************/
416
417   //options
418   TString choice(opt);
419   Bool_t sixpoints= choice.Contains("6/6");
420   //Get primary vertex
421   if(fVertexer){
422     if(fVert)delete fVert;
423     fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
424   }
425   else {
426     gAlice->GetEvent(evnumber);
427     if(!fVert){
428       Fatal("FindTracks","Vertex is missing\n");
429       return;
430     }
431   }
432   Double_t primaryVertex[3];
433   Double_t errorsprimvert[3];
434   fVert->GetXYZ(primaryVertex);
435   fVert->GetSigmaXYZ(errorsprimvert);
436   if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
437     Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
438     errorsprimvert[0]=0.0001;
439     errorsprimvert[1]=0.0001;
440   }
441   fVert->PrintStatus();
442
443   //Fill array with cluster indices for each module
444   if(!fTable){
445     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
446     fTable->FillArray(clusterTree,evnumber);   
447     fTable->FillArrayCoorAngles(clusterTree,evnumber); 
448  }
449   
450   SetEventNumber(evnumber);
451   if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree);
452
453   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
454   for(Int_t i=0;i<fGeom->GetNlayers();i++){
455     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
456   }
457   // firstmod [i] number of the first module in the ITS layer i.
458   
459     
460   AliITSlayer &layer=fgLayers[0];   // first layer
461   Int_t ntrack=0;
462
463   Int_t dim=layer.GetNumberOfClusters();
464   //loop on the different windows
465   for(Int_t nloop=0;nloop<fNloop;nloop++){
466     for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
467       ResetForFinding();
468       Int_t pflag=0;
469       AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
470       if(cl==0) continue;
471       if(cl->IsUsed()==1) continue;
472       fPhic = fTable->GetPhiCluster(0,ncl);
473       fLambdac = fTable->GetLambdaCluster(0,ncl);
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,fTable);
483           
484       nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
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,fTable);
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,fTable);
500       pflag=1;
501       if(nn[3]!=0) UpdatePoints();
502       nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable); 
503       pflag=1;
504       if(nn[4]!=0) UpdatePoints();
505       nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable); 
506           
507
508       Int_t layOK=0;
509       Int_t numberofpoints;
510       if(sixpoints) 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,opt);
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->IsUsed()==1) kl->Use();
523           }
524           continue;
525         }
526         
527         AliESDtrack outtrack;
528         outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
529         event->AddTrack(&outtrack);
530         
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           if(kl->IsUsed()==0) kl->Use();
537         }
538       } 
539       else{
540         Int_t nct = trs->GetNumberOfClustersSA();
541         while(nct--){
542           Int_t index = trs->GetClusterIndexSA(nct);
543           AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
544           if(kl->IsUsed()==1) kl->Use();
545         }
546       }
547       delete trs;
548
549     }//end loop on clusters of layer1
550
551   }//end loop2
552
553   //if 5/6 points are required, second loop starting 
554   //from second layer, to find tracks with point of 
555   //layer 1 missing
556    
557   if(!sixpoints){
558     //   counter for clusters on each layer  
559     Int_t * nn = new Int_t[fGeom->GetNlayers()-1];      
560     for(Int_t nloop=0;nloop<fNloop;nloop++){
561       AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
562       Int_t ncl2=layer2.GetNumberOfClusters();
563       while(ncl2--){ //loop starting from layer 2
564         ResetForFinding();
565         Int_t pflag=0;
566         AliITSclusterV2* cl = layer2.GetCluster(ncl2);
567         if(cl->IsUsed()==1) continue;
568         fPhic = fTable->GetPhiCluster(1,ncl2);
569         fLambdac = fTable->GetLambdaCluster(1,ncl2);
570         fPhiEstimate = fPhic;
571         AliITStrackSA* trs = new AliITStrackSA(); 
572         fPoint1[0]=primaryVertex[0];
573         fPoint1[1]=primaryVertex[1];
574         fPoint2[0]=fTable->GetXCluster(1,ncl2);
575         fPoint2[1]=fTable->GetYCluster(1,ncl2);
576    
577         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
578         for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
579           nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
580                    trs,primaryVertex[2],pflag,fTable);
581           if(nn[kk]==0)break;
582           if(kk>0){
583             UpdatePoints();
584             pflag = 1;
585           }
586         }
587         Int_t fl=0;
588         for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
589           if(nn[nnp]!=0) fl+=1;
590         }
591         if(fl>=5){  // 5/6       
592           AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
593           if(tr2==0){
594             ntrack++;
595             Int_t nct = trs->GetNumberOfClustersSA();
596             while(nct--){
597               Int_t index = trs->GetClusterIndexSA(nct);
598               AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
599               if(kl->IsUsed()==1) kl->Use();
600             }
601             continue;
602           }
603
604           AliESDtrack outtrack;
605           outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
606           event->AddTrack(&outtrack);
607
608           Int_t nct = tr2->GetNumberOfClusters();
609           while(nct--){
610             Int_t index = tr2->GetClusterIndex(nct);     
611             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
612             if(kl==0) continue;
613             if(kl->IsUsed()==0) kl->Use();
614           }
615         }       
616         else{
617           Int_t nct = trs->GetNumberOfClustersSA();
618           while(nct--){
619             Int_t index = trs->GetClusterIndexSA(nct);
620             AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
621             if(kl==0) continue;
622             if(kl->IsUsed()==1) kl->Use();
623           }
624         }
625         delete trs;
626       }//end loop on clusters of layer2
627     }
628     delete [] nn;
629   }  //end opt="5/6"  
630
631   delete [] firstmod;
632   delete fTable;   
633   UnloadClusters();
634   Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
635
636
637
638 }
639 //________________________________________________________________________
640
641 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert,char *opt){
642   //fit of the found track
643
644  
645 //options
646   TString choice(opt);
647   Bool_t sixpoints= choice.Contains("6/6");
648   
649   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
650   for(Int_t i=0;i<fGeom->GetNlayers();i++){
651     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
652   }  
653   TObjArray* listoftracks=new TObjArray(0,0);
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->IsUsed()==1) cl->Use();
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)) listoftracks->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=listoftracks->GetEntries();
843   if(dim==0){
844     delete listoftracks;
845     for(Int_t i=0;i<fGeom->GetNlayers();i++){
846       delete listlayer[i];
847     }
848     delete listlayer;
849     return 0;
850   }
851
852   AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(listoftracks,dim);
853
854   if(otrack==0) return 0;
855   Int_t * indexc = new Int_t[fGeom->GetNlayers()];
856   for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
857   for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
858     indexc[nind] = otrack->GetClusterIndex(nind);
859   }      
860   AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
861   AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);     
862   AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);     
863   AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
864   AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
865   Int_t labl[3]={-1,-1,-1};
866   if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
867     AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
868     labl[0]=cl5->GetLabel(0);
869     labl[1]=cl5->GetLabel(1);
870     labl[2]=cl5->GetLabel(2);
871   }
872   delete [] indexc;
873   if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
874     labl[0]=-1;
875     labl[1]=-1;
876     labl[2]=-1;
877   }
878   Int_t numberofpoints;
879   if(sixpoints) numberofpoints=6;
880   else numberofpoints=5;
881   Int_t label =  Label(cl0->GetLabel(0),cl1->GetLabel(0), 
882                        cl2->GetLabel(0),cl3->GetLabel(0),
883                        cl4->GetLabel(0),labl[0],
884                        cl0->GetLabel(1),cl1->GetLabel(1),
885                        cl2->GetLabel(1),cl3->GetLabel(1),
886                        cl4->GetLabel(1),labl[1],
887                        cl0->GetLabel(2),cl1->GetLabel(2),
888                        cl2->GetLabel(2),cl3->GetLabel(2),
889                        cl4->GetLabel(2),labl[2],numberofpoints);
890   
891   otrack->SetLabel(label);  
892   delete listoftracks;
893   for(Int_t i=0;i<fGeom->GetNlayers();i++){
894     delete listlayer[i];
895   }
896   delete listlayer; 
897   delete [] firstmod;
898   return otrack;
899
900 }
901
902 //_______________________________________________________________________
903 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2, TTree* clustertree){
904   // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
905   //(or AliITStrackV2 tracks found with function FindTracks of this class)
906
907   //Get primary vertex
908   if(fVertexer){
909     if(fVert)delete fVert;
910     fVert = fVertexer->FindVertexForCurrentEvent(evnum);
911   }
912   else {
913     gAlice->GetEvent(evnum);
914     if(!fVert){
915       Fatal("FindTracks","Vertex is missing\n");
916       return;
917     }
918   }
919   Double_t primaryVertex[3];
920   fVert->GetXYZ(primaryVertex);
921
922   if(!fTable){
923     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
924     fTable->FillArray(clustertree,evnum);
925     fTable->FillArrayCoorAngles(clustertree,evnum);
926   }
927   SetEventNumber(evnum);
928   if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
929   SetFlagLoadedClusters(1);
930
931   TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
932   if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
933   AliITStrackV2* ttrrt = new AliITStrackV2;
934   bra->SetAddress(&ttrrt);
935
936   for(Int_t nj=0;nj<treev2->GetEntries();nj++){
937     treev2->GetEvent(nj);
938     Int_t ncl = ttrrt->GetNumberOfClusters();
939     for(Int_t k=0;k<ncl;k++){
940       Int_t index = ttrrt->GetClusterIndex(k);
941       AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
942       if(clui->IsUsed()==0) clui->Use();  
943     }
944   }
945   delete ttrrt;
946   
947 }
948
949 //_______________________________________________________________________
950 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,AliESD *event, TTree* 
951 clustertree){
952   // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
953
954   //Get primary vertex
955   if(fVertexer){
956     if(fVert)delete fVert;
957     fVert = fVertexer->FindVertexForCurrentEvent(evnum);
958   }
959   else {
960     gAlice->GetEvent(evnum);
961     if(!fVert){
962       Fatal("FindTracks","Vertex is missing\n");
963       return;
964     }
965   }
966   Double_t primaryVertex[3];
967   fVert->GetXYZ(primaryVertex);
968
969   if(!fTable){
970     fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
971     fTable->FillArray(clustertree,evnum);
972     fTable->FillArrayCoorAngles(clustertree,evnum); 
973   }
974   SetEventNumber(evnum);
975   if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
976   SetFlagLoadedClusters(1);
977   Int_t ntracks = event->GetNumberOfTracks();
978   while (ntracks--) {
979     AliESDtrack *esd=event->GetTrack(ntracks);
980     if ((esd->GetStatus()&
981          AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue; 
982     UInt_t idx[6];
983     Int_t ncl = esd->GetITSclusters(idx);
984     for(Int_t  clu=0; clu<ncl; clu++){
985       AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
986       if(cl->IsUsed()==0) cl->Use();
987     }
988   }
989   
990   Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
991   
992
993 }
994
995 //_______________________________________________________
996 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
997   //function used to to find the clusters associated to the track
998   Int_t nc=0;
999   AliITSlayer &lay = fgLayers[layer];
1000   Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
1001   for(Int_t i=0;i<fGeom->GetNlayers();i++){
1002     firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
1003   }
1004   if(pflag==1){
1005       
1006     Float_t cx1,cx2,cy1,cy2;
1007     FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1008     Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
1009     if(fun==0) return 0;
1010     Double_t fi1 =TMath::ATan2(cy1,cx1);
1011     Double_t fi2 =TMath::ATan2(cy2,cx2);
1012     fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1013   }
1014
1015   Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1016   Double_t zed1  =  TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1017   Double_t zed2  =  TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1018   
1019   Double_t fi = fPhiEstimate;
1020   Int_t nmod  = lay.FindDetectorIndex(fi,zed)+firstmod[layer];
1021  
1022   Int_t nm[8]={0,0,0,0,0,0,0,0};
1023   nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed)+firstmod[layer];
1024   nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed)+firstmod[layer];
1025   nm[2]  = lay.FindDetectorIndex(fi,zed1)+firstmod[layer];
1026   nm[3]  = lay.FindDetectorIndex(fi,zed2)+firstmod[layer];
1027   nm[4]  = lay.FindDetectorIndex(fi+phiwindow,zed1)+firstmod[layer];
1028   nm[5]  = lay.FindDetectorIndex(fi-phiwindow,zed1)+firstmod[layer];
1029   nm[6]  = lay.FindDetectorIndex(fi+phiwindow,zed2)+firstmod[layer];
1030   nm[7]  = lay.FindDetectorIndex(fi-phiwindow,zed2)+firstmod[layer];
1031
1032  
1033    if(nmod<0) return 0;
1034   Int_t nn=0;
1035   TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1036   TArrayI* list = new TArrayI(array->GetSize());
1037   for(Int_t i=0;i<array->GetSize();i++){
1038     Int_t in=(Int_t)array->At(i);
1039     list->AddAt(in,nn);
1040     nn++;
1041   }
1042   
1043   for(Int_t ii=0;ii<8;ii++){
1044     if(nm[ii]!=nmod && nm[ii]>=0){
1045       TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]);
1046       list->Set(list->GetSize()+ar->GetSize());
1047       for(Int_t j=0;j<ar->GetSize();j++){
1048         Int_t in=(Int_t)ar->At(j);
1049         list->AddAt(in,nn);
1050         nn++;
1051       }
1052     }
1053   }
1054   for(Int_t i=0;i<list->GetSize();i++){
1055     Int_t index = (Int_t)list->At(i);
1056     AliITSclusterV2* cllay = lay.GetCluster(index);
1057     if(cllay==0) continue;
1058     if(cllay->IsUsed()==1) continue;  
1059     Double_t phi   = fTable->GetPhiCluster(layer,index);
1060     Double_t lambda= fTable->GetLambdaCluster(layer,index);
1061  
1062     if(TMath::Abs(fLambdac-lambda)<lambdawindow && 
1063        TMath::Abs(fPhiEstimate-phi)<phiwindow){
1064       nc+=1;
1065       fLambdac = lambda;
1066       if(trs->GetNumberOfClustersSA()==20){
1067         delete list;
1068         return 0;
1069       }
1070       trs->AddClusterSA(layer,index);
1071       cllay->Use();
1072       fPhiEstimate=phi;
1073       fPointc[0]=fTable->GetXCluster(layer,index);
1074       fPointc[1]=fTable->GetYCluster(layer,index);
1075     }
1076
1077   }
1078   delete list;
1079   delete [] firstmod;
1080   return nc;
1081
1082 }
1083
1084
1085
1086 //________________________________________________________________
1087 void AliITStrackerSA::UpdatePoints(){
1088   //update of points for the estimation of the curvature  
1089
1090   fPoint1[0]=fPoint2[0]; 
1091   fPoint1[1]=fPoint2[1];
1092   fPoint2[0]=fPoint3[0];
1093   fPoint2[1]=fPoint3[1];
1094   fPoint3[0]=fPointc[0];
1095   fPoint3[1]=fPointc[1];
1096
1097   
1098 }
1099
1100
1101
1102 //___________________________________________________________________
1103 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){
1104
1105    //given (x,y) of three recpoints (in global coordinates) 
1106    //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1107
1108    Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1109    if(den==0) return 0;
1110    a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1111    b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1112    c = -x1*x1-y1*y1-a*x1-b*y1;
1113    return 1;
1114  }
1115 //__________________________________________________________________________
1116  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){
1117  
1118  //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1119  //c2 is -rlayer*rlayer
1120
1121   if(a1==0) return 0;
1122  Float_t m = c2-c1; 
1123  Float_t aA = (b1*b1)/(a1*a1)+1;
1124  Float_t bB = (-2*m*b1/(a1*a1));
1125  Float_t cC = c2+(m*m)/(a1*a1);
1126  if((bB*bB-4*aA*cC)<0) return 0;
1127  
1128  y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA); 
1129  y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA); 
1130  x1 = (c2-c1-b1*y1)/a1;
1131  x2 = (c2-c1-b1*y2)/a1;
1132
1133  return 1; 
1134 }
1135 //____________________________________________________________________
1136 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t 
1137 x2,Double_t y2,Double_t x3,Double_t y3){
1138
1139   //calculates the curvature of track  
1140   Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1141   if(den==0) return 0;
1142   Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1143   Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1144   Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1145   Double_t xc=-a/2.;
1146
1147   if((a*a+b*b-4*c)<0) return 0;
1148   Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1149   if(rad==0) return 0;
1150   
1151   if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1152   if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1153   //  if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1154   // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1155   
1156   return 1/rad;
1157  
1158 }
1159 //____________________________________________________________________
1160 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1161
1162   //Returns the point closest to pp
1163
1164   Double_t diff1 = p1-pp;
1165   Double_t diff2 = p2-pp;
1166   
1167   if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1168   else fPhiEstimate=p2;  
1169   return fPhiEstimate;
1170   
1171 }
1172
1173
1174 //_________________________________________________________________
1175 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1176   // returns track with lowes chi square  
1177   if(dim==1){
1178     AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1179     return trk;
1180   }
1181   if(dim==0) return 0;
1182   Double_t * chi2 = new Double_t[dim];
1183   Int_t * index = new Int_t[dim];
1184   for(Int_t i=0;i<dim;i++){
1185     AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1186     chi2[i]=trk->GetChi2();
1187     index[i]=i;
1188   }
1189
1190   Int_t w=0;Double_t value;
1191   Int_t lp;
1192   while(w<dim){
1193     for(Int_t j=w+1;j<dim;j++){
1194       if(chi2[w]<chi2[j]){
1195         value=chi2[w];
1196         chi2[w]=chi2[j];
1197         chi2[j]=value;
1198         lp=index[w];
1199         index[w]=index[j];
1200         index[j]=lp;
1201       }
1202     }
1203     w++;
1204   }
1205
1206   AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1207   delete [] chi2;
1208   delete [] index;
1209   return trk;
1210   
1211 }
1212
1213 //__________________________________________________________
1214 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1215
1216   //function used to determine the track label
1217   
1218   Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1219   Int_t aa[6]={1,1,1,1,1,1};
1220   Int_t ff=0; 
1221   Int_t ll=0;
1222   Int_t k=0;Int_t w=0;Int_t num=6;
1223   if(lb[5]==-1) num=5;
1224   
1225   while(k<num){
1226   
1227     for(Int_t i=k+1;i<num;i++){
1228     
1229       if(lb[k]==lb[i] && aa[k]!=0){
1230       
1231         aa[k]+=1;
1232         aa[i]=0;
1233       }
1234     }
1235   k++;
1236   }
1237
1238   while(w<num){
1239   
1240     for(Int_t j=0;j<6;j++){
1241       if(aa[w]<aa[j]){
1242       ff=aa[w];
1243       aa[w]=aa[j];
1244       aa[j]=ff;
1245       ll=lb[w];
1246       lb[w]=lb[j];
1247       lb[j]=ll;
1248      }
1249     }
1250   w++;
1251   }
1252   if(num==6)  return lb[5];
1253   else return lb[4];
1254 }
1255
1256 //_____________________________________________________________________________
1257 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,
1258 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1259
1260  
1261   //function used to assign label to the found track. If track is fake, the label is negative
1262
1263   Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1264   Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1265   Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1266   Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1267   Int_t lflag=0;Int_t num=6;
1268   if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1269
1270   for(Int_t i=0;i<num;i++){
1271     if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1272   }
1273
1274   if(lflag>=numberofpoints) return ll;
1275   else return -ll;
1276
1277   
1278 }
1279
1280 //_____________________________________________________________________________
1281 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1282   // Set sizes of the phi and lambda windows used for track finding
1283   fNloop = n;
1284   if(phi){ // user defined values
1285     fPhiWin = new Double_t[fNloop];
1286     fLambdaWin = new Double_t[fNloop];
1287     for(Int_t k=0;k<fNloop;k++){
1288       fPhiWin[k]=phi[k];
1289       fLambdaWin[k]=lam[k];
1290     }
1291   }
1292   else {  // default values
1293     
1294     Double_t phid[40]   = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1295                            0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1296                            0.005,0.0053,0.0055,
1297                            0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1298                            0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1299                            0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1300                            0.014};
1301     Double_t lambdad[40] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1302                             0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1303                             0.005,0.0053,0.0055,
1304                             0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1305                             0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1306                             0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1307                             0.014};
1308     
1309     if(fNloop!=40){
1310       fNloop = 40;
1311       Warning("SetWindowSizes","Number of loop forced to the default value %d",fNloop);
1312     }
1313    
1314     fPhiWin = new Double_t[fNloop];
1315     fLambdaWin = new Double_t[fNloop];
1316    
1317     for(Int_t k=0;k<fNloop;k++){
1318       fPhiWin[k]=phid[k];
1319       fLambdaWin[k]=lambdad[k];
1320     }
1321   
1322   }
1323
1324 }
1325
1326