]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITStrackerUpgrade.cxx
fixes for memory leaks and a bug in loading clusters
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITStrackerUpgrade.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ////////////////////////////////////////////////////
19 //  ITS Upgrade Stand alone tracker class         //
20 //  Authors: A.Mastroserio C.Terrevoli            //
21 //  e-mail:annalisa.mastroserio@cern.ch           //
22 //         cristina.terrevoli@ba.infn.it          //
23 ////////////////////////////////////////////////////
24
25 #include <stdlib.h>
26
27 #include <TArrayI.h>
28 #include <TBranch.h>
29 #include <TObjArray.h>
30 #include <TTree.h>
31 #include "AliRunLoader.h"
32 #include "AliITSLoader.h"
33 #include "AliESDEvent.h"
34 #include "AliESDVertex.h"
35 #include "AliESDtrack.h"
36 #include "AliITSVertexer.h"
37 #include "AliITSclusterTable.h"
38 #include "AliITSRecPoint.h"
39 #include "AliITStrackSA.h"
40 #include "AliITStrackerMI.h"
41 #include "AliITSlayerUpgrade.h"
42 #include "AliITSsegmentationUpgrade.h"
43 #include "AliITStrackerUpgrade.h"
44 #include "AliITSReconstructor.h"
45 #include "AliLog.h"
46 #include "AliRun.h"
47
48 ClassImp(AliITStrackerUpgrade)
49
50 //____________________________________________________________________________
51   AliITStrackerUpgrade::AliITStrackerUpgrade():AliITStrackerMI(),
52                                                fNLayer(6),
53                                                fPhiEstimate(0),
54                                                fITSStandAlone(0),
55                                                fLambdac(0),
56                                                fPhic(0),
57                                                fCoef1(0),
58                                                fCoef2(0),
59                                                fCoef3(0),
60                                                fNloop(0),
61                                                fPhiWin(0),
62                                                fLambdaWin(0),
63                                                fVert(0),
64                                                fVertexer(0),
65                                                fListOfTracks(0),
66                                                fListOfSATracks(0),
67                                                fITSclusters(0),
68                                                fInwardFlag(0),
69                                                fOuterStartLayer(0),
70                                                fInnerStartLayer(5),
71                                                fMinNPoints(0),
72                                                fMinQ(0.),
73                                                fLayers(0),
74                                                fCluLayer(0),
75                                                fCluCoord(0){
76     // Default constructor
77     Init();
78  
79   }
80 //________________________________________________________________________
81 AliITStrackerUpgrade::AliITStrackerUpgrade(const AliITStrackerUpgrade& tracker):AliITStrackerMI(),
82                                                                                 fNLayer(tracker.fNLayer),
83                                                                                 fPhiEstimate(tracker.fPhiEstimate),
84                                                                                 fITSStandAlone(tracker.fITSStandAlone),
85                                                                                 fLambdac(tracker.fLambdac),
86                                                                                 fPhic(tracker.fPhic),
87                                                                                 fCoef1(tracker.fCoef1),
88                                                                                 fCoef2(tracker.fCoef2),
89                                                                                 fCoef3(tracker.fCoef3),
90                                                                                 fNloop(tracker.fNloop),
91                                                                                 fPhiWin(tracker.fPhiWin),
92                                                                                 fLambdaWin(tracker.fLambdaWin),
93                                                                                 fVert(tracker.fVert),
94                                                                                 fVertexer(tracker.fVertexer),
95                                                                                 fListOfTracks(tracker.fListOfTracks),
96                                                                                 fListOfSATracks(tracker.fListOfSATracks),
97                                                                                 fITSclusters(tracker.fITSclusters),
98                                                                                 fInwardFlag(tracker.fInwardFlag),
99                                                                                 fOuterStartLayer(tracker.fOuterStartLayer),
100                                                                                 fInnerStartLayer(tracker.fInnerStartLayer),
101                                                                                 fMinNPoints(tracker.fMinNPoints),
102                                                                                 fMinQ(tracker.fMinQ),
103                                                                                 fLayers(tracker.fLayers),
104                                                                                 fCluLayer(tracker.fCluLayer),
105                                                                                 fCluCoord(tracker.fCluCoord) {
106   // Copy constructor
107   for(Int_t i=0;i<2;i++){
108     fPoint1[i]=tracker.fPoint1[i];
109     fPoint2[i]=tracker.fPoint2[i];
110     fPoint3[i]=tracker.fPoint3[i];
111     fPointc[i]=tracker.fPointc[i];
112   }
113   if(tracker.fVertexer && tracker.fVert){
114     fVert = new AliESDVertex(*tracker.fVert);
115   }
116   else {
117     fVert = tracker.fVert;
118   }
119   for(Int_t i=0;i<6;i++){
120     fCluLayer[i] = tracker.fCluLayer[i];
121     fCluCoord[i] = tracker.fCluCoord[i];
122   }
123 }
124 //______________________________________________________________________
125 AliITStrackerUpgrade& AliITStrackerUpgrade::operator=(const AliITStrackerUpgrade& source){
126   // Assignment operator. 
127   this->~AliITStrackerUpgrade();
128   new(this) AliITStrackerUpgrade(source);
129   return *this;
130  
131 }
132
133 //____________________________________________________________________________
134 AliITStrackerUpgrade::~AliITStrackerUpgrade(){
135   // destructor
136   // if fVertexer is not null, the AliESDVertex obj. is owned by this class
137   // and is deleted here
138   if(fVertexer){
139     if(fVert)delete fVert;
140   }
141   fVert = 0;
142   fVertexer = 0;
143  
144   if(fPhiWin)delete []fPhiWin;
145   if(fLambdaWin)delete []fLambdaWin;
146   fListOfTracks->Delete();
147   delete fListOfTracks;
148   fListOfSATracks->Delete();
149   delete fListOfSATracks;
150   if(fCluLayer){
151     for(Int_t i=0;i<6;i++){
152       if(fCluLayer[i]){
153         fCluLayer[i]->Delete();
154         delete fCluLayer[i];
155       }
156     }
157     delete [] fCluLayer;
158   }
159   if(fCluCoord){
160     for(Int_t i=0;i<AliITSgeomTGeo::GetNLayers();i++){
161       if(fCluCoord[i]){
162         fCluCoord[i]->Delete();
163         delete fCluCoord[i];
164       }
165     }
166     delete [] fCluCoord;
167   }
168   
169 }
170
171 //____________________________________________________________________________
172 Int_t AliITStrackerUpgrade::Clusters2Tracks(AliESDEvent *event){
173   // This method is used to find and fit the tracks. By default the corresponding
174   // method in the parent class is invoked. In this way a combined tracking
175   // TPC+ITS is performed. If the flag fITSStandAlone is true, the tracking
176   // is done in the ITS only. In the standard reconstruction chain this option
177   // can be set via AliReconstruction::SetOption("ITS","onlyITS")
178   Int_t rc=0;
179   if(!fITSStandAlone){
180     rc=AliITStrackerMI::Clusters2Tracks(event);
181   }
182   else {
183     AliDebug(1,"Stand Alone flag set: doing tracking in ITS alone\n");
184   }
185   if(!rc){
186     rc=FindTracks(event,kFALSE);
187     if(AliITSReconstructor::GetRecoParam()->GetSAUseAllClusters()==kTRUE) {
188       rc=FindTracks(event,kTRUE);
189     }
190   }
191   return rc;
192 }
193 //_________________________________________________-
194 void AliITStrackerUpgrade::Init(){
195   //  Reset all data members
196   fPhiEstimate=0;
197   for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
198   fLambdac=0;
199   fPhic=0;
200   fCoef1=0;
201   fCoef2=0;
202   fCoef3=0;
203   fPointc[0]=0;
204   fPointc[1]=0;
205   fVert = 0;
206   fVertexer = 0;
207   Int_t nLoops=AliITSReconstructor::GetRecoParam()->GetNLoopsSA();
208   if(nLoops==33){
209     SetFixedWindowSizes();
210   }else{
211     Double_t phimin=AliITSReconstructor::GetRecoParam()->GetMinPhiSA();
212     Double_t phimax=AliITSReconstructor::GetRecoParam()->GetMaxPhiSA();
213     Double_t lambmin=AliITSReconstructor::GetRecoParam()->GetMinLambdaSA();
214     Double_t lambmax=AliITSReconstructor::GetRecoParam()->GetMaxLambdaSA();
215     SetCalculatedWindowSizes(nLoops,phimin,phimax,lambmin,lambmax);
216   }
217   fMinQ=AliITSReconstructor::GetRecoParam()->GetSAMinClusterCharge();
218   fITSclusters = 0;
219   SetOuterStartLayer(1);
220   SetSAFlag(kFALSE);
221   fListOfTracks=new TClonesArray("AliITStrackMI",100);
222   fListOfSATracks=new TClonesArray("AliITStrackSA",100);
223   fLayers=new AliITSlayerUpgrade*[6];//to be fixed
224   fCluLayer = 0;
225   fCluCoord = 0;
226   fMinNPoints = 3;
227   //AliITSsegmentationUpgrade *segmentation2 = 0x0;
228   //if(!segmentation2)
229   //  segmentation2 = new AliITSsegmentationUpgrade();
230   for(Int_t layer=0; layer<6; layer++){
231     Double_t p=0.;
232     Double_t zC= 0.;
233     fLayers[layer] = new AliITSlayerUpgrade(p,zC);
234   }
235
236         
237 }
238 //_______________________________________________________________________
239 Int_t AliITStrackerUpgrade::LoadClusters(TTree *clusTree){ 
240 //
241 // Load clusters for tracking
242 // 
243
244   TClonesArray statITSCluster("AliITSRecPoint");
245   TClonesArray *ITSCluster = &statITSCluster;
246
247     TBranch* itsClusterBranch=clusTree->GetBranch("ITSRecPoints");
248     if (!itsClusterBranch){
249       AliError("can't get the branch with the ITS clusters ! \n");
250       return 1;
251     }
252     itsClusterBranch->SetAddress(&ITSCluster);
253     clusTree->GetEvent(0);
254     Int_t nCluster = ITSCluster->GetEntriesFast();
255     for(Int_t i=0; i<nCluster; i++){
256       AliITSRecPoint *recp = (AliITSRecPoint*)ITSCluster->UncheckedAt(i);
257       fLayers[recp->GetLayer()]->InsertCluster(new AliITSRecPoint(*recp));
258     }//loop clusters
259
260
261
262    
263   SetClusterTree(clusTree);
264   return 0;
265 }
266
267 //_______________________________________________________
268 void AliITStrackerUpgrade::ResetForFinding(){
269   //  Reset data members used in all loops during track finding
270   fPhiEstimate=0;
271   for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
272   fLambdac=0;
273   fPhic=0;
274   fCoef1=0;
275   fCoef2=0;
276   fCoef3=0;
277   fPointc[0]=0;
278   fPointc[1]=0;
279   fListOfTracks->Clear();
280   fListOfSATracks->Clear();
281 }
282
283  
284
285 //______________________________________________________________________
286 Int_t AliITStrackerUpgrade::FindTracks(AliESDEvent* event,Bool_t useAllClusters){
287
288   // Track finder using the ESD object
289   AliDebug(2,Form(" field is %f",event->GetMagneticField()));
290   printf(" field is %f",event->GetMagneticField());
291   AliDebug(2,Form("SKIPPING %d %d %d %d %d %d",ForceSkippingOfLayer(0),ForceSkippingOfLayer(1),ForceSkippingOfLayer(2),ForceSkippingOfLayer(3),ForceSkippingOfLayer(4),ForceSkippingOfLayer(5)));
292   if(!fITSclusters){
293     Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
294     return -1;
295   }
296
297   //Reads event and mark clusters of traks already found, with flag kITSin
298   Int_t nentr=event->GetNumberOfTracks();
299   if(!useAllClusters) {
300     while (nentr--) {
301       AliESDtrack *track=event->GetTrack(nentr);
302       if ((track->GetStatus()&AliESDtrack::kITSin) == AliESDtrack::kITSin){
303         Int_t idx[12];
304         Int_t ncl = track->GetITSclusters(idx);
305         for(Int_t k=0;k<ncl;k++){
306           AliITSRecPoint* cll = (AliITSRecPoint*)GetCluster(idx[k]);
307           cll->SetBit(kSAflag);
308         }
309       }
310     }
311   }else{
312     while (nentr--) {
313       AliESDtrack *track=event->GetTrack(nentr);
314       if ((track->GetStatus()&AliESDtrack::kITSin) == AliESDtrack::kITSin){
315         Int_t idx[12];
316         Int_t ncl = track->GetITSclusters(idx);
317         for(Int_t k=0;k<ncl;k++){
318           AliITSRecPoint* cll = (AliITSRecPoint*)GetCluster(idx[k]);
319           cll->ResetBit(kSAflag);
320         }
321       }
322     }
323   }
324   //Get primary vertex
325   Double_t primaryVertex[3];
326   event->GetVertex()->GetXYZ(primaryVertex);
327   Int_t nclusters[/*AliITSsegmentationUpgrade::GetNLayers()*/6]={0,0,0,0,0,0};
328   Int_t dmar[/*AliITSsegmentationUpgrade::GetNLayers()*/6]={0,0,0,0,0,0};
329   if (fCluLayer == 0) {
330     fCluLayer = new TClonesArray*[6];
331     fCluCoord = new TClonesArray*[6];
332     for(Int_t i=0;i<6;i++) {
333       fCluLayer[i]=0;
334       fCluCoord[i]=0;
335     }
336   }
337   for(Int_t i=0;i<6;i++){
338     if (!ForceSkippingOfLayer(i)) {
339       for(Int_t cli=0;cli<fLayers[i]->GetNumberOfClusters();cli++){
340         AliITSRecPoint* cls = (AliITSRecPoint*)fLayers[i]->GetCluster(cli);
341         if(cls->TestBit(kSAflag)==kTRUE) continue; //clusters used by TPC prol.
342         if(cls->GetQ()==0) continue; //fake clusters dead zones
343         if(i>1 && cls->GetQ()<=fMinQ) continue; // cut on SDD and SSD cluster charge
344         nclusters[i]++;
345       }
346     } 
347     dmar[i]=0;
348     if(!fCluLayer[i]){
349       fCluLayer[i] = new TClonesArray("AliITSRecPoint",nclusters[i]);
350     }else{
351       fCluLayer[i]->Delete();
352       fCluLayer[i]->Expand(nclusters[i]);
353     }
354     if(!fCluCoord[i]){
355       fCluCoord[i] = new TClonesArray("AliITSclusterTable",nclusters[i]);
356     }else{
357       fCluCoord[i]->Delete();
358       fCluCoord[i]->Expand(nclusters[i]);
359     }
360   }
361   for(Int_t ilay=0;ilay<AliITSgeomTGeo::GetNLayers();ilay++){
362     TClonesArray &clulay = *fCluLayer[ilay];
363     TClonesArray &clucoo = *fCluCoord[ilay];
364     if (!ForceSkippingOfLayer(ilay)){
365     AliDebug(2,Form("number of clusters in layer %i : %i",ilay,fLayers[ilay]->GetNumberOfClusters()));
366       for(Int_t cli=0;cli<fLayers[ilay]->GetNumberOfClusters();cli++){
367         AliITSRecPoint* cls = (AliITSRecPoint*)fLayers[ilay]->GetCluster(cli);
368         if(cls->TestBit(kSAflag)==kTRUE) continue;
369         if(cls->GetQ()==0) continue;
370         Double_t phi=0;Double_t lambda=0;
371         Double_t x=0;Double_t y=0;Double_t z=0;
372         Float_t sx=0;Float_t sy=0;Float_t sz=0;
373         GetCoorAngles(cls,phi,lambda,x,y,z,primaryVertex);
374         GetCoorErrors(cls,sx,sy,sz);
375         new (clulay[dmar[ilay]]) AliITSRecPoint(*cls);
376         new (clucoo[dmar[ilay]]) AliITSclusterTable(x,y,z,sx,sy,sz,phi,lambda,cli);
377         dmar[ilay]++;
378       }
379     }
380   }
381
382
383   // track counter
384   Int_t ntrack=0;
385
386   static Int_t nClusLay[6];//counter for clusters on each layer
387   Int_t startLayForSeed=0;
388   Int_t lastLayForSeed=fOuterStartLayer;
389   Int_t nSeedSteps=lastLayForSeed-startLayForSeed;
390   Int_t seedStep=1;
391   if(fInwardFlag){
392     startLayForSeed=6-1;
393     lastLayForSeed=fInnerStartLayer;
394     nSeedSteps=startLayForSeed-lastLayForSeed;
395     seedStep=-1;
396   }
397   // loop on minimum number of points
398   for(Int_t iMinNPoints=6; iMinNPoints>=fMinNPoints; iMinNPoints--) {
399     // loop on starting layer for track finding
400     for(Int_t iSeedLay=0; iSeedLay<=nSeedSteps; iSeedLay++) {
401       Int_t theLay=startLayForSeed+iSeedLay*seedStep;
402       if(ForceSkippingOfLayer(theLay)) continue;
403       Int_t minNPoints=iMinNPoints-theLay;
404       if(fInwardFlag) minNPoints=iMinNPoints-(6-1-theLay);
405       for(Int_t i=theLay+1;i<6;i++)
406         if(ForceSkippingOfLayer(i))
407           minNPoints--;
408       if(minNPoints<fMinNPoints) continue;
409
410       // loop on phi and lambda window size
411       for(Int_t nloop=0;nloop<fNloop;nloop++){
412         Int_t nclTheLay=fCluLayer[theLay]->GetEntries();
413         while(nclTheLay--){
414           ResetForFinding();
415           Bool_t useRP=SetFirstPoint(theLay,nclTheLay,primaryVertex);
416           if(!useRP) continue;
417           AliITStrackSA trs;
418
419           Int_t pflag=0;
420           Int_t kk;
421           for(kk=0;kk<6;kk++) nClusLay[kk] = 0;
422
423           kk=0;
424           nClusLay[kk] = SearchClusters(theLay,fPhiWin[nloop],fLambdaWin[nloop],
425                                         &trs,primaryVertex[2],pflag);
426           Int_t nextLay=theLay+seedStep;
427           Bool_t goon=kTRUE;
428           while(goon){
429
430             kk++;
431             nClusLay[kk] = SearchClusters(nextLay,fPhiWin[nloop],fLambdaWin[nloop],
432                                           &trs,primaryVertex[2],pflag);
433             if(nClusLay[kk]!=0){
434
435               pflag=1;
436               if(kk==1) {
437                 fPoint3[0]=fPointc[0];
438                 fPoint3[1]=fPointc[1];
439               } else {
440                 UpdatePoints();
441               }
442             }
443             nextLay+=seedStep;
444             if(nextLay<0 || nextLay==6) goon=kFALSE;
445           }
446
447
448           Int_t layOK=0;
449           if(!fInwardFlag){
450             for(Int_t nnp=0;nnp<6-theLay;nnp++){
451               if(nClusLay[nnp]!=0) layOK+=1;
452             }
453           }else{
454             for(Int_t nnp=theLay; nnp>=0; nnp--){
455               if(nClusLay[nnp]!=0) layOK+=1;
456             }
457           }
458           if(layOK>=minNPoints){
459             AliDebug(2,Form("---NPOINTS: %d; MAP: %d %d %d %d %d %d\n",layOK,nClusLay[0],nClusLay[1],nClusLay[2],nClusLay[3],nClusLay[4],nClusLay[5]));
460             AliITStrackV2* tr2 = 0;
461             Bool_t onePoint = kFALSE;
462             tr2 = FitTrack(&trs,primaryVertex,onePoint);
463             if(!tr2){
464               continue;
465             }
466             AliDebug(2,Form("---NPOINTS fit: %d\n",tr2->GetNumberOfClusters()));
467         
468             StoreTrack(tr2,event,useAllClusters);
469             ntrack++;
470
471           }
472
473         }//end loop on clusters of theLay
474       } //end loop on window sizes
475     } //end loop on theLay
476   }//end loop on min points
477   if(!useAllClusters) AliInfo(Form("Number of found tracks: %d",event->GetNumberOfTracks()));
478
479   ResetForFinding();
480   return 0;
481
482 }
483  
484 //________________________________________________________________________
485
486 AliITStrackV2* AliITStrackerUpgrade::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Bool_t onePoint) {
487
488   
489   const Int_t kMaxClu=AliITStrackSA::kMaxNumberOfClusters;
490
491   
492   static Int_t clind[6][kMaxClu];
493   static Int_t clmark[6][kMaxClu];
494   static Int_t end[6];
495   static AliITSRecPoint *listlayer[6][kMaxClu];
496
497   for(Int_t i=0;i<6;i++) {
498     end[i]=0;
499     for(Int_t j=0;j<kMaxClu; j++){
500       clind[i][j]=0;
501       clmark[i][j]=0;
502       listlayer[i][j]=0;
503     }
504   }
505   Int_t inx[6]; 
506   for (Int_t k=0; k<6; k++) inx[k]=-1;
507   Int_t nclusters = tr->GetNumberOfClustersSA();
508   for(Int_t ncl=0;ncl<nclusters;ncl++){
509     Int_t index = tr->GetClusterIndexSA(ncl); 
510     Int_t lay = (index & 0xf0000000) >> 28;
511     Int_t cli = index&0x0fffffff;
512     AliITSRecPoint* cl = (AliITSRecPoint*)fLayers[lay]->GetCluster(cli);
513     if(cl->TestBit(kSAflag)==kTRUE) cl->ResetBit(kSAflag);  
514     inx[lay]=index;
515  
516     Int_t nInLay=end[lay];
517     listlayer[lay][nInLay]=cl;
518     clind[lay][nInLay]=index;
519     end[lay]++;
520   }
521
522   for(Int_t nlay=0;nlay<6;nlay++){
523     for(Int_t ncl=0;ncl<tr->GetNumberOfMarked(nlay);ncl++){
524       Int_t mark = tr->GetClusterMark(nlay,ncl);
525       clmark[nlay][ncl]=mark;
526     }
527   }
528
529
530   Int_t firstLay=-1,secondLay=-1;
531   for(Int_t i=0;i<6;i++) {
532     if(end[i]==0) {
533       end[i]=1;
534     }else{
535       if(firstLay==-1) {
536         firstLay=i;
537       } else if(secondLay==-1) {
538         secondLay=i;
539       }
540     }
541   }
542
543   if(firstLay==-1 || (secondLay==-1 && !onePoint)) return 0;
544   TClonesArray &arrMI= *fListOfTracks;
545   TClonesArray &arrSA= *fListOfSATracks;
546   Int_t nFoundTracks=0;
547
548   for(Int_t l0=0;l0<end[0];l0++){ //loop on layer 1
549     AliITSRecPoint* cl0 = (AliITSRecPoint*)listlayer[0][l0];
550     for(Int_t l1=0;l1<end[1];l1++){ //loop on layer 2
551       AliITSRecPoint* cl1 = (AliITSRecPoint*)listlayer[1][l1];
552       for(Int_t l2=0;l2<end[2];l2++){  //loop on layer 3
553         AliITSRecPoint* cl2 = (AliITSRecPoint*)listlayer[2][l2];
554         for(Int_t l3=0;l3<end[3];l3++){ //loop on layer 4   
555           AliITSRecPoint* cl3 = (AliITSRecPoint*)listlayer[3][l3];
556           for(Int_t l4=0;l4<end[4];l4++){ //loop on layer 5
557             AliITSRecPoint* cl4 = (AliITSRecPoint*)listlayer[4][l4];
558             for(Int_t l5=0;l5<end[5];l5++){ //loop on layer 6  
559               AliITSRecPoint* cl5 = (AliITSRecPoint*)listlayer[5][l5];
560
561
562               Double_t x1,y1,z1,sx1,sy1,sz1;
563               Double_t x2,y2,z2,sx2,sy2,sz2;
564               AliITSRecPoint* p1=0;
565               AliITSRecPoint* p2=0;
566               Int_t index1=0,index2=0;
567               Int_t mrk1=0,mrk2=0;
568
569               switch(firstLay) {
570               case 0:
571                 p1=cl0;
572                 index1=clind[0][l0];mrk1=clmark[0][l0];
573                 break;
574               case 1:
575                 p1=cl1;
576                 index1=clind[1][l1];mrk1=clmark[1][l1];
577                 break;
578               case 2:
579                 p1=cl2;
580                 index1=clind[2][l2];mrk1=clmark[2][l2];
581                 break;
582               case 3:
583                 p1=cl3;
584                 index1=clind[3][l3];mrk1=clmark[3][l3];
585                 break;
586               case 4:
587                 p1=cl4;
588                 index1=clind[4][l4];mrk1=clmark[4][l4];
589                 break;
590               }
591
592               switch(secondLay) {
593               case 1:
594                 p2=cl1;
595                 index2=clind[1][l1];mrk2=clmark[1][l1];
596                 break;
597               case 2:
598                 p2=cl2;
599                 index2=clind[2][l2];mrk2=clmark[2][l2];
600                 break;
601               case 3:
602                 p2=cl3;
603                 index2=clind[3][l3];mrk2=clmark[3][l3];
604                 break;
605               case 4:
606                 p2=cl4;
607                 index2=clind[4][l4];mrk2=clmark[4][l4];
608                 break;
609               case 5:
610                 p2=cl5;
611                 index2=clind[5][l5];mrk2=clmark[5][l5];
612                 break;
613               default:
614                 p2=0;
615                 index2=-1;mrk2=-1;
616                 break;
617               }
618
619               Int_t layer = 0;
620               Double_t radius = 0.;
621               Double_t yclu1 = p1->GetY();
622               Double_t zclu1 = p1->GetZ();
623               layer=p1->GetLayer();
624               AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade();
625               radius=segmentation->GetRadius(layer);
626               Double_t cv=0.,tgl2=0.,phi2=0.;
627               Int_t cln1=mrk1;
628               AliITSclusterTable* arr1 = (AliITSclusterTable*)GetClusterCoord(firstLay,cln1);
629               x1 = arr1->GetX();
630               y1 = arr1->GetY();
631               z1 = arr1->GetZ();
632               sx1 = arr1->GetSx();
633               sy1 = arr1->GetSy();
634               sz1 = arr1->GetSz();
635
636               if(secondLay>0) {
637                 Int_t cln2=mrk2;
638                 AliITSclusterTable* arr2 = (AliITSclusterTable*)GetClusterCoord(secondLay,cln2);
639                 x2 = arr2->GetX();
640                 y2 = arr2->GetY();
641                 z2 = arr2->GetZ();
642                 sx2 = arr2->GetSx();
643                 sy2 = arr2->GetSy();
644                 sz2 = arr2->GetSz();
645                 cv = Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
646                 tgl2 = (z2-z1)/TMath::Sqrt((x2-x1)*(x2-x1)+(y2-y1)*(y2-y1));
647                 phi2 = TMath::ATan2((y2-y1),(x2-x1));
648               }
649               Double_t x,y,z = 0.;
650               Double_t xz[2];
651               xz[0]= p1->GetDetLocalX(); 
652               xz[1]= p1->GetDetLocalZ();   
653               Bool_t check2;
654               check2 = segmentation->DetToGlobal(layer,xz[0], xz[1],x,y,z);
655
656
657               Double_t phiclrad, phicldeg;
658               phiclrad=TMath::ATan2(y,x); 
659               if (phiclrad<0) phiclrad+=TMath::TwoPi();
660               else if (phiclrad>=TMath::TwoPi()) phiclrad-=TMath::TwoPi();
661
662               phicldeg=180.*phiclrad/TMath::Pi();
663               Int_t ladder=(Int_t)(phicldeg/18.);//virtual segmentation
664               Double_t alpha= (ladder*18.+9.)*TMath::Pi()/180.;
665
666
667               AliITStrackSA trac(alpha,radius,yclu1,zclu1,phi2,tgl2,cv,1);
668
669
670               if(cl5!=0) {
671                 trac.AddClusterV2(5,(clind[5][l5] & 0x0fffffff)>>0);
672                 trac.AddClusterMark(5,clmark[5][l5]);
673               }
674               if(cl4!=0){
675                 trac.AddClusterV2(4,(clind[4][l4] & 0x0fffffff)>>0);
676                 trac.AddClusterMark(4,clmark[4][l4]);
677               }
678               if(cl3!=0){
679                 trac.AddClusterV2(3,(clind[3][l3] & 0x0fffffff)>>0);
680                 trac.AddClusterMark(3,clmark[3][l3]);
681               }
682               if(cl2!=0){
683                 trac.AddClusterV2(2,(clind[2][l2] & 0x0fffffff)>>0);
684                 trac.AddClusterMark(2,clmark[2][l2]);
685               }
686               if(cl1!=0){
687                 trac.AddClusterV2(1,(clind[1][l1] & 0x0fffffff)>>0);
688                 trac.AddClusterMark(1,clmark[1][l1]);
689               }
690               if(cl0!=0){
691                 trac.AddClusterV2(0,(clind[0][l0] & 0x0fffffff)>>0);
692                 trac.AddClusterMark(0,clmark[0][l0]);
693               }
694   
695               //fit with Kalman filter using AliITStrackerMI::RefitAtBase()
696               AliITStrackSA ot(trac);
697               ot.ResetCovariance(10.);
698               ot.ResetClusters();
699               
700               // Propagate inside the innermost layer with a cluster 
701               if(ot.Propagate(ot.GetX()-0.1*ot.GetX())) {
702                 Double_t rt=0.;
703                 rt=segmentation->GetRadius(5);
704                 if(RefitAtBase(rt,&ot,inx)){ //fit from layer 1 to layer 6
705                   AliITStrackMI otrack2(ot);
706                   otrack2.ResetCovariance(10.); 
707                   otrack2.ResetClusters();
708                   
709                   //fit from layer 6 to layer 1
710                   if(RefitAtBase(/*AliITSRecoParam::GetrInsideSPD1()*/segmentation->GetRadius(0),&otrack2,inx)) {//check clind
711                     new(arrMI[nFoundTracks]) AliITStrackMI(otrack2);
712                     new(arrSA[nFoundTracks]) AliITStrackSA(trac);
713                     ++nFoundTracks;
714                   }
715                               
716                 }       
717               }
718             }//end loop layer 6
719           }//end loop layer 5
720         }//end loop layer 4        
721       }//end loop layer 3
722     }//end loop layer 2 
723   }//end loop layer 1
724
725
726
727
728   if(fListOfTracks->GetEntries()==0) return 0;
729
730   Int_t lowchi2 = FindTrackLowChiSquare();
731   AliITStrackV2* otrack =(AliITStrackV2*)fListOfTracks->At(lowchi2);
732   AliITStrackSA* trsa = (AliITStrackSA*)fListOfSATracks->At(lowchi2);
733  
734   if(otrack==0) {
735     return 0;
736   }
737   Int_t indexc[6];
738   for(Int_t i=0;i<6;i++) indexc[i]=0;
739   for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
740     indexc[nind] = otrack->GetClusterIndex(nind);
741   }      
742  
743   //remove clusters of found track
744   for(Int_t nlay=0;nlay<6;nlay++){
745     for(Int_t cln=0;cln<trsa->GetNumberOfMarked(nlay);cln++){
746       Int_t index = trsa->GetClusterMark(nlay,cln);
747       fCluLayer[nlay]->RemoveAt(index);
748       RemoveClusterCoord(nlay,index);
749       fCluLayer[nlay]->Compress();
750     }    
751   }
752   return otrack;
753
754 }
755
756 //_______________________________________________________
757 void AliITStrackerUpgrade::StoreTrack(AliITStrackV2 *t,AliESDEvent *event, Bool_t pureSA) /*const*/ 
758 {
759   //
760   // Add new track to the ESD
761   //
762   AliESDtrack outtrack;
763   outtrack.UpdateTrackParams(t,AliESDtrack::kITSin);
764   if(pureSA) outtrack.SetStatus(AliESDtrack::kITSpureSA);
765   outtrack.SetStatus(AliESDtrack::kITSout);
766   outtrack.SetStatus(AliESDtrack::kITSrefit);
767   Double_t sdedx[4]={0.,0.,0.,0.};
768   for(Int_t i=0; i<4; i++) sdedx[i]=t->GetSampledEdx(i);
769   outtrack.SetITSdEdxSamples(sdedx);
770   event->AddTrack(&outtrack);
771   return;
772 }
773
774
775 //_______________________________________________________
776 Int_t AliITStrackerUpgrade::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t /*zvertex*/,Int_t pflag){
777   //function used to to find the clusters associated to the track
778
779    AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade();
780
781   AliDebug(2,"Starting...");
782   if(ForceSkippingOfLayer(layer)) {
783   AliDebug(2,Form("Forcing skipping of layer %i. Exiting",layer));
784   return 0;
785   }
786
787   Int_t nc=0;
788   Double_t r=segmentation->GetRadius(layer);
789   if(pflag==1){      
790     Float_t cx1,cx2,cy1,cy2;
791     FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
792     if (FindIntersection(fCoef1,fCoef2,fCoef3,-r*r,cx1,cy1,cx2,cy2)==0)
793       return 0;
794     Double_t fi1=TMath::ATan2(cy1-fPoint1[1],cx1-fPoint1[0]);
795     Double_t fi2=TMath::ATan2(cy2-fPoint1[1],cx2-fPoint1[0]);
796     fPhiEstimate=ChoosePoint(fi1,fi2,fPhic);
797   }
798
799  
800   Int_t ncl = fCluLayer[layer]->GetEntries();
801   AliDebug(2,Form(" Number of clusters %i in layer %i.",ncl,layer));
802   for (Int_t index=0; index<ncl; index++) {
803     AliITSRecPoint *c = (AliITSRecPoint*)fCluLayer[layer]->At(index);
804     if (!c) continue;
805     if (c->GetQ()<=0) continue;
806     if(layer>1 && c->GetQ()<=fMinQ) continue;
807     
808     AliITSclusterTable* arr = (AliITSclusterTable*)GetClusterCoord(layer,index);
809     
810     Double_t phi;
811     Double_t xT,yT;
812     xT=arr->GetX();
813     yT=arr->GetY();
814     phi=arr->GetPhi();
815     if (TMath::Abs(phi-fPhiEstimate)>phiwindow)
816       {
817         continue;
818       }
819
820     Double_t lambda = arr->GetLambda();
821     if (TMath::Abs(lambda-fLambdac)>lambdawindow) continue;
822
823     if(trs->GetNumberOfClustersSA()==trs->GetMaxNumberOfClusters()) return 0;
824     if(trs->GetNumberOfMarked(layer)==trs->GetMaxNMarkedPerLayer()) return 0;
825     Int_t orind = arr->GetOrInd();
826     trs->AddClusterSA(layer,orind);
827     trs->AddClusterMark(layer,index);
828        
829     nc++;
830     fLambdac=lambda;
831     fPhiEstimate=phi;
832
833     fPointc[0]=arr->GetX();
834     fPointc[1]=arr->GetY();
835   }
836   delete segmentation;
837   return nc;
838 }
839
840 //________________________________________________________________
841 Bool_t AliITStrackerUpgrade::SetFirstPoint(Int_t lay, Int_t clu, Double_t* primaryVertex){
842   // Sets the first point (seed) for tracking
843
844   AliITSRecPoint* cl = (AliITSRecPoint*)fCluLayer[lay]->At(clu);
845   if(!cl) return kFALSE;
846   if (cl->GetQ()<=0) return kFALSE;
847   if(lay>1 && cl->GetQ()<=fMinQ) return kFALSE;
848   AliITSclusterTable* arr = (AliITSclusterTable*)GetClusterCoord(lay,clu);
849   fPhic = arr->GetPhi();
850   fLambdac = arr->GetLambda();
851   fPhiEstimate = fPhic;
852   fPoint1[0]=primaryVertex[0];
853   fPoint1[1]=primaryVertex[1];
854   fPoint2[0]=arr->GetX();
855   fPoint2[1]=arr->GetY();
856   return kTRUE; 
857 }
858
859 //________________________________________________________________
860 void AliITStrackerUpgrade::UpdatePoints(){
861   //update of points for the estimation of the curvature  
862
863   fPoint2[0]=fPoint3[0];
864   fPoint2[1]=fPoint3[1];
865   fPoint3[0]=fPointc[0];
866   fPoint3[1]=fPointc[1];
867
868   
869 }
870
871 //___________________________________________________________________
872 Int_t AliITStrackerUpgrade::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){
873
874   //given (x,y) of three recpoints (in global coordinates) 
875   //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
876
877   Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
878   if(den==0) return 0;
879   a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
880   b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
881   c = -x1*x1-y1*y1-a*x1-b*y1;
882   return 1;
883 }
884 //__________________________________________________________________________
885 Int_t AliITStrackerUpgrade::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){
886  
887   //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
888   //c2 is -rlayer*rlayer
889
890   if(a1==0) return 0;
891   Double_t m = c2-c1; 
892   Double_t aA = (b1*b1)/(a1*a1)+1;
893   Double_t bB = (-2*m*b1/(a1*a1));
894   Double_t cC = c2+(m*m)/(a1*a1);
895   Double_t dD = bB*bB-4*aA*cC;
896   if(dD<0) return 0;
897  
898   y1 = (-bB+TMath::Sqrt(dD))/(2*aA); 
899   y2 = (-bB-TMath::Sqrt(dD))/(2*aA); 
900   x1 = (c2-c1-b1*y1)/a1;
901   x2 = (c2-c1-b1*y2)/a1;
902
903   return 1; 
904 }
905 //____________________________________________________________________
906 Double_t AliITStrackerUpgrade::Curvature(Double_t x1,Double_t y1,Double_t 
907                                          x2,Double_t y2,Double_t x3,Double_t y3){
908
909   //calculates the curvature of track  
910   Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
911   if(den==0) return 0;
912   Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
913   Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
914   Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
915   Double_t xc=-a/2.;
916
917   if((a*a+b*b-4*c)<0) return 0;
918   Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
919   if(rad==0) return 0;
920   
921   if((x1>0 && y1>0 && x1<xc)) rad*=-1;
922   if((x1<0 && y1>0 && x1<xc)) rad*=-1;
923   //  if((x1<0 && y1<0 && x1<xc)) rad*=-1;
924   // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
925   
926   return 1/rad;
927  
928 }
929
930
931 //____________________________________________________________________
932 Double_t AliITStrackerUpgrade::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
933
934   //Returns the point closest to pp
935
936   Double_t diff1 = p1-pp;
937   Double_t diff2 = p2-pp;
938   
939   if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
940   else fPhiEstimate=p2;  
941   return fPhiEstimate;
942   
943 }
944
945
946 //_________________________________________________________________
947 Int_t AliITStrackerUpgrade::FindTrackLowChiSquare() const {
948   // returns track with lowest chi square  
949   Int_t dim=fListOfTracks->GetEntries();
950   if(dim<=1) return 0;
951   AliITStrackV2* trk = (AliITStrackV2*)fListOfTracks->At(0);
952   Double_t minChi2=trk->GetChi2();
953   Int_t index=0;
954   for(Int_t i=1;i<dim;i++){
955     trk = (AliITStrackV2*)fListOfTracks->At(i);
956     Double_t chi2=trk->GetChi2();
957     if(chi2<minChi2){
958       minChi2=chi2;
959       index=i;
960     }
961   }
962   return index;
963 }
964
965 //__________________________________________________________
966 /*Int_t AliITStrackerUpgrade::FindLabel(Int_t l0, Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5){
967
968 //function used to determine the track label
969   
970 Int_t lb[6] = {l0,l1,l2,l3,l4,l5};
971 Int_t aa[6]={1,1,1,1,1,1};
972 Int_t ff=0; 
973 Int_t ll=0;
974 Int_t k=0;Int_t w=0;Int_t num=6;
975 for(Int_t i=5;i>=0;i--) if(lb[i]==-1) num=i;
976   
977 while(k<num){
978   
979 for(Int_t i=k+1;i<num;i++){
980     
981 if(lb[k]==lb[i] && aa[k]!=0){
982       
983 aa[k]+=1;
984 aa[i]=0;
985 }
986 }
987 k++;
988 }
989
990 while(w<num){
991  
992 for(Int_t j=0;j<6;j++){
993 if(aa[w]<aa[j]){
994 ff=aa[w];
995 aa[w]=aa[j];
996 aa[j]=ff;
997 ll=lb[w];
998 lb[w]=lb[j];
999 lb[j]=ll;
1000 }
1001 }
1002 w++;
1003 }
1004   
1005 if(num<1) return -1; 
1006 return lb[num-1];
1007 }
1008 */
1009 //_____________________________________________________________________________
1010 /*Int_t AliITStrackerUpgrade::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,
1011   Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t minNPoints){
1012
1013  
1014   //function used to assign label to the found track. If track is fake, the label is negative
1015
1016   Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1017   Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1018   Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1019   Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1020   Int_t lflag=0;Int_t num=6;
1021   if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1022
1023   for(Int_t i=0;i<num;i++){
1024   if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1025   }
1026
1027   if(lflag>=minNPoints) return ll;
1028   else return -ll;
1029
1030   
1031   }
1032 */
1033 //_____________________________________________________________________________
1034 void AliITStrackerUpgrade::SetCalculatedWindowSizes(Int_t n, Float_t phimin, Float_t phimax, Float_t lambdamin, Float_t lambdamax){
1035   // Set sizes of the phi and lambda windows used for track finding
1036   fNloop = n;
1037   if(fPhiWin) delete [] fPhiWin;
1038   if(fLambdaWin) delete [] fLambdaWin;
1039   fPhiWin = new Double_t[fNloop];
1040   fLambdaWin = new Double_t[fNloop];
1041   Float_t stepPhi=(phimax-phimin)/(Float_t)(fNloop-1);
1042   Float_t stepLambda=(lambdamax-lambdamin)/(Float_t)(fNloop-1);
1043   for(Int_t k=0;k<fNloop;k++){
1044     Float_t phi=phimin+k*stepPhi;
1045     Float_t lam=lambdamin+k*stepLambda;
1046     fPhiWin[k]=phi;
1047     fLambdaWin[k]=lam;
1048   }
1049 }
1050 //_____________________________________________________________________________
1051 void AliITStrackerUpgrade::SetFixedWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1052   // Set sizes of the phi and lambda windows used for track finding
1053   fNloop = n;
1054   if(phi){ // user defined values
1055     fPhiWin = new Double_t[fNloop];
1056     fLambdaWin = new Double_t[fNloop];
1057     for(Int_t k=0;k<fNloop;k++){
1058       fPhiWin[k]=phi[k];
1059       fLambdaWin[k]=lam[k];
1060     }
1061   }
1062   else {  // default values
1063             
1064     Double_t phid[33]   = {0.002,0.003,0.004,0.0045,0.0047,
1065                            0.005,0.0053,0.0055,
1066                            0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1067                            0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1068                            0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,0.0140,0.0145};
1069     Double_t lambdad[33] = {0.003,0.004,0.005,0.005,0.005,
1070                             0.005,0.005,0.006,
1071                             0.006,0.006,0.006,0.007,0.007,0.007,0.007,
1072                             0.007,0.007,0.007,0.007,0.007,0.007,0.007,
1073                             0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008,0.008};
1074     
1075     if(fNloop!=33){
1076       fNloop = 33;
1077     }
1078     
1079     
1080     fPhiWin = new Double_t[fNloop];
1081     fLambdaWin = new Double_t[fNloop];
1082
1083     Double_t factor=AliITSReconstructor::GetRecoParam()->GetFactorSAWindowSizes(); 
1084   
1085     for(Int_t k=0;k<fNloop;k++){
1086       fPhiWin[k]=phid[k]*factor;
1087       fLambdaWin[k]=lambdad[k]*factor;
1088     }
1089   
1090   }
1091
1092 }
1093 //_______________________________________________________________________
1094 void AliITStrackerUpgrade::GetCoorAngles(AliITSRecPoint* cl,Double_t &phi,Double_t &lambda, Double_t &x, Double_t &y,Double_t &z,Double_t* vertex){
1095   //Returns values of phi (azimuthal) and lambda angles for a given cluster
1096   /*  
1097       Double_t rot[9];     fGeom->GetRotMatrix(module,rot);
1098       Int_t lay,lad,det; fGeom->GetModuleId(module,lay,lad,det);
1099       Float_t tx,ty,tz;  fGeom->GetTrans(lay,lad,det,tx,ty,tz);     
1100
1101       Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
1102       Double_t phi1=TMath::Pi()/2+alpha;
1103       if (lay==1) phi1+=TMath::Pi();
1104
1105       Float_t cp=TMath::Cos(phi1), sp=TMath::Sin(phi1);
1106       Float_t r=tx*cp+ty*sp;
1107
1108       xyz= r*cp - cl->GetY()*sp;
1109       y= r*sp + cl->GetY()*cp;
1110       z=cl->GetZ();
1111   */
1112   
1113   Double_t xz[2];
1114   xz[0]= cl->GetDetLocalX(); 
1115   xz[1]= cl->GetDetLocalZ() ; 
1116   AliITSsegmentationUpgrade *segmentation2 = new AliITSsegmentationUpgrade();
1117   Bool_t check2;
1118   Int_t ilayer;
1119   ilayer = cl->GetLayer();
1120   check2 = segmentation2->DetToGlobal(ilayer,xz[0], xz[1],x,y,z);
1121
1122   if(x!=0 && y!=0)  
1123     phi=TMath::ATan2(y-vertex[1],x-vertex[0]);
1124   lambda=TMath::ATan2(z-vertex[2],TMath::Sqrt((x-vertex[0])*(x-vertex[0])+(y-vertex[1])*(y-vertex[1])));
1125   if(segmentation2) delete segmentation2;
1126 }
1127
1128 //________________________________________________________________________
1129 void AliITStrackerUpgrade::GetCoorErrors(AliITSRecPoint* cl,Float_t &sx,Float_t &sy, Float_t &sz){
1130
1131   //returns sigmax, y, z of cluster in global coordinates
1132   /*
1133     Double_t rot[9];     fGeom->GetRotMatrix(module,rot);
1134     Int_t lay,lad,det; 
1135     AliITSgeomTGeo::GetModuleId(module,lay,lad,det);
1136  
1137     Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
1138     Double_t phi=TMath::Pi()/2+alpha;
1139     if (lay==1) phi+=TMath::Pi();
1140
1141     Float_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
1142   */
1143   Float_t covm[6];
1144   cl->GetGlobalCov(covm);
1145   sx=TMath::Sqrt(covm[0]);
1146   sy=TMath::Sqrt(covm[3]);
1147   sz=TMath::Sqrt(covm[5]);
1148
1149
1150
1151   //  Float_t covm[6];
1152   //  cl->GetGlobalCov(covm);
1153   //  sx=12*1e-04;//TMath::Sqrt(covm[0]);
1154   //  sy=12*1e-04;//TMath::Sqrt(covm[3]);
1155   //  sz=40*1e-04;//TMath::Sqrt(covm[5]);
1156   /*
1157     sx = TMath::Sqrt(sp*sp*cl->GetSigmaY2());
1158     sy = TMath::Sqrt(cp*cp*cl->GetSigmaY2());
1159     sz = TMath::Sqrt(cl->GetSigmaZ2());
1160   */
1161 }
1162 //_____________________________________________________________
1163 void AliITStrackerUpgrade::UnloadClusters() {
1164   //--------------------------------------------------------------------
1165   //This function unloads ITS clusters
1166   //--------------------------------------------------------------------
1167   for (Int_t i=0; i<6; i++) fLayers[i]->ResetClusters();
1168 }
1169 //______________________________________________________________________
1170 Bool_t AliITStrackerUpgrade::RefitAtBase(Double_t xx,AliITStrackMI *track,
1171                                          const Int_t *clusters)
1172 {
1173   //--------------------------------------------------------------------
1174   // Simplified version for ITS upgrade studies -- does not use module info
1175   //--------------------------------------------------------------------
1176   Int_t index[AliITSgeomTGeo::kNLayers];
1177   Int_t k;
1178   for (k=0; k<AliITSgeomTGeo::GetNLayers(); k++) index[k]=-1;
1179   //
1180   for (k=0; k<AliITSgeomTGeo::GetNLayers(); k++) {
1181     index[k]=clusters[k];
1182   }
1183   AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade();
1184
1185   ULong_t trStatus=0;
1186   if(track->GetESDtrack()) trStatus=track->GetStatus();
1187   Int_t innermostlayer=0;
1188   if(trStatus&AliESDtrack::kTPCin)  {
1189     innermostlayer=5;
1190     Double_t drphi = TMath::Abs(track->GetD(0.,0.));
1191     for(innermostlayer=0; innermostlayer<AliITSgeomTGeo::GetNLayers(); innermostlayer++) {
1192       if( (drphi < (segmentation->GetRadius(innermostlayer)+1.)) ||
1193           index[innermostlayer] >= 0 ) break;
1194     }
1195     AliDebug(2,Form(" drphi  %f  innermost %d",drphi,innermostlayer));
1196   }
1197
1198   Int_t from, to, step;
1199   if (xx > track->GetX()) {
1200     from=innermostlayer; to=6;
1201     step=+1;
1202   } else {
1203     from=AliITSgeomTGeo::GetNLayers()-1; to=innermostlayer-1;
1204     step=-1;
1205   }
1206   TString dir = (step>0 ? "outward" : "inward");
1207
1208   for (Int_t ilayer = from; ilayer != to; ilayer += step) {
1209     Double_t r=0.;
1210     r=segmentation->GetRadius(ilayer);
1211
1212     if (step<0 && xx>r){
1213       break;
1214     }
1215
1216     // material between SSD and SDD, SDD and SPD
1217     //      Double_t hI=ilayer-0.5*step; 
1218     //      if (TMath::Abs(hI-3.5)<0.01) // SDDouter
1219     //        if(!CorrectForShieldMaterial(track,"SDD",dir)) return kFALSE;
1220     //      if (TMath::Abs(hI-1.5)<0.01) // SPDouter
1221     //        if(!CorrectForShieldMaterial(track,"SPD",dir)) return kFALSE;
1222
1223
1224     Double_t oldGlobXYZ[3];
1225
1226     if (!track->GetXYZ(oldGlobXYZ)) {
1227       if(segmentation) delete segmentation;
1228       return kFALSE;
1229     }
1230     // continue if we are already beyond this layer
1231     Double_t oldGlobR = TMath::Sqrt(oldGlobXYZ[0]*oldGlobXYZ[0]+oldGlobXYZ[1]*oldGlobXYZ[1]);
1232     if(step>0 && oldGlobR > r){
1233       continue; // going outward
1234     }
1235     if(step<0 && oldGlobR < r) {
1236       continue; // going inward
1237     }
1238     Double_t phi,z;
1239     if (!track->GetPhiZat(r,phi,z)){
1240       if(segmentation) delete segmentation;
1241       return kFALSE;
1242     }
1243     // only for ITS-SA tracks refit
1244     track->SetCheckInvariant(kFALSE);
1245       
1246     if (phi<0) phi+=TMath::TwoPi();//from 0 to 360 (rad) 
1247     else if (phi>=TMath::TwoPi()) phi-=TMath::TwoPi();//
1248     Double_t phideg=0.;  
1249     phideg=180.*phi/TMath::Pi();// in deg  
1250     Int_t ladder=(Int_t)(phideg/18.);// virtual segmentation 
1251     Double_t alpha= (ladder*18.+9.)*TMath::Pi()/180.;
1252     
1253     if (!track->Propagate(alpha,r)) {
1254       if(segmentation) delete segmentation;
1255       return kFALSE;
1256     }
1257
1258     const AliITSRecPoint *clAcc=0;
1259     Double_t maxchi2=1000.*AliITSReconstructor::GetRecoParam()->GetMaxChi2();
1260
1261     Int_t idx=index[ilayer];
1262     if (idx>=0) { // cluster in this layer   
1263       Int_t cli = idx&0x0fffffff;
1264       const AliITSRecPoint *cl=(AliITSRecPoint *)fLayers[ilayer]->GetCluster(cli);
1265
1266       if (cl) {                                                                  
1267         Int_t cllayer = (idx & 0xf0000000) >> 28;;                               
1268         Double_t chi2=GetPredictedChi2MI(track,cl,cllayer);                      
1269         if (chi2<maxchi2) {                                                      
1270           clAcc=cl;                                                              
1271           maxchi2=chi2;                                                          
1272         } else {                                                                 
1273            if(segmentation) delete segmentation;
1274            return kFALSE;                                                         
1275         }                                                                        
1276       }                                                                          
1277     }              
1278
1279     if (clAcc) {
1280       if (!UpdateMI(track,clAcc,maxchi2,idx)){
1281         if(segmentation) delete segmentation;
1282         return kFALSE;
1283       }
1284       track->SetSampledEdx(clAcc->GetQ(),ilayer-2);
1285     }
1286
1287
1288     // Correct for material of the current layer
1289     // cross material
1290     // add time if going outward
1291     if(!CorrectForLayerMaterial(track,ilayer,oldGlobXYZ,dir)){
1292       if(segmentation) delete segmentation;
1293       return kFALSE;
1294     }
1295     track->SetCheckInvariant(kFALSE);
1296   } // end loop on layers
1297
1298   if (!track->PropagateTo(xx,0.,0.)){
1299     if(segmentation) delete segmentation;
1300     return kFALSE;
1301   } 
1302   if(segmentation) delete segmentation;
1303   return kTRUE;
1304 }
1305
1306 //_____________________________________________________________________
1307 Int_t AliITStrackerUpgrade::UpdateMI(AliITStrackMI* track, const AliITSRecPoint* cl,Double_t chi2,Int_t index) const
1308 {
1309   //
1310   // Update ITS track
1311   //
1312   Int_t layer = (index & 0xf0000000) >> 28;
1313   track->SetClIndex(layer, index);
1314   if (TMath::Abs(cl->GetQ())<1.e-13) return 0;  // ingore the "virtual" clusters
1315
1316
1317   // Take into account the mis-alignment (bring track to cluster plane)
1318   Double_t xTrOrig=track->GetX();
1319   Double_t trxyz[3]; track->GetXYZ(trxyz);
1320   Float_t xclu1 = cl->GetX();//upgrade clusters tracking coordinate
1321   //Float_t yclu1 = cl->GetY();
1322   //Float_t zclu1 = cl->GetZ();
1323   if (!track->Propagate(xTrOrig+xclu1)){
1324     return 0;
1325   }
1326   AliCluster c(*cl);
1327   c.SetSigmaY2(track->GetSigmaY(layer)*track->GetSigmaY(layer));
1328   c.SetSigmaZ2(track->GetSigmaZ(layer)*track->GetSigmaZ(layer));
1329   c.SetSigmaYZ(track->GetSigmaYZ(layer));
1330
1331   // Bring the track back to detector plane in ideal geometry  
1332   Int_t updated = track->UpdateMI(&c,chi2,index);
1333   if (!track->Propagate(xTrOrig)) {
1334     return 0;
1335   }
1336   if(!updated){
1337     AliDebug(2,"update failed");
1338   }
1339   return updated;
1340 }
1341 //____________________________________________________________________-
1342 Int_t AliITStrackerUpgrade::CorrectForLayerMaterial(AliITStrackMI *t,
1343                                                     Int_t layerindex,
1344                                                     Double_t oldGlobXYZ[3],
1345                                                     TString direction) {
1346   //-------------------------------------------------------------------
1347   // Propagate beyond layer and correct for material
1348   // (material budget in different ways according to fUseTGeo value)
1349   // Add time if going outward (PropagateTo or PropagateToTGeo)
1350   //-------------------------------------------------------------------
1351
1352   // Define budget mode:
1353   // 0: material from AliITSRecoParam (hard coded)
1354   // 1: material from TGeo in stepsof X cm (on the fly)
1355   //    X = AliITSRecoParam::GetStepSizeTGeo()
1356   // 2: material from lut
1357   // 3: material from TGeo in one step (same for all hypotheses)
1358   Int_t mode;
1359   switch(fUseTGeo) {
1360   case 0:
1361     mode=0;
1362     break;
1363   case 1:
1364     mode=1;
1365     break;
1366   case 2:
1367     mode=2;
1368     break;
1369   case 3:
1370     if(fTrackingPhase.Contains("Clusters2Tracks"))
1371       { mode=3; } else { mode=1; }
1372     break;
1373   case 4:
1374     if(fTrackingPhase.Contains("Clusters2Tracks"))
1375       { mode=3; } else { mode=2; }
1376     break;
1377   default:
1378     mode=0;
1379     break;
1380   }
1381   if(fTrackingPhase.Contains("Default")) {
1382     mode=0;
1383   }
1384   Float_t  dir = (direction.Contains("inward") ? 1. : -1.);
1385   Double_t r = 0.;
1386   AliITSsegmentationUpgrade *segmentation = 0x0;
1387   if(!segmentation)
1388     segmentation = new AliITSsegmentationUpgrade();
1389   r=segmentation->GetRadius(layerindex);
1390   Double_t deltar=(layerindex<2 ? 0.10*r : 0.05*r);
1391   Double_t rToGo=TMath::Sqrt(t->GetX()*t->GetX()+t->GetY()*t->GetY())-deltar*dir;
1392   Double_t xToGo;
1393   if (!t->GetLocalXat(rToGo,xToGo)) {
1394     return 0;
1395   }
1396   Int_t index=6*fCurrentEsdTrack+layerindex;
1397
1398
1399   Double_t xOverX0=0.0,x0=0.0,lengthTimesMeanDensity=0.0;
1400   Int_t nsteps=1;
1401   // back before material (no correction)
1402   Double_t rOld,xOld;
1403   rOld=TMath::Sqrt(oldGlobXYZ[0]*oldGlobXYZ[0]+oldGlobXYZ[1]*oldGlobXYZ[1]);
1404   if (!t->GetLocalXat(rOld,xOld)) return 0;
1405   if (!t->Propagate(xOld)) return 0;
1406
1407   switch(mode) {
1408   case 1:
1409     x0=21.82;
1410     xOverX0 = segmentation->GetThickness(layerindex)/x0;
1411
1412     lengthTimesMeanDensity = xOverX0*x0;
1413     lengthTimesMeanDensity *= dir;
1414     // Bring the track beyond the material
1415     if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)){
1416       return 0;
1417     }
1418     break;
1419   case 0:
1420     nsteps = (Int_t)(TMath::Abs(xOld-xToGo)/AliITSReconstructor::GetRecoParam()->GetStepSizeTGeo())+1;
1421     if (!t->PropagateToTGeo(xToGo,nsteps)) return 0; // cross the material and apply correction
1422     break;
1423   case 2:
1424     if(fxOverX0Layer[layerindex]<0) BuildMaterialLUT("Layers");
1425     xOverX0 = fxOverX0Layer[layerindex];
1426     lengthTimesMeanDensity = fxTimesRhoLayer[layerindex];
1427     lengthTimesMeanDensity *= dir;
1428     // Bring the track beyond the material
1429     if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)) return 0;
1430     break;
1431   case 3:
1432     if(!fxOverX0LayerTrks || index<0 || index>=6*fNtracks) Error("CorrectForLayerMaterial","Incorrect usage of UseTGeo option!\n");
1433     if(fxOverX0LayerTrks[index]<0) {
1434       nsteps = (Int_t)(TMath::Abs(xOld-xToGo)/AliITSReconstructor::GetRecoParam()->GetStepSizeTGeo())+1;
1435       if (!t->PropagateToTGeo(xToGo,nsteps,xOverX0,lengthTimesMeanDensity)) return 0;
1436       Double_t angle=TMath::Sqrt((1.+t->GetTgl()*t->GetTgl())/
1437                                  ((1.-t->GetSnp())*(1.+t->GetSnp())));
1438       fxOverX0LayerTrks[index] = TMath::Abs(xOverX0)/angle;
1439       fxTimesRhoLayerTrks[index] = TMath::Abs(lengthTimesMeanDensity)/angle;
1440       return 1;
1441     }
1442     xOverX0 = fxOverX0LayerTrks[index];
1443     lengthTimesMeanDensity = fxTimesRhoLayerTrks[index];
1444     lengthTimesMeanDensity *= dir;
1445     // Bring the track beyond the material
1446     if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)){
1447       return 0;
1448     }
1449     break;
1450   }
1451
1452
1453   return 1;
1454 }
1455
1456
1457 //_____________________________________________________________________________
1458 Double_t AliITStrackerUpgrade::GetPredictedChi2MI(AliITStrackMI* track, const AliITSRecPoint *cluster,Int_t layer)
1459 {
1460   //
1461   // Compute predicted chi2
1462   //
1463   Float_t erry,errz,covyz;
1464   Float_t theta = track->GetTgl();
1465   Float_t phi   = track->GetSnp();
1466   phi = TMath::Abs(phi)*TMath::Sqrt(1./((1.-phi)*(1.+phi)));
1467   GetError(layer,cluster,theta,phi,track->GetExpQ(),erry,errz,covyz);
1468   AliDebug(2,Form(" chi2: tr-cl   %f  %f   tr X %f cl X %f",track->GetY()-cluster->GetY(),track->GetZ()-cluster->GetZ(),track->GetX(),cluster->GetX()));
1469   // Take into account the mis-alignment (bring track to cluster plane)
1470   Double_t xTrOrig=track->GetX();
1471   if (!track->Propagate(xTrOrig+cluster->GetX())) return 1000.;
1472   AliDebug(2,Form(" chi2: tr-cl   %f  %f   tr X %f cl X %f",track->GetY()-cluster->GetY(),track->GetZ()-cluster->GetZ(),track->GetX(),cluster->GetX()));
1473   Double_t chi2 = track->GetPredictedChi2MI(cluster->GetY(),cluster->GetZ(),erry,errz,covyz);
1474   // Bring the track back to detector plane in ideal geometry
1475   // [mis-alignment will be accounted for in UpdateMI()]
1476   if (!track->Propagate(xTrOrig)) return 1000.;
1477   Float_t ny=0.,nz=0.;
1478   //GetNTeor(layer,cluster,theta,phi,ny,nz);
1479   Double_t delta = cluster->GetNy()+cluster->GetNz()-nz-ny;
1480   if (delta>1){
1481     chi2+=0.5*TMath::Min(delta/2,2.);
1482     chi2+=2.*cluster->GetDeltaProbability();
1483   }
1484   //
1485   track->SetNy(layer,ny);
1486   track->SetNz(layer,nz);
1487   track->SetSigmaY(layer,erry);
1488   track->SetSigmaZ(layer, errz);
1489   track->SetSigmaYZ(layer,covyz);
1490   //track->fNormQ[layer] = cluster->GetQ()/TMath::Sqrt(1+theta*theta+phi*phi);
1491   track->SetNormQ(layer,cluster->GetQ()/TMath::Sqrt((1.+ track->GetTgl()*track->GetTgl())/((1.-track->GetSnp())*(1.+track->GetSnp()))));
1492   return chi2;
1493 }
1494 //________________________________________________________________
1495 Int_t AliITStrackerUpgrade::GetError(Int_t /*layer*/,
1496                                      const AliITSRecPoint *cl,
1497                                      Float_t /*tgl*/,Float_t /*tgphitr*/,Float_t /*expQ*/,
1498                                      Float_t &erry,Float_t &errz,Float_t &covyz,
1499                                      Bool_t /*addMisalErr*/)
1500 {
1501   //
1502   // Calculate cluster position error for the upgrade
1503   //
1504   Int_t retval=0;
1505   covyz=0.;
1506   retval = GetErrorOrigRecPoint(cl,erry,errz,covyz);
1507
1508
1509
1510   return retval;
1511                                                                                    
1512 }
1513 //____________________________________________________________________-
1514
1515 Int_t AliITStrackerUpgrade::GetErrorOrigRecPoint(const AliITSRecPoint*cl,
1516                                                  Float_t &erry,Float_t &errz,Float_t &covyz)
1517 {
1518   //
1519   // Calculate cluster position error (just take error from AliITSRecPoint)
1520   //
1521   erry   = TMath::Sqrt(cl->GetSigmaY2());
1522   errz   = TMath::Sqrt(cl->GetSigmaZ2());
1523   covyz  = cl->GetSigmaYZ();
1524   return 1;
1525 }
1526 //__________________________
1527 void AliITStrackerUpgrade::GetNTeor(Int_t layer,const AliITSRecPoint* /*cl*/,
1528                                     Float_t tgl,Float_t tgphitr,
1529                                     Float_t &ny,Float_t &nz)
1530 {
1531   //
1532   // Get "mean shape" (original parametrization from AliITStrackerMI)
1533   //
1534   tgl = TMath::Abs(tgl);
1535   tgphitr = TMath::Abs(tgphitr);
1536
1537   // SPD
1538   if (layer==0) {
1539     ny = 1.+tgphitr*3.2;
1540     nz = 1.+tgl*0.34;
1541     return;
1542   }
1543   if (layer==1) {
1544     ny = 1.+tgphitr*3.2;
1545     nz = 1.+tgl*0.28;
1546     return;
1547   }
1548   // SSD
1549   if (layer==4 || layer==5) {
1550     ny = 2.02+tgphitr*1.95;
1551     nz = 2.02+tgphitr*2.35;
1552     return;
1553   }
1554   // SDD
1555   ny  = 6.6-2.7*tgphitr;
1556   nz  = 2.8-3.11*tgphitr+0.45*tgl;
1557   return;
1558 }
1559 //_________________________________________________________________
1560 Int_t AliITStrackerUpgrade::PropagateBack(AliESDEvent *event) {
1561   //--------------------------------------------------------------------
1562   // This functions propagates reconstructed ITS tracks back
1563   // The clusters must be loaded !
1564   //--------------------------------------------------------------------
1565   fTrackingPhase="PropagateBack";
1566   Int_t nentr=event->GetNumberOfTracks();
1567   Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
1568   Int_t ntrk=0;
1569   for (Int_t i=0; i<nentr; i++) {
1570     AliESDtrack *esd=event->GetTrack(i);
1571
1572     if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
1573     if (esd->GetStatus()&AliESDtrack::kITSout) continue;
1574
1575     AliITStrackMI *t=0;
1576     try {
1577       t=new AliITStrackMI(*esd);
1578     } catch (const Char_t *msg) {
1579       //Warning("PropagateBack",msg);
1580       delete t;
1581       continue;
1582     }
1583     t->SetExpQ(TMath::Max(0.8*t->GetESDtrack()->GetTPCsignal(),30.));
1584
1585     ResetTrackToFollow(*t);
1586
1587     /*
1588     // propagate to vertex [SR, GSI 17.02.2003]
1589     // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
1590     if (CorrectForPipeMaterial(&fTrackToFollow,"inward")) {
1591     if (fTrackToFollow.PropagateToVertex(event->GetVertex()))
1592     fTrackToFollow.StartTimeIntegral();
1593     // from vertex to outside pipe
1594     CorrectForPipeMaterial(&fTrackToFollow,"outward");
1595     }*/
1596     // Start time integral and add distance from current position to vertex
1597     Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
1598     fTrackToFollow.GetXYZ(xyzTrk);
1599     Double_t dst2 = 0.;
1600     for (Int_t icoord=0; icoord<3; icoord++) {
1601       Double_t di = xyzTrk[icoord] - xyzVtx[icoord];
1602       dst2 += di*di;
1603     }
1604     fTrackToFollow.StartTimeIntegral();
1605     fTrackToFollow.AddTimeStep(TMath::Sqrt(dst2));
1606
1607     fTrackToFollow.ResetCovariance(10.); 
1608     fTrackToFollow.ResetClusters();
1609     //
1610     Int_t inx[6];
1611     for (Int_t k=0; k<6; k++) inx[k]=-1;
1612     Int_t nclusters = t->GetNumberOfClusters();
1613     for(Int_t ncl=0;ncl<nclusters;ncl++){
1614       Int_t index = t-> GetClIndex(ncl);
1615       Int_t lay = (index & 0xf0000000) >> 28;
1616       inx[lay]=index;
1617     }
1618     //   
1619     if (RefitAtBase(43.6,&fTrackToFollow,inx)) {
1620       //fTrackToFollow.SetLabel(t->GetLabel());//
1621       //fTrackToFollow.CookdEdx();
1622       //CookLabel(&fTrackToFollow,0.); //For comparison only
1623       fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
1624       UseClusters(&fTrackToFollow);
1625       ntrk++;
1626     }
1627     delete t;
1628   }
1629
1630   AliInfo(Form("Number of back propagated ITS tracks: %d out of %d ESD tracks",ntrk,nentr));
1631
1632   fTrackingPhase="Default";
1633
1634   return 0;
1635 }
1636 //_________________________________________________________________________
1637 AliCluster *AliITStrackerUpgrade::GetCluster(Int_t index) const {
1638   //--------------------------------------------------------------------
1639   //       Return pointer to a given cluster
1640   //--------------------------------------------------------------------
1641   Int_t l=(index & 0xf0000000) >> 28;
1642   AliInfo(Form("index %i  cluster index %i layer %i", index,index & 0x0fffffff,index & 0xf0000000));
1643   return fLayers[l]->GetCluster(index);
1644 }
1645 //______________________________________________________________________________
1646 Int_t AliITStrackerUpgrade::CorrectForPipeMaterial(AliITStrackMI *t, TString direction) {
1647 //-------------------------------------------------------------------
1648 // Propagate beyond beam pipe and correct for material
1649 // (material budget in different ways according to fUseTGeo value)
1650 // Add time if going outward (PropagateTo or PropagateToTGeo)
1651 //-------------------------------------------------------------------
1652
1653 // Define budget mode:
1654 // 0: material from AliITSRecoParam (hard coded)
1655 // 1: material from TGeo in one step (on the fly)
1656 // 2: material from lut
1657 // 3: material from TGeo in one step (same for all hypotheses)
1658 Int_t mode;
1659 switch(fUseTGeo) {
1660  case 0:
1661    mode=0;
1662    break;
1663  case 1:
1664    mode=1;
1665    break;
1666  case 2:
1667    mode=2;
1668    break;
1669  case 3:
1670    if(fTrackingPhase.Contains("Clusters2Tracks"))
1671      { mode=3; } else { mode=1; }
1672    break;
1673  case 4:
1674    if(fTrackingPhase.Contains("Clusters2Tracks"))
1675      { mode=3; } else { mode=2; }
1676    break;
1677  default:
1678    mode=0;
1679    break;
1680  }
1681 if(fTrackingPhase.Contains("Default")) mode=0;
1682 Int_t index=fCurrentEsdTrack;
1683
1684 Float_t  dir = (direction.Contains("inward") ? 1. : -1.);
1685 Double_t rToGo=(dir>0 ? AliITSRecoParam::GetrInsidePipe() : AliITSRecoParam::GetrOutsidePipe());
1686 Double_t xToGo;
1687 if (!t->GetLocalXat(rToGo,xToGo)) return 0;
1688
1689 Double_t xOverX0,x0,lengthTimesMeanDensity;
1690
1691 switch(mode) {
1692  case 1:
1693    xOverX0 = AliITSRecoParam::GetdPipe();
1694    x0 = AliITSRecoParam::GetX0Be();
1695    lengthTimesMeanDensity = xOverX0*x0;
1696    lengthTimesMeanDensity *= dir;
1697    if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)) return 0;
1698    break;
1699  case 0:
1700    if (!t->PropagateToTGeo(xToGo,1)) return 0;
1701    break;
1702  case 2:
1703    if(fxOverX0Pipe<0) BuildMaterialLUT("Pipe");
1704    xOverX0 = fxOverX0Pipe;
1705    lengthTimesMeanDensity = fxTimesRhoPipe;
1706    lengthTimesMeanDensity *= dir;
1707    if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)) return 0;
1708    break;
1709  case 3:
1710    if(!fxOverX0PipeTrks || index<0 || index>=fNtracks) Error("CorrectForPipeMaterial","Incorrect usage of UseTGeo option!\n");
1711    if(fxOverX0PipeTrks[index]<0) {
1712      if (!t->PropagateToTGeo(xToGo,1,xOverX0,lengthTimesMeanDensity)) return 0;
1713      Double_t angle=TMath::Sqrt((1.+t->GetTgl()*t->GetTgl())/
1714                                 ((1.-t->GetSnp())*(1.+t->GetSnp())));
1715      fxOverX0PipeTrks[index] = TMath::Abs(xOverX0)/angle;
1716      fxTimesRhoPipeTrks[index] = TMath::Abs(lengthTimesMeanDensity)/angle;
1717      return 1;
1718    }
1719    xOverX0 = fxOverX0PipeTrks[index];
1720    lengthTimesMeanDensity = fxTimesRhoPipeTrks[index];
1721    lengthTimesMeanDensity *= dir;
1722    if (!t->PropagateTo(xToGo,xOverX0,lengthTimesMeanDensity/xOverX0)) return 0;
1723    break;
1724  }
1725
1726 return 1;
1727 }
1728 //__________________________________________________________________________
1729 Int_t AliITStrackerUpgrade::RefitInward(AliESDEvent *event) {
1730   //--------------------------------------------------------------------
1731   // This functions refits ITS tracks using the
1732   // "inward propagated" TPC tracks
1733   // The clusters must be loaded !
1734   //--------------------------------------------------------------------
1735   fTrackingPhase="RefitInward";
1736   Int_t nentr=event->GetNumberOfTracks();
1737   for (Int_t i=0; i<nentr; i++) {
1738     AliESDtrack *esd=event->GetTrack(i);
1739
1740     // if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
1741     // if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
1742     //if (esd->GetStatus()&AliESDtrack::kTPCout)
1743     //  if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
1744     AliITStrackMI *t=0;
1745     try {
1746       t=new AliITStrackMI(*esd);
1747     } catch (const Char_t *msg) {
1748       //Warning("RefitInward",msg);
1749       delete t;
1750       continue;
1751     }
1752
1753     ResetTrackToFollow(*t);
1754     // fTrackToFollow.ResetClusters();
1755     //  if ((esd->GetStatus()&AliESDtrack::kTPCin)==0)
1756     //  fTrackToFollow.ResetCovariance(10.);
1757     //Refitting...
1758     //The beam pipe
1759     // if (CorrectForPipeMaterial(&fTrackToFollow,"inward")) {
1760     fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
1761     AliESDtrack  *esdTrack =fTrackToFollow.GetESDtrack();
1762     Double_t r[3]={0.,0.,0.};
1763     Double_t maxD=3.;
1764     esdTrack->RelateToVertex(event->GetVertex(),GetBz(r),maxD);
1765     // }
1766     delete t;
1767   }
1768
1769   fTrackingPhase="Default";
1770   return 0;
1771 }
1772
1773
1774