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