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