]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUTrackerSA.cxx
update of the tracker
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUTrackerSA.cxx
1 //-------------------------------------------------------------------------
2 //               Implementation of the ITS tracker class
3 //    It reads AliITSUClusterPix clusters and and fills the ESD with tracks
4 //-------------------------------------------------------------------------
5
6 #include <TBranch.h>
7 #include <TMath.h>
8 using TMath::Abs;
9 using TMath::Sort;
10 using TMath::Sqrt;
11 #include <TTree.h>
12 #include <algorithm>
13 using std::sort;
14
15 // Vc library
16 //#include "Vc/Vc"
17 //#include "AliITSUTrackerSAauxVc.h" // Structs and other stuff using Vc library  
18 #include "AliLog.h"
19 #include "AliESDEvent.h"
20 #include "AliITSUClusterPix.h"
21 #include "AliITSUTrackerSA.h"
22 #include "AliITSUReconstructor.h"
23 #include "AliITSURecoDet.h"
24 #include "AliESDtrack.h"
25
26 #include <Riostream.h>
27
28
29
30 //#include "AliITSUtrackSA.h"      // Some dedicated SA track class ?  
31
32 ClassImp(AliITSUTrackerSA)
33
34 const Double_t AliITSUTrackerSA::fgkToler =  1e-6;// tolerance for layer on-surface check
35
36 //________________________________________________________________________________
37 AliITSUTrackerSA::AliITSUTrackerSA(AliITSUReconstructor* rec) :  
38   fReconstructor(rec)
39   ,fITS(0)
40   ,fMatLUT(0)
41   ,fUseMatLUT(kFALSE)
42   ,fCurrMass(0.14)
43   //
44   ,fClusters(),
45   fClustersTC(),
46   fDoublets(),
47   fIndex(),
48   fNClusters(),
49   fNDoublets(),
50   fPhiCut(0.05),
51   fRPhiCut(0.03),
52   fZCut(0.01)
53 {
54   //--------------------------------------------------------------------
55   // This default constructor needs to be provided
56   //--------------------------------------------------------------------
57   for(Int_t i=0;i<7;++i) {
58     fClusters[i].reserve(5000);
59   }
60   if (rec) Init(rec);
61 }
62
63 //________________________________________________________________________________
64 AliITSUTrackerSA::AliITSUTrackerSA(const AliITSUTrackerSA &t): 
65   AliTracker(t),
66   fReconstructor(t.fReconstructor),
67   fITS(t.fITS),
68   fMatLUT(t.fMatLUT),
69   fUseMatLUT(t.fUseMatLUT),
70   fCurrMass(t.fCurrMass),
71   //
72   fClusters(),
73   fClustersTC(),
74   fIndex(),
75   fNClusters(),
76   fNDoublets(),
77   fPhiCut(),
78   fRPhiCut(),
79   fZCut()
80 {
81   //--------------------------------------------------------------------
82   // The copy constructor is protected
83   //--------------------------------------------------------------------
84 }
85
86 //________________________________________________________________________________
87 AliITSUTrackerSA::~AliITSUTrackerSA() 
88 {
89   // d-tor
90   delete fMatLUT;
91 }
92
93
94 //_________________________________________________________________________
95 void AliITSUTrackerSA::Init(AliITSUReconstructor* rec)
96 {
97   // init with external reconstructor
98   //
99   fITS = rec->GetITSInterface();
100   //
101   // create material lookup table
102   const int kNTest = 1000;
103   const double kStepsPerCM=5;
104   fMatLUT  = new AliITSUMatLUT(fITS->GetRMin(),fITS->GetRMax(),Nint(kStepsPerCM*(fITS->GetRMax()-fITS->GetRMin())));
105   double zmn = 1e6;
106   for (int ilr=fITS->GetNLayers();ilr--;) {
107     AliITSURecoLayer* lr = fITS->GetLayer(ilr);
108     if (zmn>Abs(lr->GetZMin())) zmn = Abs(lr->GetZMin());
109     if (zmn>Abs(lr->GetZMax())) zmn = Abs(lr->GetZMax());
110   }
111   fMatLUT->FillData(kNTest,-zmn,zmn);
112   //
113 }
114
115 //________________________________________________________________________________
116 Int_t AliITSUTrackerSA::Clusters2Tracks(AliESDEvent *event) {
117   //--------------------------------------------------------------------
118   // This is the main tracking function
119   // The clusters must already be loaded
120   //--------------------------------------------------------------------
121
122   // Possibly, create the track "seeds" (combinatorial)
123
124   // Possibly, increment the seeds with additional clusters (Kalman)
125
126   // Possibly, (re)fit the found tracks 
127
128   // Three iterations: 
129   // - High momentum first;
130   // - Low momentum with vertex constraint; 
131   // - Everything else. 
132
133   MakeDoublets();       // To be implemented
134   //MakeTriplets();       // Are triplets really necessary? MFT does not use them. 
135   CASelection(event);       
136   //GlobalFit();          // To be implemented
137   //ChiSquareSelection(); // To be implemented
138
139   return 0;
140 }
141
142 //________________________________________________________________________________
143 Int_t AliITSUTrackerSA::PropagateBack(AliESDEvent */*event*/) {
144   //--------------------------------------------------------------------
145   // Here, we implement the Kalman smoother ?
146   // The clusters must be already loaded
147   //--------------------------------------------------------------------
148
149   return 0;
150 }
151
152 //________________________________________________________________________________
153 Int_t AliITSUTrackerSA::RefitInward(AliESDEvent */*event*/) {
154   //--------------------------------------------------------------------
155   // Some final refit, after the outliers get removed by the smoother ?  
156   // The clusters must be loaded
157   //--------------------------------------------------------------------
158
159   return 0;
160 }
161
162 Int_t AliITSUTrackerSA::LoadClusters(TTree *cluTree) {
163   //--------------------------------------------------------------------
164   // This function reads the ITSU clusters from the tree,
165   // sort them, distribute over the internal tracker arrays, etc
166   //--------------------------------------------------------------------
167   fITS->LoadClusters(cluTree);
168   fITS->ProcessClusters();
169   //
170   for(int iL=0; iL<7; ++iL) {
171     fClustersTC[iL]=*fITS->GetLayerActive(iL)->GetClustersAddress();
172     TClonesArray *clCont=fClustersTC[iL];
173     fNClusters[iL]=clCont->GetEntriesFast();
174     Float_t phi[fNClusters[iL]];
175     fIndex[iL] = new Int_t[fNClusters[iL]];
176     for(int iC=0;iC<fNClusters[iL];++iC) {
177       const AliITSUClusterPix *cl = (AliITSUClusterPix*)clCont->At(iC);
178       float pos[3];
179       cl->GetGlobalXYZ(pos);
180       phi[iC] = pos[0]==0.f ? TMath::PiOver2() : TMath::ATan2(pos[1]-GetY(),pos[0]-GetX());
181       float angle=0.f; // TO BE UNDERSTOOD: DO I STILL NEED THE STATION ANGLE IF I USE THE GLOBAL COVARIANCE MATRIX?
182       fClusters[iL].push_back(itsCluster(pos[0],pos[1],pos[2],cl->GetSigmaY2(),cl->GetSigmaZ2(),cl->GetSigmaYZ(),phi[iC],angle));
183     }
184     TMath::Sort(fNClusters[iL],phi,fIndex[iL],kFALSE);
185   }
186 // #ifdef __DEBUG__
187 //   PrintInfo("clusters");
188 // #endif
189   return 0;
190 }
191
192 //________________________________________________________________________________
193 void AliITSUTrackerSA::UnloadClusters() {
194   //--------------------------------------------------------------------
195   // This function unloads ITSU clusters from the RAM
196   //--------------------------------------------------------------------
197   for(int i=0;i<7;++i) {
198     fClusters[i].clear();
199     fNClusters[i]=0;
200     delete fIndex[i];
201   }
202   for(int i=0;i<6;++i) fDoublets[i].clear();
203 }
204
205 //________________________________________________________________________________
206 AliCluster *AliITSUTrackerSA::GetCluster(Int_t /*index*/) const {
207   //--------------------------------------------------------------------
208   //       Return pointer to a given cluster
209   //--------------------------------------------------------------------
210   return 0;  // replace with an actual pointer 
211 }
212
213 //________________________________________________________________________________
214 void AliITSUTrackerSA::CASelection(AliESDEvent *event) {
215   // Here it's implemented the Cellular Automaton routine
216   // Firstly the level of each doublet is set according to the level of 
217   // the neighbour doublets. 
218   // Doublet are considered to be neighbour if they share one point and the 
219   // phi and theta direction difference of the two is below a cut value.
220
221   cout << "Begin of the CA selection" << endl;
222   for( int iL = 1; iL < 6; ++iL ) {
223     
224     const itsCluster* clusters1 = &fClusters[iL-1][0];
225     const itsCluster* clusters2 = &fClusters[iL][0];
226     const itsCluster* clusters3 = &fClusters[iL+1][0];    
227
228     const nPlets* doublets1 = &fDoublets[iL-1][0];
229     nPlets* doublets2 = &fDoublets[iL][0];
230
231     for ( int iD2 = 0; iD2 < fNDoublets[iL]; ++iD2 ) {
232       for ( int iD1 = 0; iD1 < fNDoublets[iL-1]; ++iD1 ) {
233         const int id1 = doublets1[iD1].id0;
234         const int id2 = doublets2[iD2].id0;
235         if ( id1 == id2 ) {
236           if ( doublets2[iD2].level <= ( doublets1[iD1].level + 1 ) ) {
237             const int id3 = doublets2[iD2].id1;
238             const float r3 = Sqrt( clusters3[id3].x * clusters3[id3].x + clusters3[id3].y * clusters3[id3].y );
239             const float r2 = Sqrt( clusters2[id2].x * clusters2[id2].x + clusters2[id2].y * clusters2[id2].y );
240             const float extrZ3 = doublets1[iD1].tanLambda * ( r3 - r2 ) + clusters2[id2].z ;
241             //cout << extrZ3 << " " << clusters3[id3].z << " " << Abs ( extrZ3 - clusters3[id3].z ) << endl;
242             if ( Abs ( extrZ3 - clusters3[id3].z ) < fZCut ) {
243               //cout << "OK Z doublets: "<< iL-1 << "," << iD1 << "\t" << iL << "," <<iD2 << endl;
244               const float det = (clusters1[id1].x - GetX())*(clusters2[id2].y - GetY()) - (clusters2[id2].x-GetX() )*(clusters1[id1].y - GetY()); // (GetX() - clusters2[id2].x)*(clusters1[id1].y - clusters2[id2].y) - (GetY() - clusters2[id2].y)*(clusters1[id1].x - clusters2[id2].x);
245               //cout << det << endl;
246               if ( Abs(det) <= 1e-12 ) {
247                 // linear extrapolation to next layer
248                 const float dsq = ( doublets1[iD1].tanPhi * (clusters3[id3].x + clusters2[id2].x) + clusters3[id3].y - clusters2[id2].y ) * 
249                   ( doublets1[iD1].tanPhi * (clusters3[id3].x + clusters2[id2].x) + clusters3[id3].y - clusters2[id2].y ) / (1 + doublets1[iD1].tanPhi * doublets1[iD1].tanPhi );
250                 if ( dsq < fRPhiCut*fRPhiCut )  {
251                   doublets2[iD2].level += doublets1[iD1].level;
252                   doublets2[iD2].neighbours.push_back(iD1);
253                 }
254               } else {
255                 const float r1sq = clusters1[id1].x * clusters1[id1].x + clusters1[id1].y * clusters1[id1].y ;
256                 const float rvsq = GetX() * GetX() + GetY() * GetY();
257                 const float deta = (rvsq - r1sq) * (clusters2[id2].y - GetY()) - (rvsq - r2*r2) * (clusters1[id1].y - GetY());
258                 const float detb = - (rvsq - r1sq) * (clusters2[id2].x - GetX()) + (rvsq - r2*r2) * (clusters1[id1].x - GetX()) ;
259                 const float a = deta/det ;
260                 const float b = detb/det ;
261                 const float c = -rvsq - a * GetX() - b * GetY();
262                 const float rc = Sqrt( a*a/4.f + b*b/4.f - c );
263                 const float d = Sqrt( (a/2.f + clusters3[id3].x) * (a/2.f + clusters3[id3].x) + (b/2.f + clusters3[id3].y) * (b/2.f + clusters3[id3].y) );
264                 //cout << d << " " << rc << " " << d - rc << endl;
265                 if ( Abs( d - rc ) < fRPhiCut ) {
266                   doublets2[iD2].level += doublets1[iD1].level;
267                   doublets2[iD2].neighbours.push_back(iD1);
268                 }
269               }
270             }
271           }
272         }
273       }
274     }
275   }
276
277   /*#ifdef __DEBUG__
278   PrintInfo("doublets");
279   #endif*/
280
281   // Hic sunt leones: the following code could be optimised to be iterative. But now I don't have time.
282   for ( int level = 6; level >= 2 ; --level ) {
283     vector<trackC> candidates;
284     int nCand=0;
285     const int limit = level > 5 ? 5 : level;
286     for ( int doubl = 5; doubl >= limit; --doubl ) {
287       for ( int iD = 0; iD < fNDoublets[doubl]; ++iD ) {
288         if ( fDoublets[doubl][iD].level == level ) {
289           cout << "Pushing new candidate" << endl;
290           candidates.push_back(trackC());
291           candidates[nCand].fPoints[doubl*2+2] = fDoublets[doubl][iD].id1;
292           candidates[nCand].fPoints[doubl*2] = fDoublets[doubl][iD].id0;
293           CandidatesTreeTraversal(candidates, iD, doubl, nCand);
294         }
295       }
296     }
297
298     nCand = candidates.size();
299     Double_t rDest = fITS->GetRMax();
300     cout << "Number of candidates at level " << level << ": " << candidates.size() << endl;
301     int index[nCand];
302     for( int i = 0; i < nCand; ++i ) index[i] = i;
303     for ( int cand = 0; cand < nCand; ++cand ) {
304       int count = 0;
305       for ( int j = 0; j < 14; ++j ) {
306         if( candidates[cand].fPoints[j] > -1 ) ++count; 
307       }
308       cout << "Candidates " << cand << ", number of points: " << count << endl;
309       //      candidates[cand].ResetCovariance(50.);
310       InitTrackParams(candidates[cand]);
311       cout << "Fit cnd: " << cand << " " << candidates[cand] << endl;
312       candidates[cand].fChi2 = RefitTrack( (AliExternalTrackParam*)&candidates[cand], candidates[cand].fPoints, rDest ,-1);
313     }
314     CompDesc comp(&candidates);
315     sort(index,index+nCand,comp);
316     for ( int cand = 0; cand < nCand; ++cand ) {
317       const int ii = index[cand];
318 #ifdef __DEBUG__
319       //cout << ii << " " << candidates[ii] << endl;
320 #endif
321       if ( candidates[ii].fChi2 < 0. ) break;
322       bool goodTrack = true;
323       for ( unsigned int point = 0; point < 14; ++point ) {
324         if ( candidates[ii].fPoints[point] != -1 ) {
325           if( !(fClusters[ point/2 ][ candidates[ii].fPoints[point] ].isUsed ) ) {
326             fClusters[ point/2 ][ candidates[ii].fPoints[point] ].isUsed = true;
327           } else {
328             goodTrack = false;
329           }
330         }
331       }
332       if ( goodTrack ) {
333         AliESDtrack outTrack;
334         outTrack.SetOuterParam((AliExternalTrackParam*)&candidates[ii],AliESDtrack::kITSpureSA);
335         event->AddTrack(&outTrack);
336       }
337     }
338
339   }
340
341   // A first try... 
342   // for ( int level = 6; level >=3 ; --level ) {
343   //   vector<int> points[7];
344   //   for ( int doubl = 5; doubl >= 2; --doubl ) {
345   //     if( ( doubl + 1 - level ) < 0 ) break;
346   //     for ( int iD = 0; iD < fNDoublets[doubl]; ++iD ) {
347   //    if ( fDoublets[doubl][iD].level == level ) {
348   //      points[doubl+1].push_back( fDoublets[doubl][iD].id1 );
349   //      points[doubl].push_back( fDoublets[doubl][iD].id0 );
350      
351   //      vector<int> currentVector = (fDoublets[doubl][iD].neighbours);
352   //      for( int iL = 1; iL <= doubl; ++iL ) {
353   //        const nPlets * doublets = &fDoublets[doubl-iL][0];
354   //        vector<int> tmp;
355   //        for ( unsigned int iV = 0 ; iV < currentVector.size() ; ++iV ) {
356   //          points[doubl-iL].push_back( doublets[ currentVector[ iV ] ].id0 ) ;
357   //          for ( unsigned int iN = 0 ; iN < doublets[ currentVector[ iV ] ].neighbours.size(); ++iN ) {
358   //            tmp.push_back( doublets[ currentVector[ iV ] ].neighbours[ iN ] );
359   //          }
360   //        }
361   //        currentVector.swap( tmp );
362   //      }
363   //    }
364   //     }
365   //   }
366   // }
367 }
368
369 //________________________________________________________________________________
370 void AliITSUTrackerSA::MakeDoublets() {
371   // Make associations between two points on adjacent layers within an azimuthal window. 
372   // Under consideration:
373   // - track parameter estimation using the primary vertex position
374   // To do:
375   // - last iteration
376
377   cout << "Vertex of used by the tracker: " << GetX() << " " << GetY() << " " << GetZ() << endl;
378
379   for( int iL = 0 ; iL < 6 ; ++iL ) {
380     fNDoublets[iL] = 0; 
381     const itsCluster* clusters1 = &fClusters[iL][0];
382     const itsCluster* clusters2 = &fClusters[iL+1][0];
383
384     // 0 - 2Pi junction treatment (part I)
385     for ( int iCC1 = 0 ; iCC1 < fNClusters[iL] ; ++iCC1 ) {
386       bool flag = true;
387       const int iC1 = fIndex[iL][iCC1];
388       for ( int iCC2 = fNClusters[iL+1]-1; iCC2 > -1 ; --iCC2 ) {
389         const int iC2 = fIndex[iL+1][iCC2];
390         if( (TMath::TwoPi() - (clusters2[iC2].phi-clusters1[iC1].phi) ) < fPhiCut ) {
391           fDoublets[iL].push_back(nPlets(iC1,iC2));
392           fDoublets[iL][fNDoublets[iL]].tanPhi = (clusters1[iC1].y-clusters2[iC2].y)/(clusters1[iC1].x-clusters2[iC2].x);
393           float r1  = Sqrt(clusters1[iC1].x * clusters1[iC1].x + clusters1[iC1].y * clusters1[iC1].y);
394           cout << clusters2[iC2].x * clusters2[iC2].x + clusters2[iC2].y * clusters2[iC2].y << flush << endl;
395           float r2  = Sqrt(clusters2[iC2].x * clusters2[iC2].x + clusters2[iC2].y * clusters2[iC2].y);
396           fDoublets[iL][fNDoublets[iL]].tanLambda = (clusters1[iC1].z-clusters2[iC2].z)/(r1-r2);
397           ++fNDoublets[iL];
398           flag = false;
399         } else break;
400
401       }
402       if (flag) break;
403     }
404
405     
406     // "Central" points 
407     for ( int iCC1 = 0 ; iCC1 < fNClusters[iL] ; ++iCC1 ) {
408       const int iC1 = fIndex[iL][iCC1];
409       for ( int iCC2 = 0; iCC2 < fNClusters[iL+1] ; ++iCC2 ) {
410         const int iC2 = fIndex[iL+1][iCC2];
411         if( Abs( clusters1[iC1].phi - clusters2[iC2].phi ) < fPhiCut ) {
412           fDoublets[iL].push_back(nPlets(iC1,iC2));
413           fDoublets[iL][fNDoublets[iL]].tanPhi = (clusters1[iC1].y-clusters2[iC2].y)/(clusters1[iC1].x-clusters2[iC2].x);
414           float r1  = Sqrt(clusters1[iC1].x * clusters1[iC1].x + clusters1[iC1].y * clusters1[iC1].y);
415           float r2  = Sqrt(clusters2[iC2].x * clusters2[iC2].x + clusters2[iC2].y * clusters2[iC2].y);
416           fDoublets[iL][fNDoublets[iL]].tanLambda = (clusters1[iC1].z-clusters2[iC2].z)/(r1-r2);
417           ++fNDoublets[iL];
418         } else if( clusters2[iC2].phi - clusters1[iC1].phi > fPhiCut ) break;
419       
420       }
421
422     }
423
424     // 0 - 2Pi junction treatment (part II)
425     for ( int iCC1 = fNClusters[iL]-1; iCC1 > -1 ; --iCC1 ) {
426       bool flag = true;
427       const int iC1 = fIndex[iL][iCC1];
428       for ( int iCC2 = 0; iCC2 < fNClusters[iL+1] ; ++iCC2 ) {
429         const int iC2 = fIndex[iL+1][iCC2];
430         if( (TMath::TwoPi() - (clusters1[iC1].phi-clusters2[iC2].phi) ) < fPhiCut ) { 
431           fDoublets[iL].push_back(nPlets(iC1,iC2));
432           fDoublets[iL][fNDoublets[iL]].tanPhi = (clusters1[iC1].y-clusters2[iC2].y)/(clusters1[iC1].x-clusters2[iC2].x);
433           float r1  = Sqrt(clusters1[iC1].x * clusters1[iC1].x + clusters1[iC1].y * clusters1[iC1].y);
434           float r2  = Sqrt(clusters2[iC2].x * clusters2[iC2].x + clusters2[iC2].y * clusters2[iC2].y);
435           fDoublets[iL][fNDoublets[iL]].tanLambda = (clusters1[iC1].z-clusters2[iC2].z)/(r1-r2);
436           ++fNDoublets[iL];
437           flag = false;
438         } else break;
439
440       }
441
442       if (flag) break;
443     }
444   
445   }
446   // #ifdef __DEBUG__
447   // PrintInfo("doublets");
448   // #endif
449 }
450 //______________________________________________________________________________
451 Bool_t AliITSUTrackerSA::InitTrackParams(trackC &track)
452 {
453   // Set the initial guess on track kinematics for propagation.
454   // Assume at least 3 points available
455   int lrOcc[AliITSUAux::kMaxLayers], nCl=0;
456   //
457   // we will need endpoints and middle layer
458   for (int i=fITS->GetNLayersActive();i--;) if (track.fPoints[i<<0x1]>-1) lrOcc[nCl++] = i;
459   if (nCl<3) {
460     AliError(Form("Cannot estimate momentum of tracks with %d clusters",nCl));
461     cout << track << endl;
462     return kFALSE;
463   }
464   //
465   int lr0   = lrOcc[0];
466   int lr1   = lrOcc[nCl/2];
467   int lr2   = lrOcc[nCl-1];
468   //
469   const itsCluster& cl0 = fClusters[lr0][ track.fPoints[lr0<<0x1] ];
470   const itsCluster& cl1 = fClusters[lr1][ track.fPoints[lr1<<0x1] ];
471   const itsCluster& cl2 = fClusters[lr2][ track.fPoints[lr2<<0x1] ];
472   double cv = Curvature(cl0.x,cl0.y, cl1.x,cl1.y, cl2.x,cl2.y);
473   double tgl = (cl2.z-cl0.z)/TMath::Sqrt((cl2.x-cl0.x)*(cl2.x-cl0.x)+(cl2.y-cl0.y)*(cl2.y-cl0.y));
474   //  double phi = TMath::ATan2((cl2.y-cl1.y),(cl2.x-cl1.x));
475   //
476   AliITSUClusterPix* clus = (AliITSUClusterPix*)fClustersTC[ lr0 ]->At( track.fPoints[lr0<<0x1] );
477   AliITSURecoLayer* lr = fITS->GetLayerActive(lr0);
478   AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
479   double x = sens->GetXTF() + clus->GetX();
480   double alp = sens->GetPhiTF();
481   //  printf("Alp: %f phi: %f\n",alp,phi);
482   double par[5] = {clus->GetY(),clus->GetZ(),0,tgl,cv};
483   double cov[15] = {
484     5*5,
485     0, 5*5,
486     0, 0, 0.7*0.7,
487     0,0,0,0.7*0.7,
488     0,0,0,0,10
489   };
490   track.Set(x,alp,par,cov);
491   return kTRUE;
492 }
493
494 //______________________________________________________________________________
495 void AliITSUTrackerSA::CandidatesTreeTraversal(vector<trackC> &track,  int &iD, int &doubl, int &nCand) {
496   
497   if ( doubl < 1 ) {
498 #ifdef __DEBUG__
499     cout << "ERROR IN CandidatesTreeTraversal" << endl;
500 #endif
501     return;
502   }
503
504   int currD = doubl-1;
505   track[nCand].fPoints[ 2 * currD ] = fDoublets[ currD ] [ fDoublets[ doubl ][ iD ].neighbours[0] ].id0;
506
507
508   if ( fDoublets[ currD ] [ fDoublets[ doubl ][ iD ].neighbours[0] ].level == 1 ) {
509
510     cout << "Setting track " << nCand << " information using its first cluster: " << endl;
511     cout << "Layer " << currD << ", id " << track[nCand].fPoints[ 2 * currD ] <<endl;
512     AliITSUClusterPix* clus = (AliITSUClusterPix*)fClustersTC[ currD ]->At( track[nCand].fPoints[ 2 * currD ] ); // assign your 1st cluster (from which the fit is starting)
513     AliITSURecoLayer* lr = fITS->GetLayerActive(currD) ; // assign the layer which the cluster belongs to
514     
515     AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
516     double x = sens->GetXTF() + clus->GetX();
517     double alp = sens->GetPhiTF();
518     double par[5] = {clus->GetY(),clus->GetZ(),0,0,0.01};
519     double cov[15] = {
520       5*5,
521       0, 5*5,
522       0, 0, 0.7*0.7,
523       0,0,0,0.7*0.7,
524       0,0,0,0,10
525     };
526     track[nCand].Set(x,alp,par,cov);
527     //    cout << "after set " << track[nCand] << endl;
528     return;
529   } else {
530     CandidatesTreeTraversal(track,fDoublets[doubl][iD].neighbours[0],currD,nCand);
531   }
532   
533   for( unsigned int iD2 = 1 ; iD2 < fDoublets[ doubl ][ iD ].neighbours.size() ; ++iD2 ) {
534     track.push_back(track[nCand++]);
535     track[nCand].fPoints[ 2 * currD ] = fDoublets[ currD ] [ fDoublets[ doubl ][ iD ].neighbours[ iD2 ] ].id0;
536     if ( fDoublets[ currD ] [ fDoublets[ doubl ][ iD ].neighbours[iD2] ].level == 1 ) {
537       AliITSUClusterPix* clus = (AliITSUClusterPix*)fClustersTC[ currD ]->At( track[nCand].fPoints[ 2 * currD ] ); // assign your 1st cluster (from which the fit is starting)
538       AliITSURecoLayer* lr = fITS->GetLayerActive(currD) ; // assign the layer which the cluster belongs to
539       
540       AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
541       double x = sens->GetXTF() + clus->GetX();
542       double alp = sens->GetPhiTF();
543       double par[5] = {clus->GetY(),clus->GetZ(),0,0,0};
544       double cov[15] = {
545         5*5,
546         0, 5*5,
547         0, 0, 0.7*0.7,
548         0,0,0,0.7*0.7,
549         0,0,0,0,10
550       };
551       track[nCand].Set(x,alp,par,cov);
552       return;
553     } else {
554       CandidatesTreeTraversal(track,fDoublets[doubl][iD].neighbours[iD2],currD,nCand);
555     }
556   }
557 }
558
559 //______________________________________________________________________________
560 Double_t AliITSUTrackerSA::RefitTrack(AliExternalTrackParam* trc, 
561                                       Int_t clInfo[2*AliITSUAux::kMaxLayers],
562                                       Double_t rDest, Int_t stopCond)
563 {
564   // refit track till radius rDest. 
565   // if stopCond<0 : propagate till last cluster then stop
566   // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
567   // if stopCond>0 : rDest must be reached
568   //
569   // The clList should provide the indices of clusters at corresponding layer (as stored in the layer
570   // TClonesArray, with convention (allowing for up to 2 clusters per layer due to the overlaps):
571   // if there is a cluster on given layer I, then it should be stored at clInfo[2*I-1]
572   // if there is an additional cluster on this layer, it goes to clInfo[2*I]
573   // -1 means no cluster
574   //
575   double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
576   int dir,lrStart,lrStop;
577   //
578   dir = rCurr<rDest ? 1 : -1;
579   lrStart = fITS->FindFirstLayerID(rCurr,dir);
580   lrStop  = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
581   cout << "Start/End layers and direction " << lrStart << " " << lrStop << " " <<  dir << endl;
582   //
583   if (lrStop<0 || lrStart<0) AliFatal(Form("Failed to find start(%d) or last(%d) layers. "
584                                            "Track from %.3f to %.3f",lrStart,lrStop,rCurr,rDest));
585   //
586   int nCl = 0;
587   for (int i=2*fITS->GetNLayersActive();i--;) {if (clInfo[i]<0) continue; nCl++;}
588   //
589   AliExternalTrackParam tmpTr(*trc);
590   double chi2 = 0;
591   int iclLr[2],nclLr;
592   int nclFit = 0;
593   //
594   int lrStop1 = lrStop+dir;
595   for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
596     AliITSURecoLayer* lr = fITS->GetLayer(ilr);
597     if ( dir*(rCurr-lr->GetR(dir))>0) {
598       cout << ilr << " passed!" << endl;
599       continue;
600     } // this layer is already passed
601     int ilrA2,ilrA = lr->GetActiveID();
602     // passive layer or active w/o hits will be traversed on the way to next cluster
603     if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) {
604       cout << ilr << " is inactive or without cluster for current candidates" << endl;
605       continue;
606     } 
607     cout << "OK layer " << ilr << endl;
608     //
609     // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
610     nclLr=0;
611     if (dir>0) { // clusters are stored in increasing radius order
612       iclLr[nclLr++]=clInfo[ilrA2++];
613       if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
614     }
615     else {
616       if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
617       iclLr[nclLr++]=clInfo[ilrA2];
618     }
619     //
620     Bool_t transportedToLayer = kFALSE;
621     for (int icl=0;icl<nclLr;icl++) {
622       AliITSUClusterPix* clus =  (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
623       AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
624       if (!tmpTr.Rotate(sens->GetPhiTF())) { cout << "failed rotation" << endl; return -1; }
625       //
626       double xClus = sens->GetXTF()+clus->GetX();
627       if (!transportedToLayer) {
628         if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
629           cout << "failed transport to the entrance" << endl; return -1; } // go to the entrance to the layer
630         lrStart = ilr;
631         transportedToLayer = kTRUE;
632       }
633       //
634       if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) { cout << "failed propagation of the seed X:" << xClus << endl; tmpTr.Print(); return -1; }
635       //
636       Double_t p[2]={clus->GetY(), clus->GetZ()};
637       Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
638       double chi2cl = tmpTr.GetPredictedChi2(p,cov);
639       chi2 += chi2cl;
640       //
641       if ( !tmpTr.Update(p,cov) ) { cout << "failed update of the covariance" << endl; return -1; }
642       if (++nclFit==nCl && stopCond<0) {
643         *trc = tmpTr;
644         //      printf("Fit chi2: %f for %d clusters\n",chi2,nclFit);
645         return chi2; // it was requested to not propagate after last update
646       }
647     }
648     //
649   }
650   // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
651   // Still, try to go as close as possible to rDest.
652   //
653   //  printf("Fit chi2: %f for %d clusters\n",chi2,nclFit);
654   //
655   if (lrStart!=lrStop) {
656     if (!TransportToLayer(&tmpTr,lrStart,lrStop)) return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
657     if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
658   }
659   // go to the destination radius. Note that here we don't select direction to avoid precision problems
660   if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
661     return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
662   }
663   *trc = tmpTr;
664
665   return chi2;
666 }
667
668 //______________________________________________________________________________
669 Bool_t AliITSUTrackerSA::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
670 {
671   // propagate seed to given x applying material correction if requested
672   const Double_t kEpsilon = 1e-5;
673   Double_t xpos     = seed->GetX();
674   Int_t dir         = (xpos<xToGo) ? 1:-1;
675   Double_t xyz0[3],xyz1[3];
676   //
677   Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
678   if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
679   while ( (xToGo-xpos)*dir > kEpsilon){
680     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
681     Double_t x    = xpos+step;
682     Double_t bz=GetBz();   // getting the local Bz
683     if (!seed->PropagateTo(x,bz)) return kFALSE;
684     double ds = 0;
685     if (matCorr || updTime) {
686       xyz0[0]=xyz1[0]; // global pos at the beginning of step
687       xyz0[1]=xyz1[1];
688       xyz0[2]=xyz1[2];
689       seed->GetXYZ(xyz1);    //  // global pos at the end of step
690       //
691       if (matCorr) {
692         Double_t xrho,xx0;
693         ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);     
694         if (dir>0) xrho = -xrho; // outward should be negative
695         if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
696       }
697       else { // matCorr is not requested but time integral is
698         double d0 = xyz1[0]-xyz0[0];
699         double d1 = xyz1[1]-xyz0[1];
700         double d2 = xyz1[2]-xyz0[2];    
701         ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
702       }
703     }
704     if (updTime) seed->AddTimeStep(ds);
705     //
706     xpos = seed->GetX();
707   }
708   return kTRUE;
709 }
710
711 //_________________________________________________________________________
712 Bool_t AliITSUTrackerSA::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
713 {
714   // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
715   //  
716   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
717   //
718   int dir = lTo > lFrom ? 1:-1;
719   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
720   Bool_t checkFirst = kTRUE;
721   Bool_t limReached = kFALSE;
722   while(lFrom!=lTo) {
723     if (lrFr) {
724       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
725       checkFirst = kFALSE;
726     }
727     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
728     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
729     //
730     // go the entrance of the layer, assuming no materials in between
731     double xToGo = lrTo->GetR(-dir);
732     if (rLim>0) {
733       if (dir>0) {
734         if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
735       }
736       else {
737         if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
738       }
739     }
740     //    double xts = xToGo;
741     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
742       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
743       //      seed->Print("etp");
744       return kFALSE;
745     }
746     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
747       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
748       //seed->Print("etp");
749       return kFALSE;
750     }
751     lrFr = lrTo;
752     if (limReached) break;
753   }
754   return kTRUE;
755   //
756 }
757
758 //_________________________________________________________________________
759 Bool_t AliITSUTrackerSA::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
760 {
761   // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
762   //  
763   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
764   //
765   int dir = lTo > lFrom ? 1:-1;
766   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
767   Bool_t checkFirst = kTRUE;
768   while(lFrom!=lTo) {
769     if (lrFr) {
770       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
771       checkFirst = kFALSE;
772     }
773     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
774     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
775     //
776     // go the entrance of the layer, assuming no materials in between
777     double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
778     //
779     //    double xts = xToGo;
780     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
781       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
782       //      seed->Print("etp");
783       return kFALSE;
784     }
785     if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
786     //
787 #ifdef _ITSU_DEBUG_
788     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
789 #endif
790     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
791       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
792       //seed->Print("etp");
793       return kFALSE;
794     }
795     lrFr = lrTo;
796   }
797   return kTRUE;
798   //
799 }
800
801 //_________________________________________________________________________
802 Bool_t AliITSUTrackerSA::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
803 {
804   // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
805   // If check is requested, do this only provided the track has not exited the layer already
806   double xToGo = lr->GetR(dir);
807   if (check) { // do we need to track till the surface of the current layer ?
808     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
809     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
810     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
811   }
812   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
813   // go via layer to its boundary, applying material correction.
814   if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
815   //
816   return kTRUE;
817   //
818 }
819
820 //_________________________________________________________________________
821 Bool_t AliITSUTrackerSA::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
822 {
823   // go to the entrance of lr in direction dir, w/o applying material corrections.
824   // If check is requested, do this only provided the track did not reach the layer already
825   double xToGo = lr->GetR(-dir);
826   if (check) { // do we need to track till the surface of the current layer ?
827     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
828     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
829     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
830   }
831   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
832   // go via layer to its boundary, applying material correction.
833   if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
834   return kTRUE;
835   //
836 }
837
838 //____________________________________________________
839 Double_t AliITSUTrackerSA::GetMaterialBudget(const double* pnt0,const double* pnt1, double& x2x0, double& rhol) const
840 {
841   double par[7];
842   if (fUseMatLUT && fMatLUT) {
843     double d = fMatLUT->GetMatBudget(pnt0,pnt1,par);
844     x2x0 = par[AliITSUMatLUT::kParX2X0];
845     rhol = par[AliITSUMatLUT::kParRhoL];    
846     return d;
847   }
848   else {
849     MeanMaterialBudget(pnt0,pnt1,par);
850     x2x0 = par[1];
851     rhol = par[0]*par[4];    
852     return par[4];
853   }
854 }
855
856 //____________________________________________________________________
857 Double_t AliITSUTrackerSA::Curvature(Double_t x1,Double_t y1,Double_t 
858 x2,Double_t y2,Double_t x3,Double_t y3)
859 {
860
861   //calculates the curvature of track  
862   Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
863   if(den==0) return 0;
864   Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
865   Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
866   Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
867   Double_t xc=-a/2.;
868
869   if((a*a+b*b-4*c)<0) return 0;
870   Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
871   if(rad==0) return 0;
872   
873   if((x1>0 && y1>0 && x1<xc)) rad*=-1;
874   if((x1<0 && y1>0 && x1<xc)) rad*=-1;
875   //  if((x1<0 && y1<0 && x1<xc)) rad*=-1;
876   // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
877   
878   return 1/rad;
879  
880 }
881
882 #ifdef __DEBUG__ 
883 //____________________________________________________
884 void AliITSUTrackerSA::PrintInfo(TString what) {
885   //
886   if( what.Contains("clusters") ) {
887     cout << "Dumping clusters info" << endl;
888     for ( int i = 0; i < 7; ++i ) {
889       cout << "**** Layer " << i << " ****" << endl;
890       for ( int c = 0; c < fNClusters[i]; ++c ) {
891         cout << "*** Cluster " << c << " ***" <<endl;
892         cout << fClusters[i][c] << endl;
893       }
894     }
895   }
896   //
897   if( what.Contains("doublets") ) {
898     cout << "Dumping doublets info" << endl;
899     for ( int i = 0; i < 6; ++i ) {
900       cout << "**** Doublets array " << i << " ****" << endl;
901       for ( int c = 0; c < fNDoublets[i]; ++c ) {
902         cout << "*** Doublet " << c << " ***" <<endl;
903         cout << fDoublets[i][c] << endl;
904       }
905     }
906   }
907 }
908
909 #endif