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