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