]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUTrackerGlo.cxx
Stupid bug fix in new superlight mode (from Zurich airport)
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUTrackerGlo.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /*
17   TO CHECK
18   GetBz usage
19
20  */
21
22 //-------------------------------------------------------------------------
23 //  Implementation of the ITS Upgrade tracker mother class.              
24 //-------------------------------------------------------------------------
25 #include <TTree.h>
26 #include <Riostream.h> 
27 #include <TMath.h>
28 #include <TFile.h>
29
30 #include "AliITSUTrackerGlo.h"
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
33 #include "AliITSURecoDet.h"
34 #include "AliITSURecoSens.h"
35 #include "AliITSUReconstructor.h"
36 #include "AliITSReconstructor.h"
37 #include "AliITSUSeed.h"
38 #include "AliITSUClusterPix.h"
39 #include "AliITSUGeomTGeo.h"
40 #include "AliCodeTimer.h"
41 #include "AliRefArray.h"
42 using namespace AliITSUAux;
43 using namespace TMath;
44
45 //----------------- tmp stuff -----------------
46
47 ClassImp(AliITSUTrackerGlo)
48
49 const Double_t AliITSUTrackerGlo::fgkToler =  1e-6;// tolerance for layer on-surface check
50
51
52 //_________________________________________________________________________
53 AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
54 :  fReconstructor(rec)
55   ,fITS(0)
56   ,fMatLUT(0)
57   ,fCurrESDtrack(0)
58   ,fCurrESDtrMClb(kDummyLabel)
59   ,fCurrMass(kPionMass)
60   ,fCountProlongationTrials(0)
61   ,fCountITSin(0)
62   ,fCountITSout(0)
63   ,fCountITSrefit(0)
64   ,fNTracksESD(0)
65   ,fHypStore(100)
66   ,fLayerMaxCandidates(1000)
67   ,fLayerCandidates(0)
68   ,fNBranchesAdded(0)
69   ,fNCandidatesAdded(0)
70   ,fCurrHyp(0)
71   ,fWorkHyp(0)
72   ,fSeedsPool("AliITSUSeed",0)
73   ,fFreeSeedsID(0)
74   ,fESDIndex(0)
75   ,fNFreeSeeds(0)
76   ,fLastSeedID(0)
77   ,fNLrActive(0)
78   ,fDefTrackConds(0)
79   ,fCurrTrackCond(0)
80   ,fCurrActLrID(-1)
81   ,fCurrLayer(0)
82   ,fTrackPhaseID(-1)
83   ,fCurrPassID(-1)
84   ,fUseMatLUT(kFALSE)
85 #ifdef  _ITSU_TUNING_MODE_
86   ,fCHistoArrCorr(0)
87   ,fCHistoArrFake(0)
88 #endif
89 {
90   // Default constructor
91   if (rec) Init(rec);
92 }
93
94 //_________________________________________________________________________
95 AliITSUTrackerGlo::~AliITSUTrackerGlo()
96 {
97  // Default destructor
98  //  
99   delete[] fLayerCandidates;
100   if (fWorkHyp) fWorkHyp->SetTPCSeed(0); // this hypothesis does not own the seed
101   delete fWorkHyp;
102   delete fMatLUT;
103   //
104 #ifdef  _ITSU_TUNING_MODE_
105   if (fCHistoArrCorr || fCHistoArrFake) {
106     TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
107     ctrOut->cd();
108     AliInfo("Storing control histos");
109     //    ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
110     if (fCHistoArrCorr) {
111       fCHistoArrCorr->Write();
112       delete fCHistoArrCorr;
113     }
114     if (fCHistoArrFake) {
115       fCHistoArrFake->Write();
116       delete fCHistoArrFake;
117     }
118     ctrOut->Close();
119     delete ctrOut;
120     fCHistoArrCorr = 0;
121     fCHistoArrFake = 0;    
122   }
123 #endif 
124   //
125 }
126
127 //_________________________________________________________________________
128 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
129 {
130   // init with external reconstructor
131   //
132   fITS = rec->GetITSInterface();
133   fNLrActive = fITS->GetNLayersActive();
134   fWorkHyp = new AliITSUTrackHyp(fNLrActive);
135   // create material lookup table
136   const int kNTest = 1000;
137   const double kStepsPerCM=5;
138   fMatLUT  = new AliITSUMatLUT(fITS->GetRMin(),fITS->GetRMax(),Nint(kStepsPerCM*(fITS->GetRMax()-fITS->GetRMin())));
139   double zmn = 1e6;
140   for (int ilr=fITS->GetNLayers();ilr--;) {
141     AliITSURecoLayer* lr = fITS->GetLayer(ilr);
142     if (zmn>Abs(lr->GetZMin())) zmn = Abs(lr->GetZMin());
143     if (zmn>Abs(lr->GetZMax())) zmn = Abs(lr->GetZMax());
144   }
145   fMatLUT->FillData(kNTest,-zmn,zmn);
146   //
147   if (fLayerMaxCandidates<1) fLayerMaxCandidates = 1000;
148   fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
149   fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
150   fFreeSeedsID.Set(1000);
151   fESDIndex.Set(1000);
152   //
153 }
154
155 //_________________________________________________________________________
156 void AliITSUTrackerGlo::CreateDefaultTrackCond()
157 {
158   // creates default tracking conditions to be used when recoparam does not provide them
159
160   AliITSUTrackCond* cond = new AliITSUTrackCond();
161   //
162   cond->SetNLayers(fNLrActive);
163   cond->AddNewCondition(fNLrActive);
164   cond->AddGroupPattern( 0xffff, fNLrActive ); // require all layers hit
165   cond->Init();
166   //
167   fDefTrackConds.AddLast(cond);
168   //
169   AliInfo("Created conditions: ");
170   for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
171   //
172 }
173
174
175 //_________________________________________________________________________
176 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
177 {
178   //
179   AliCodeTimerAuto("",0);
180   SetTrackingPhase(kClus2Tracks);
181   //
182 #ifdef  _ITSU_TUNING_MODE_
183   if (!fCHistoArrCorr) fCHistoArrCorr = BookControlHistos("Corr");
184   if (!fCHistoArrFake) fCHistoArrFake = BookControlHistos("Fake");
185 #endif
186   static int evID = 0;
187   static TArrayF esdTrPt(fESDIndex.GetSize()); 
188   //
189   TObjArray *trackConds = 0;
190   //
191   fCountProlongationTrials = 0;
192   fCountITSin = 0;
193   fCountITSout = 0;
194   fCountITSrefit = 0;
195   //
196   ResetSeedsPool();
197   fNTracksESD = esdEv->GetNumberOfTracks();
198   AliInfo(Form("Will try to find prolongations for %d tracks",fNTracksESD));
199   int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
200   fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
201   if (nTrackCond<1) {
202     if (!fDefTrackConds.GetEntriesFast()) {
203       AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
204       CreateDefaultTrackCond();
205     }
206     trackConds = &fDefTrackConds;
207     nTrackCond = trackConds->GetEntriesFast();
208   } 
209   else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
210   //
211   static Bool_t first = kTRUE;
212   if (first) {
213     AliITSUReconstructor::GetRecoParam()->Print();
214     first = kFALSE;
215   }
216   fHypStore.Delete();
217   if (fHypStore.GetSize()<fNTracksESD) fHypStore.Expand(fNTracksESD+100);
218   //
219   fITS->ProcessClusters();
220   //
221 #ifdef  _ITSU_TUNING_MODE_
222   FlagSplitClusters(); // tmp RS
223 #endif
224   //
225   // the tracks will be reconstructed in decreasing pt order, sort them
226   if (fESDIndex.GetSize()<fNTracksESD) {
227     fESDIndex.Set(fNTracksESD+200);
228     esdTrPt.Set(fNTracksESD+200);
229   }
230   int* esdTrackIndex = fESDIndex.GetArray();
231   float* trPt = esdTrPt.GetArray();
232   for (int itr=fNTracksESD;itr--;) {
233     AliESDtrack* esdTr = esdEv->GetTrack(itr);
234     esdTr->SetStatus(AliESDtrack::kITSupg);   // Flag all tracks as participating in the ITSup reco
235     trPt[itr] = esdTr->Pt();
236   }
237   Sort(fNTracksESD,trPt,esdTrackIndex,kTRUE);    
238   //
239   for (int icnd=0;icnd<nTrackCond;icnd++) {
240     fCurrPassID = icnd;
241     fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
242     if (!fCurrTrackCond->IsInitDone()) fCurrTrackCond->Init();
243     // select ESD tracks to propagate
244     for (int itr=0;itr<fNTracksESD;itr++) {
245       int trID = esdTrackIndex[itr];
246       fCurrESDtrack = esdEv->GetTrack(trID);
247       fCurrESDtrMClb = fCurrESDtrack->GetLabel();
248       //
249       
250       if (!NeedToProlong(fCurrESDtrack,trID)) continue;  // are we interested in this track?
251       /*
252         // if specific tracks should be checked, create a global array int wtc[] = {evITS*10000+trID, ...};
253       Bool_t dbg = kFALSE;
254       int nwtc = sizeof(wtc)/sizeof(int);
255
256       for (int k=0;k<nwtc;k++) {
257         if (wtc[k]==evID*10000+trID) {
258           dbg = kTRUE;
259           printf("At esdTr: %d MC: %d\n",wtc[k],fCurrESDtrMClb);
260           break;
261         }
262       }
263       AliLog::SetClassDebugLevel("AliITSUTrackerGlo",dbg ? 10:0);
264       */
265 #ifdef _ITSU_DEBUG_
266       AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMassForTracking(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
267 #endif
268       FindTrack(fCurrESDtrack, trID);
269     }   
270     //
271 #ifdef _ITSU_DEBUG_
272     if (AliDebugLevelClass()>+2) {
273       AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
274       fHypStore.Print();
275     }
276 #endif
277     FinalizeHypotheses();
278 #ifdef  _ITSU_TUNING_MODE_
279     CheckClusterUsage(); //!!RS
280 #endif
281     //
282     //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);     // in case wtc array was defined, uncomment this
283   }
284   //
285   AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,fNTracksESD));
286   evID++;
287   return 0;
288 }
289
290 //_________________________________________________________________________
291 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
292 {
293   //
294   // Do outward fits in ITS
295   AliCodeTimerAuto("",0);
296   //
297   SetTrackingPhase(kPropBack);
298   fNTracksESD = esdEv->GetNumberOfTracks();
299   fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
300 #ifdef _ITSU_DEBUG_
301   AliDebug(1,Form("Will propagate back %d tracks",fNTracksESD));
302 #endif
303   //
304   double bz0 = GetBz();
305   Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
306   AliITSUTrackHyp dummyTr;
307   const double kWatchStep=10.; // for larger steps watch arc vs segment difference
308   Double_t times[AliPID::kSPECIESC];
309   //
310   for (int itr=0;itr<fNTracksESD;itr++) {
311     fCurrESDtrack = esdEv->GetTrack(itr);
312     fCurrESDtrMClb = fCurrESDtrack->GetLabel();
313     // Start time integral and add distance from current position to vertex 
314     if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
315     if (!fCurrESDtrack->IsOn(AliESDtrack::kTPCin)) continue; // skip ITS s.a.
316     //
317     fCurrESDtrack->GetXYZ(xyzTrk); 
318     Double_t dst = 0.;     // set initial track lenght, tof
319     {
320       double dxs = xyzTrk[0] - xyzVtx[0];
321       double dys = xyzTrk[1] - xyzVtx[1];
322       //      double dzs = xyzTrk[2] - xyzVtx[2];
323       // RS: for large segment steps d use approximation of cicrular arc s by
324       // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
325       // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
326       // Hence s^2/d^2 = (1+1/6 p^2)^2
327       dst = dxs*dxs + dys*dys;
328       if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
329         double crv = Abs(fCurrESDtrack->GetC(bz0));
330         double fcarc = 1.+crv*crv*dst/6.;
331         dst *= fcarc*fcarc;
332       }
333       // RS: temporary hack since we don't have SPD vertex:      
334       //      dst += dzs*dzs;
335       dst *= 1+fCurrESDtrack->GetTgl()*fCurrESDtrack->GetTgl();
336       dst = Sqrt(dst); 
337     }
338     //    
339     fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
340     //
341     if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
342       dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
343       dummyTr.StartTimeIntegral();
344       dummyTr.AddTimeStep(dst);
345       dummyTr.GetIntegratedTimes(times,AliPID::kSPECIESC); 
346       fCurrESDtrack->SetIntegratedTimes(times);
347       fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
348       continue;
349     }
350     //
351     fCurrHyp = GetTrackHyp(itr);
352     fCurrMass = fCurrHyp->GetMass();
353     fCurrHyp->StartTimeIntegral();
354     fCurrHyp->AddTimeStep(dst);
355     fCurrHyp->ResetCovariance(10000);
356     int nclFit = 0;
357     double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax(),nclFit);
358     if (chi2>0) { // propagate to exit from the ITS/TPC screen
359       int ndf = nclFit*2-5;
360       if (ndf>0) chi2 /= ndf;
361       fCurrHyp->SetChi2(chi2);
362       UpdateESDTrack(fCurrHyp,AliESDtrack::kITSout);
363     }
364     else {
365 #ifdef _ITSU_DEBUG_
366       AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
367 #endif
368       //fCurrHyp->AliExternalTrackParam::Print();
369       //fCurrHyp->GetWinner()->Print();
370     }
371     //
372   }
373   //
374   AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
375                fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
376   //
377   return 0;
378 }
379
380 //_________________________________________________________________________
381 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
382 {
383   //
384   // refit the tracks inward, using current cov.matrix
385   AliCodeTimerAuto("",0);
386   //
387   SetTrackingPhase(kRefitInw);
388   fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
389   fNTracksESD = esdEv->GetNumberOfTracks();
390   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
391
392   AliDebug(1,Form("Will refit inward %d tracks",fNTracksESD));
393   //  Bool_t uselogMS = AliExternalTrackParam::GetUseLogTermMS();
394   //  AliExternalTrackParam::SetUseLogTermMS(kTRUE);
395   //
396   for (int itr=0;itr<fNTracksESD;itr++) {
397     fCurrESDtrack = esdEv->GetTrack(itr);
398     fCurrESDtrMClb = fCurrESDtrack->GetLabel();
399     // Start time integral and add distance from current position to vertex 
400     UInt_t trStat = fCurrESDtrack->GetStatus();
401     if ( !(trStat & AliESDtrack::kITSout) ) continue;
402     if (   trStat & AliESDtrack::kITSrefit ) continue; // already done
403     if ( !(trStat & AliESDtrack::kTPCin)   ) continue; // skip ITS s.a.
404     if (  (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
405     //
406     fCurrHyp  = GetTrackHyp(itr);
407     fCurrMass = fCurrHyp->GetMass();
408     *fCurrHyp = *fCurrESDtrack;  // fetch current ESDtrack kinematics
409     //
410     int nclFit = 0;
411     //    fCurrHyp->ResetCovariance(1000); // in case we want to do SA refit
412     double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin(),nclFit);
413     if (chi2>0) { // propagate up to inside radius of the beam pipe      
414       int ndf = nclFit*2-5;
415       if (ndf>0) chi2 /= ndf;
416       fCurrHyp->SetChi2(chi2);
417       UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
418     }
419     else {
420 #ifdef _ITSU_DEBUG_
421       AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
422 #endif
423     }
424   }    
425   //  AliExternalTrackParam::SetUseLogTermMS(uselogMS);
426   //
427   AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
428                fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
429   //
430   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
431   return 0;
432 }
433
434 //_________________________________________________________________________
435 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
436 {
437   // read from tree (if pointer provided) or directly from the ITS reco interface
438   //
439   return fReconstructor->LoadClusters(treeRP);
440
441
442 //_________________________________________________________________________
443 void AliITSUTrackerGlo::UnloadClusters()
444 {
445   //
446   // Remove clusters from the memory 
447   //
448   /*
449   AliITSURecoDet *det=fReconstructor->GetITSInterface();
450   Int_t nlayers=det->GetNLayersActive();
451   for (Int_t i=0; i<nlayers; i++) {
452       TClonesArray *clusters=*(det->GetLayerActive(i)->GetClustersAddress());
453       clusters->Delete();
454   }
455   */
456
457
458 //_________________________________________________________________________
459 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
460 {
461   //
462   // To be implemented 
463   //
464   
465   Info("GetCluster","To be implemented");
466   return 0x0;
467
468
469 //_________________________________________________________________________
470 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr, Int_t esdID)
471 {
472   // do we need to match this track to ITS?
473   //
474   static double bz = GetBz();
475   //
476   // is it already prolonged
477   if (esdTr->IsOn(AliESDtrack::kITSin)) return kFALSE;
478   AliITSUTrackHyp* hyp = GetTrackHyp(esdID); // in case there is unfinished hypothesis from prev.pass, clean it
479   if (hyp) {
480     CleanHypothesis(hyp);
481     if (hyp->GetSkip()) return kFALSE; // need to skip this
482     return kTRUE;
483   }
484   //
485   if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
486       esdTr->IsOn(AliESDtrack::kTPCout) ||
487       esdTr->IsOn(AliESDtrack::kITSin)  ||
488       esdTr->GetKinkIndex(0)>0) return kFALSE;
489   //
490   if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
491   //
492   float dtz[2];
493   esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz); 
494   // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
495   if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation()) 
496        && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
497            Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
498   //
499   //  if (esdTr->Pt()<3) return kFALSE;//RS
500   return kTRUE;
501 }
502
503 //_________________________________________________________________________
504 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
505 {
506   // find prolongaion candidates finding for single seed
507   //
508   AliITSUSeed seedUC;  // copy of the seed from the upper layer
509   AliITSUSeed seedT;   // transient seed between the seedUC and new prolongation hypothesis
510   //
511   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
512   AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID);  // initialize prolongations hypotheses tree
513   if (!hypTr || hypTr->GetSkip()) return;
514   //
515   double bz = GetBz();
516   fCurrHyp = fWorkHyp;
517   fCurrHyp->InitFrom(hypTr);
518   //
519   AliITSURecoSens *hitSens[AliITSURecoLayer::kMaxSensMatching];
520   //
521   int ilaUp = fNLrActive;     // previous active layer really checked (some may be excluded!)
522   //
523   for (int ila=fNLrActive;ila--;ilaUp=fCurrActLrID) {
524     if (fCurrTrackCond->IsLayerExcluded(ila)) continue;
525     fCurrActLrID = ila;
526     fCurrLayer = fITS->GetLayerActive(ila);
527     Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
528     //
529     // for the outermost layer the seed is created from the ESD track
530     int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
531     int maxNBranches   = fCurrTrackCond->GetMaxBranches(ila);
532     int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
533     //
534     for (int isd=0;isd<nSeedsUp;isd++) {
535       AliITSUSeed* seedU;
536       if (ilaUp==fNLrActive) {
537         seedU = 0;
538         seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
539       }
540       else {
541         seedU = fCurrHyp->GetSeed(ilaUp,isd);  // seed on prev.active layer to prolong  
542         if (seedU->IsKilled()) continue;
543         seedUC = *seedU;                       // its copy will be prolonged
544         seedUC.SetParent(seedU);        
545       }
546       seedUC.ResetFMatrix();                    // reset the matrix for propagation to next layer
547       // go till next active layer
548 #ifdef _ITSU_DEBUG_
549       AliDebug(2,Form("working on Lr:%d Seed:%d of %d for esdID=%d (MClb:%d) | pT=%.3f",ila,isd,nSeedsUp,esdID,fCurrESDtrMClb,seedUC.Pt()));
550 #endif
551       if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
552         //
553 #ifdef _ITSU_DEBUG_
554         AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
555 #endif
556         // Check if the seed satisfies to track definition
557         if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
558         continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
559       }
560       if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
561         if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
562         continue;
563       }
564       /*
565       //RS toremove
566       int mcquest = -1;
567       if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
568         mcquest = fCurrESDtrMClb;
569         seedUC.Print("etp");
570         printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
571       }
572       //
573       int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest);  // find detectors which may be hit by the track
574       */
575       int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens);  // find detectors which may be hit by the track
576 #ifdef _ITSU_DEBUG_
577       AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
578 #endif
579       //
580       seedUC.SetLr(ila);
581       //
582       for (int isn=nsens;isn--;) {
583         AliITSURecoSens* sens = hitSens[isn];
584         int ncl = sens->GetNClusters();
585         if (!ncl) continue;
586         seedT = seedUC;
587         //
588         // We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
589         // since the transport matrix should be defined in this frame.
590         double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
591         if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
592 #ifdef _ITSU_DEBUG_
593           if (AliDebugLevelClass()>2) {
594             printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
595             seedT.Print("etp");
596           }
597 #endif
598           continue;
599         }
600         if (xs<seedT.GetX()) {
601           if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
602         }
603         else { // some low precision tracks may hit the sensor plane outside of the layer radius
604 #ifdef _ITSU_DEBUG_
605           if (AliDebugLevelClass()>2) {
606             if (!seedT.ContainsFake()) {
607               printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
608               sens->Print();
609               seedT.Print("etp");           
610             }
611           }
612 #endif
613           if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
614         }
615         //      if (!seedT.PropagateToX(xs,bz)) continue;
616         //      if (!seedT.Rotate(sens->GetPhiTF())) continue;
617         if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
618         //
619         int clID0 = sens->GetFirstClusterId();
620         for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
621           //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
622           int clID = clID0 + icl;
623           if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
624           int res = CheckCluster(&seedT,ila,clID0+icl);
625           //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
626           //
627           if (res==kStopSearchOnSensor) break;     // stop looking on this sensor
628           if (res==kClusterNotMatching) continue;  // cluster does not match
629           // cluster is matching and it was added to the hypotheses tree
630         }
631       }
632       // cluster search is done. Do we need to have a version of this seed skipping current layer
633       if (!NeedToKill(&seedUC,kMissingCluster)) {
634         AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
635         if (nsens) { // to do: make penalty to account for probability to miss the cluster for good reason
636           double penalty = -fCurrTrackCond->GetMissPenalty(ila);
637           seedSkp->SetChi2Cl(penalty);
638         }
639         if (seedSkp->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(ila)) {
640           MarkSeedFree(seedSkp);
641         }
642         else AddSeedBranch(seedSkp);
643       }
644       // transfer the new branches of the seed to the hypothesis container
645       if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
646       //
647     }
648     if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
649     //    ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
650 #ifdef _ITSU_DEBUG_
651     if (AliDebugLevelClass()>2) { //RS
652       printf(">>> All hypotheses on lr %d: \n",ila);
653       for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
654     }
655 #endif
656     //
657     /*
658     for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
659       if (ila!=0) continue;
660       double vecL[5] = {0};
661       double matL[15] = {0};
662       AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
663       while(sp->GetParent()) {
664         sp->Smooth(vecL,matL);
665         if (sp->GetLayerID()>=fNLrActive-1) break;
666         sp = (AliITSUSeed*)sp->GetParent();
667       }
668       */
669   }
670   //
671   SaveReducedHypothesesTree(hypTr);
672 #ifdef _ITSU_DEBUG_
673   if (AliDebugLevelClass()>1) {
674     printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
675     hypTr->Print("l");
676   }
677 #endif
678   fCurrHyp = 0;
679   //  AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
680   //
681 }
682
683 //_________________________________________________________________________
684 AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
685 {
686   // init prolongaion candidates finding for single seed
687   AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
688   if (hyp) return hyp;
689   //
690   fCountProlongationTrials++;
691   //
692   fCurrMass = esdTr->GetMassForTracking();
693   //
694   hyp = new AliITSUTrackHyp(fNLrActive);
695   hyp->SetESDTrack(esdTr);
696   hyp->SetUniqueID(esdID);
697   hyp->SetMass(fCurrMass);
698   hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
699   SetTrackHyp(hyp,esdID);
700   Bool_t lut = fUseMatLUT; // for propagation from TPC use TGeo material budget always
701   fUseMatLUT = kFALSE;
702   if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
703   fUseMatLUT = lut;
704   //
705   return hyp;
706   // TO DO
707 }
708
709 //_________________________________________________________________________
710 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
711 {
712   // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
713   //  
714   //
715   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
716   //
717   int dir = lTo > lFrom ? 1:-1;
718   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
719   Bool_t checkFirst = kTRUE;
720   Bool_t limReached = kFALSE;
721   while(lFrom!=lTo) {
722     if (lrFr) {
723       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
724       checkFirst = kFALSE;
725     }
726     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
727     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
728     //
729     // go the entrance of the layer, assuming no materials in between
730     double xToGo = lrTo->GetR(-dir);
731     if (rLim>0) {
732       if (dir>0) {
733         if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
734       }
735       else {
736         if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
737       }
738     }
739     //    double xts = xToGo;
740     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
741       //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
742       //seed->Print("etp");
743       //      
744       return kFALSE;
745     }
746 #ifdef _ITSU_DEBUG_
747     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
748 #endif
749     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
750       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
751       //seed->Print("etp");
752       return kFALSE;
753     }
754     lrFr = lrTo;
755     if (limReached) break;
756   }
757   return kTRUE;
758   //
759 }
760
761 //_________________________________________________________________________
762 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
763 {
764   // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
765   //  
766   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
767   //
768   int dir = lTo > lFrom ? 1:-1;
769   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
770   Bool_t checkFirst = kTRUE;
771   Bool_t limReached = kFALSE;
772   while(lFrom!=lTo) {
773     if (lrFr) {
774       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
775       checkFirst = kFALSE;
776     }
777     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
778     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
779     //
780     // go the entrance of the layer, assuming no materials in between
781     double xToGo = lrTo->GetR(-dir);
782     if (rLim>0) {
783       if (dir>0) {
784         if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
785       }
786       else {
787         if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
788       }
789     }
790     //    double xts = xToGo;
791     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
792       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
793       //      seed->Print("etp");
794       return kFALSE;
795     }
796 #ifdef _ITSU_DEBUG_
797     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
798 #endif
799     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
800       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
801       //seed->Print("etp");
802       return kFALSE;
803     }
804     lrFr = lrTo;
805     if (limReached) break;
806   }
807   return kTRUE;
808   //
809 }
810
811 //_________________________________________________________________________
812 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
813 {
814   // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
815   //  
816   if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
817   //
818   int dir = lTo > lFrom ? 1:-1;
819   AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
820   Bool_t checkFirst = kTRUE;
821   while(lFrom!=lTo) {
822     if (lrFr) {
823       if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
824       checkFirst = kFALSE;
825     }
826     AliITSURecoLayer* lrTo =  fITS->GetLayer( (lFrom+=dir) );
827     if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
828     //
829     // go the entrance of the layer, assuming no materials in between
830     double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
831     //
832     //    double xts = xToGo;
833     if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
834       //      printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
835       //      seed->Print("etp");
836       return kFALSE;
837     }
838     if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
839     //
840 #ifdef _ITSU_DEBUG_
841     AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
842 #endif
843     if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
844       //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
845       //seed->Print("etp");
846       return kFALSE;
847     }
848     lrFr = lrTo;
849   }
850   return kTRUE;
851   //
852 }
853
854 //_________________________________________________________________________
855 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
856 {
857   // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
858   // If check is requested, do this only provided the track has not exited the layer already
859   double xToGo = lr->GetR(dir);
860   if (check) { // do we need to track till the surface of the current layer ?
861     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
862 #ifdef _ITSU_DEBUG_
863     AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
864 #endif
865     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
866     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
867   }
868   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
869   // go via layer to its boundary, applying material correction.
870   if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
871   //
872   return kTRUE;
873   //
874 }
875
876 //_________________________________________________________________________
877 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
878 {
879   // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
880   // If check is requested, do this only provided the track has not exited the layer already
881   double xToGo = lr->GetR(dir);
882   if (check) { // do we need to track till the surface of the current layer ?
883     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
884     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
885     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
886   }
887 #ifdef _ITSU_DEBUG_
888   AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
889 #endif
890   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
891   // go via layer to its boundary, applying material correction.
892   if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
893   //
894   return kTRUE;
895   //
896 }
897
898
899 //_________________________________________________________________________
900 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
901 {
902   // go to the entrance of lr in direction dir, w/o applying material corrections.
903   // If check is requested, do this only provided the track did not reach the layer already
904   double xToGo = lr->GetR(-dir);
905   if (check) { // do we need to track till the surface of the current layer ?
906     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
907     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
908     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
909   }
910   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
911   // go via layer to its boundary, applying material correction.
912   if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
913   return kTRUE;
914   //
915 }
916
917 //_________________________________________________________________________
918 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
919 {
920   // go to the entrance of lr in direction dir, w/o applying material corrections.
921   // If check is requested, do this only provided the track did not reach the layer already
922   double xToGo = lr->GetR(-dir);
923   if (check) { // do we need to track till the surface of the current layer ?
924     double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
925     if      (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
926     else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
927   }
928   if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
929   // go via layer to its boundary, applying material correction.
930   if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
931   return kTRUE;
932   //
933 }
934
935 //_________________________________________________________________________
936 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
937 {
938   // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
939   // as well as some aux info
940   double bz = GetBz();
941   AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
942   seed->GetXYZ(&fTrImpData[kTrXIn]);    // lab position at the entrance from above
943   static AliExternalTrackParam sc;   // seed copy for manipulations
944   sc = *seed;
945   //
946   fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
947   if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
948   BringTo02Pi(fTrImpData[kTrPhiIn]);
949   double dr  = lrA->GetDR();                              // approximate X dist at the inner radius
950   if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
951     // special case: track does not reach inner radius, might be tangential
952     double r = sc.GetD(0,0,bz);
953     double x;
954     if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
955     dr = Abs(sc.GetX() - x);
956     if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
957   }
958   //
959   fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
960   BringTo02Pi(fTrImpData[kTrPhiOut]);
961   double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
962   double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
963   sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
964   sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
965   double phi0  = MeanPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
966   double dphi0 = DeltaPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
967   //
968   fTrImpData[kTrPhi0] = phi0;
969   fTrImpData[kTrZ0]   = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
970   dphi0 += sgy/lrA->GetR();
971   fTrImpData[kTrDPhi] =  dphi0<PiOver2() ? dphi0 : PiOver2();
972   fTrImpData[kTrDZ]   = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn])   + sgz;
973   //  
974   return kTRUE;
975 }
976
977 //________________________________________
978 void AliITSUTrackerGlo::ResetSeedsPool()
979 {
980   // mark all seeds in the pool as unused
981   AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
982   fNFreeSeeds = 0;
983   fSeedsPool.Clear(); // seeds don't allocate memory
984 }
985
986
987 //________________________________________
988 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd) 
989 {
990   // account that this seed is "deleted" 
991   int id = sd->GetPoolID();
992   if (id<0) {
993     AliError(Form("Freeing of seed %p NOT from the pool is requested",sd)); 
994     return;
995   }
996   //  AliInfo(Form("%d %p",id, seed));
997   fSeedsPool.RemoveAt(id);
998   if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
999   fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
1000 }
1001
1002 //_________________________________________________________________________
1003 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID) 
1004 {
1005   // Check if the cluster (in tracking frame!) is matching to track. 
1006   // The track must be already propagated to sensor tracking frame.
1007   // Returns:  kStopSearchOnSensor if the search on given sensor should be stopped, 
1008   //           kClusterMatching    if the cluster is matching
1009   //           kClusterMatching    otherwise
1010   //
1011   const double kTolerX = 5e-4;
1012   //  AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1013   AliCluster *cl = fCurrLayer->GetCluster(clID);
1014   //
1015   Bool_t goodCl = kFALSE;
1016   int currLabel = Abs(fCurrESDtrMClb);
1017   //
1018   if (cl->GetLabel(0)>=0) {
1019     for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
1020   }
1021   //
1022   if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
1023     if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
1024       //
1025 #ifdef _ITSU_DEBUG_
1026       if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1027         AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
1028         track->Print("etp");
1029         cl->Print();
1030       }
1031 #endif
1032       return kStopSearchOnSensor; // propagation failed, seedT is intact
1033     }
1034   }
1035   double dy = cl->GetY()-track->GetY();
1036   double dz = cl->GetZ()-track->GetZ();
1037   //
1038   double dy2 = dy*dy;
1039   double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
1040     fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
1041   if (dy2>tol2) {                          // the clusters are sorted in Z(col) then in Y(row). 
1042     //
1043 #ifdef _ITSU_DEBUG_
1044     if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {    
1045       AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
1046       track->Print("etp");
1047       cl->Print();
1048       PrintSeedClusters(track);
1049     }
1050 #endif
1051     // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
1052     // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
1053     if (dy<0) return kStopSearchOnSensor;  // No chance that other cluster of this sensor will match (all Y's will be even larger)
1054     else      return kClusterNotMatching;   // Other clusters may match
1055   }
1056   double dz2 = dz*dz;
1057   tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
1058     fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
1059   if (dz2>tol2) {
1060 #ifdef _ITSU_DEBUG_
1061     if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
1062       AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
1063       track->Print("etp");
1064       cl->Print();
1065       PrintSeedClusters(track);
1066     }
1067 #endif
1068     return kClusterNotMatching; // Other clusters may match
1069   }
1070   //
1071   // check chi2
1072   Double_t p[2]={cl->GetY(), cl->GetZ()};
1073   Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
1074   double chi2 = track->GetPredictedChi2(p,cov);
1075   //
1076   if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
1077 #ifdef _ITSU_DEBUG_
1078     if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
1079       AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1080                       lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
1081       track->Print("etp");
1082       cl->Print("");
1083       PrintSeedClusters(track);
1084     }
1085 #endif
1086     return kClusterNotMatching;
1087   }
1088   //
1089   track = NewSeedFromPool(track);  // input track will be reused, use its clone for updates
1090   if (!track->Update()) {
1091 #ifdef _ITSU_DEBUG_
1092     if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
1093       AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
1094       track->Print("etp");
1095       cl->Print("");
1096       PrintSeedClusters(track);
1097     }
1098 #endif
1099     MarkSeedFree(track);
1100     return kClusterNotMatching;
1101   }
1102   track->SetChi2Cl(chi2);
1103   track->SetLrClusterID(lr,clID);
1104   //
1105   if (track->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(lr)) {
1106 #ifdef _ITSU_DEBUG_
1107     if (AliDebugLevelClass()>2 && goodCl &&  !track->ContainsFake()) {
1108       AliDebug(2,Form("Lost good cl on L:%d : Chi2Glo=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1109                       lr,track->GetChi2GloNrm(),fCurrTrackCond->GetMaxChi2GloNrm(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb)); 
1110       track->Print("etp");
1111       cl->Print("");
1112       PrintSeedClusters(track);
1113     }
1114 #endif
1115     MarkSeedFree(track);
1116     return kClusterNotMatching;
1117   }
1118   //  cl->IncreaseClusterUsage(); // do this only for winners
1119   //
1120   track->SetFake(!goodCl);
1121   //
1122 #ifdef _ITSU_DEBUG_
1123   AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
1124                goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
1125 #endif
1126   //
1127   AddSeedBranch(track);
1128   //
1129   return kClusterMatching;
1130 }
1131
1132 //_________________________________________________________________________
1133 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
1134 {
1135   // check if the seed should not be discarded
1136   const UShort_t kMask = 0xffff;
1137   if (flag==kMissingCluster) {
1138     int lastChecked = seed->GetLayerID();
1139     UShort_t patt = seed->GetHitsPattern();
1140     if (lastChecked) patt |= (~(kMask<<lastChecked))&fCurrTrackCond->GetAllowedLayers(); // not all layers were checked, complete unchecked ones by potential hits
1141     Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
1142     return !seedOK;
1143   }
1144   return kTRUE;
1145 }
1146
1147 //______________________________________________________________________________
1148 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
1149 {
1150   // propagate seed to given x applying material correction if requested
1151   const Double_t kEpsilon = 1e-5;
1152   Double_t xpos     = seed->GetX();
1153   Int_t dir         = (xpos<xToGo) ? 1:-1;
1154   Double_t xyz0[3],xyz1[3];
1155   //
1156   Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1157   if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
1158   while ( (xToGo-xpos)*dir > kEpsilon){
1159     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1160     Double_t x    = xpos+step;
1161     Double_t bz=GetBz();   // getting the local Bz
1162     if (!seed->PropagateToX(x,bz)) return kFALSE;
1163     double ds = 0;
1164     if (matCorr || updTime) {
1165       xyz0[0]=xyz1[0]; // global pos at the beginning of step
1166       xyz0[1]=xyz1[1];
1167       xyz0[2]=xyz1[2];
1168       seed->GetXYZ(xyz1);    //  // global pos at the end of step
1169       if (matCorr) {
1170         Double_t xrho,xx0;
1171         ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);     
1172         if (dir>0) xrho = -xrho; // outward should be negative
1173         if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
1174       }
1175        else { // matCorr is not requested but time integral is
1176         double d0 = xyz1[0]-xyz0[0];
1177         double d1 = xyz1[1]-xyz0[1];
1178         double d2 = xyz1[2]-xyz0[2];    
1179         ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1180       }     
1181     }
1182     if (updTime) seed->AddTimeStep(ds);
1183     xpos = seed->GetX();
1184   }
1185   return kTRUE;
1186 }
1187
1188 //______________________________________________________________________________
1189 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr) 
1190 {
1191   // propagate seed to given x applying material correction if requested
1192   const Double_t kEpsilon = 1e-5;
1193   Double_t xpos     = seed->GetX();
1194   Int_t dir         = (xpos<xToGo) ? 1:-1;
1195   Double_t xyz0[3],xyz1[3];
1196   //
1197 #ifdef _ITSU_DEBUG_
1198   if (AliDebugLevelClass()>1) {
1199     AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
1200     seed->AliExternalTrackParam::Print();
1201   }
1202 #endif
1203   Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1204   if (matCorr || updTime) seed->GetXYZ(xyz1);   //starting global position
1205   while ( (xToGo-xpos)*dir > kEpsilon){
1206     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1207     Double_t x    = xpos+step;
1208     Double_t bz=GetBz();   // getting the local Bz
1209     if (!seed->PropagateTo(x,bz))  return kFALSE;
1210     double ds = 0;
1211     if (matCorr || updTime) {
1212       xyz0[0]=xyz1[0]; // global pos at the beginning of step
1213       xyz0[1]=xyz1[1];
1214       xyz0[2]=xyz1[2];
1215       seed->GetXYZ(xyz1);    //  // global pos at the end of step
1216       //
1217       if (matCorr) {
1218         Double_t xrho,xx0;
1219         ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);     
1220         if (dir>0) xrho = -xrho; // outward should be negative
1221         if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1222       }
1223       else { // matCorr is not requested but time integral is
1224         double d0 = xyz1[0]-xyz0[0];
1225         double d1 = xyz1[1]-xyz0[1];
1226         double d2 = xyz1[2]-xyz0[2];    
1227         ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1228       }
1229     }
1230     if (updTime) seed->AddTimeStep(ds);
1231     //
1232     xpos = seed->GetX();
1233   }
1234 #ifdef _ITSU_DEBUG_
1235   if (AliDebugLevelClass()>1) {
1236     AliDebug(2,Form("After Propagate to X=%f",xToGo));
1237     seed->AliExternalTrackParam::Print();
1238   }
1239 #endif
1240   return kTRUE;
1241 }
1242
1243 //______________________________________________________________________________
1244 Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
1245 {
1246   // finalize single TPC track hypothesis
1247   if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
1248   AliITSUSeed* winner = 0;
1249   fCurrHyp = hyp;
1250   fCurrMass = hyp->GetMass();
1251   if (!(winner=hyp->DefineWinner(fCurrTrackCond->GetActiveLrInner()))) return kFALSE;
1252   CookMCLabel(hyp);
1253   //
1254 #ifdef  _ITSU_TUNING_MODE_  // fill tuning histos
1255   TObjArray* dest = hyp->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1256   if (dest) {
1257     AliITSUSeed* sd = winner;
1258     double htrPt = hyp->Pt();
1259     do {
1260       int clID,lrID;
1261       if ( (clID=sd->GetLrCluster(lrID))<0 ) continue;
1262       ((TH2*)dest->At( GetHistoID(lrID,kHChi2Nrm     ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2GloNrm());
1263       ((TH2*)dest->At( GetHistoID(lrID,kHBestInBranch,fCurrPassID) ))->Fill(htrPt,sd->GetOrdBranch());
1264       ((TH2*)dest->At( GetHistoID(lrID,kHBestInCand  ,fCurrPassID) ))->Fill(htrPt,sd->GetOrdCand());
1265       //
1266       if (dest==fCHistoArrFake && !sd->IsFake()) continue; // for the fake seeds fill only fake clusters part
1267       //
1268       ((TH2*)dest->At( GetHistoID(lrID,kHResY        ,fCurrPassID) ))->Fill(htrPt,sd->GetResidY());
1269       ((TH2*)dest->At( GetHistoID(lrID,kHResZ        ,fCurrPassID) ))->Fill(htrPt,sd->GetResidZ());
1270       ((TH2*)dest->At( GetHistoID(lrID,kHResYP       ,fCurrPassID) ))->Fill(htrPt,sd->GetPullY());
1271       ((TH2*)dest->At( GetHistoID(lrID,kHResZP       ,fCurrPassID) ))->Fill(htrPt,sd->GetPullZ());
1272       ((TH2*)dest->At( GetHistoID(lrID,kHChi2Cl      ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2Cl());
1273       //
1274     } while((sd=(AliITSUSeed*)sd->GetParent()));
1275     //
1276     ((TH2*)dest->At( GetHistoID(-1,kHChiMatch,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSTPC());
1277     ((TH2*)dest->At( GetHistoID(-1,kHChiITSSA,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSSA());
1278   }
1279   //
1280 #endif
1281   //
1282   CheckClusterSharingConflicts(hyp);
1283   return hyp->GetWinner();  // winner might change of disappear after resolving conflicts
1284   //
1285 }
1286
1287 //______________________________________________________________________________
1288 void AliITSUTrackerGlo::FinalizeHypotheses()
1289 {
1290   // select winner for each hypothesis, remove cl. sharing conflicts
1291   AliCodeTimerAuto("",0); 
1292   //
1293   int* esdTrackIndex = fESDIndex.GetArray();
1294   for (int ih=0;ih<fNTracksESD;ih++) FinalizeHypothesis(GetTrackHyp(esdTrackIndex[ih])); // finlize in decreasing pt order
1295   //
1296   for (int ih=fNTracksESD;ih--;) UpdateESDTrack(GetTrackHyp(esdTrackIndex[ih]),AliESDtrack::kITSin); 
1297   //
1298 }
1299
1300 //______________________________________________________________________________
1301 void AliITSUTrackerGlo::CheckClusterSharingConflicts(AliITSUTrackHyp* hyp)
1302 {
1303   // remove eventual cluster sharing conflicts
1304   AliITSUSeed* winner = 0;
1305   if (!(winner=hyp->GetWinner())) return;
1306   UShort_t idH = (UShort_t)hyp->GetUniqueID()+1;
1307   int lrID,clID;
1308   int inLr = fCurrTrackCond->GetActiveLrInner();
1309   AliITSUSeed* winSD=winner;
1310   do {
1311     if ( (clID=winSD->GetLrCluster(lrID))<0 ) continue;
1312     AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1313     Int_t refID = cl->GetRecoInfo();                // was it already referred by some track (with refID-1 if >0)
1314     if (!refID) {cl->SetRecoInfo(idH); continue;}   // if not, refer from track to cluster
1315     if (refID==idH) continue;                          // ignore reference to itself
1316     // the cluster is already used by other track, need to resolve the conflict
1317     AliITSUTrackHyp* hypC = GetTrackHyp(refID-1);   // competitor
1318     //    printf("Ref to %d (%p) from %d Cl %d %d\n",refID-1,hypC, idH-1,clID,lrID);
1319     AliITSUSeed* winnerC = hypC->GetWinner();
1320     if (!winnerC) {
1321       printf("Missing winner, ERROR\n");
1322       continue;
1323     }
1324     //
1325 #ifdef _ITSU_DEBUG_
1326     if (AliDebugLevelClass()>1) { //------------------------ DEBUG -------------------
1327       AliInfo(Form("Shared cl#%4d on Lr:%d: Hyp%3d/Q:%6.3f vs Hyp%3d/Q:%6.3f",
1328                    clID,lrID,idH-1,winner->GetQualityVar(),refID-1,winnerC->GetQualityVar()));
1329       // dump winner of hyp
1330       printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1331              idH-1,winner->Pt(),winner->GetChi2GloNrm(),winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),
1332              winner->GetNLayersHit(),hyp->GetITSLabel(),hyp->GetESDTrack()->GetTPCLabel());
1333       int prevL=-1;
1334       AliITSUSeed* sd = winner;
1335       do {
1336         int lrs;
1337         int clIDt = sd->GetLrCluster(lrs);
1338         if (clIDt<0) continue;
1339         while( lrs>++prevL ) printf("%4s        ","----");
1340         printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1341       } while ((sd=(AliITSUSeed*)sd->GetParent()));
1342       printf("|\n");
1343       // dump winner of hypC
1344       printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1345              refID-1,winnerC->Pt(),winnerC->GetChi2GloNrm(),winnerC->GetChi2ITSTPC(),winnerC->GetChi2ITSSA(),
1346              winnerC->GetNLayersHit(),hypC->GetITSLabel(),hypC->GetESDTrack()->GetTPCLabel());
1347       prevL=-1;
1348       sd = winnerC;
1349       do {
1350         int lrs;
1351         int clIDt = sd->GetLrCluster(lrs);
1352         if (clIDt<0) continue;
1353         while( lrs>++prevL ) printf("%4s        ","----");
1354         printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1355       } while ((sd=(AliITSUSeed*)sd->GetParent()));
1356       printf("|\n");
1357     }  //-------------------------------------------------- DEBUG -------------------
1358 #endif
1359     //
1360     if (winner->GetQualityVar()<winnerC->GetQualityVar()) { // current tracks is better then competitor track
1361       FlagSeedClusters(winnerC,kFALSE,refID); // unflag cluster usage by loser
1362       cl->SetRecoInfo(idH);
1363       //
1364       //if ( hypC->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",refID-1));
1365       winnerC->Kill();
1366       if (hypC->DefineWinner(inLr)) {             // find new winner instead of suppressed candidate
1367         CookMCLabel(hypC);
1368         CheckClusterSharingConflicts(hypC);   // and check its sharing conflicts again
1369         if (winner->IsKilled()) break;        // the current winner might have been killed during check of new winner of hypC!
1370       }
1371     }    
1372     else { // competitor hypC is better than the hyp
1373       FlagSeedClusters(winner,kFALSE,idH); // unflag cluster usage by loser
1374       winner->Kill();
1375       //if ( hyp->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",idH-1));
1376       if (hyp->DefineWinner(inLr)) {
1377         CookMCLabel(hyp);
1378         CheckClusterSharingConflicts(hyp);
1379       }
1380       break;
1381     }
1382     //   
1383   } while ((winSD=(AliITSUSeed*)winSD->GetParent()));
1384   //
1385   return;
1386 }
1387
1388 //______________________________________________________________________________
1389 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1390 {
1391   // update ESD track with current best hypothesis
1392   if (!hyp) return;
1393   AliESDtrack* esdTr = hyp->GetESDTrack();
1394   if (!esdTr || esdTr->IsOn(flag)) return;
1395   AliITSUSeed* win = hyp->GetWinner();
1396   if (!win || win->IsKilled()) return;
1397   double chiSav;
1398   //
1399   switch (flag) {
1400   case AliESDtrack::kITSin: 
1401     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1402     fCountITSin++;
1403     // set fakes cluster info
1404     {
1405       UShort_t clfake = 0;
1406       Int_t    clSplit = 0;
1407       AliITSUSeed* sd = win;
1408       int ip=0;
1409       do {
1410         int lr, clID = sd->GetLrCluster(lr);
1411         if (sd->IsFake()) clfake |= 0x1<<lr;
1412         if (clID>=0) {
1413           esdTr->SetITSModuleIndex(ip++, sd->GetLrClusterID());
1414           AliITSUClusterPix *cl = (AliITSUClusterPix*)fITS->GetLayerActive(lr)->GetCluster(clID);
1415 #ifdef  _ITSU_TUNING_MODE_
1416           if (cl->IsSplit()) clSplit |= 0x1<<lr;
1417 #endif
1418         }
1419       } while ((sd=(AliITSUSeed*)sd->GetParent()));
1420       //
1421       // RS: Temporary set special flag for tracks from the afterburner
1422       if (fCurrPassID>0) clfake |= 0x1<<7;
1423       //
1424       esdTr->SetITSSharedMap(clfake);
1425       esdTr->SetITSModuleIndex(10,clSplit);
1426     }
1427     // TEMPORARY: store iteration id
1428     esdTr->SetITSModuleIndex(11,fCurrPassID);
1429     break;
1430     //
1431   case AliESDtrack::kITSout: 
1432     // here the stored chi2 will correspond to backward ITS-SA fit
1433     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1434     fCountITSout++;
1435     // TODO: avoid setting friend
1436     break;
1437     //
1438   case AliESDtrack::kITSrefit: 
1439
1440     // at the moment fill as a chi2 the TPC-ITS matching chi2
1441     chiSav = hyp->GetChi2();
1442     hyp->SetChi2(win->GetChi2ITSTPC());
1443     esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1444     hyp->SetChi2(chiSav);
1445     fCountITSrefit++;
1446     // TODO: avoid setting cluster info
1447     break;
1448   default:
1449     AliFatal(Form("Unknown flag %d",flag));
1450   }
1451   //
1452   esdTr->SetITSLabel(hyp->GetITSLabel());
1453   // transfer chip indices
1454   // TODO
1455 }
1456
1457 //______________________________________________________________________________
1458 Double_t AliITSUTrackerGlo::RefitTrack(AliITSUTrackHyp* trc, Double_t rDest, Int_t &nclFit, Int_t stopCond)
1459 {
1460   // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
1461   // if stopCond<0 : propagate till last cluster then stop
1462   // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
1463   // if stopCond>0 : rDest must be reached
1464   //
1465   double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1466   int dir,lrStart,lrStop;
1467   //
1468   dir = rCurr<rDest ? 1 : -1;
1469   lrStart = fITS->FindFirstLayerID(rCurr,dir);
1470   lrStop  = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
1471   //
1472 #ifdef _ITSU_DEBUG_
1473   if (AliDebugLevelClass()>2) {
1474     printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1475     printf("Before refit: "); trc->Print();
1476   }
1477 #endif
1478   if (lrStop<0 || lrStart<0) AliFatal(Form("Failed to find start(%d) or last(%d) layers. Track from %.3f to %.3f",lrStart,lrStop,rCurr,rDest));
1479   //
1480   Int_t clInfo[2*AliITSUAux::kMaxLayers];
1481   Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
1482   fCurrMass = fCurrHyp->GetMass();
1483   AliITSUTrackHyp tmpTr(*(AliKalmanTrack*)trc);
1484   double chi2 = 0;
1485   int iclLr[2],nclLr;
1486   nclFit = 0;
1487   //
1488   int lrStop1 = lrStop+dir;
1489   for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
1490     AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1491     if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1492     int ilrA2,ilrA = lr->GetActiveID();
1493     // passive layer or active w/o hits will be traversed on the way to next cluster
1494     if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue; 
1495     //
1496     // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1497     nclLr=0;
1498     if (dir>0) { // clusters are stored in increasing radius order
1499       iclLr[nclLr++]=clInfo[ilrA2++];
1500       if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
1501     }
1502     else {
1503       if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
1504       iclLr[nclLr++]=clInfo[ilrA2];
1505     }
1506     //
1507     Bool_t transportedToLayer = kFALSE;
1508     for (int icl=0;icl<nclLr;icl++) {
1509       AliITSUClusterPix* clus =  (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1510       AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1511       if (!tmpTr.Rotate(sens->GetPhiTF())) {
1512 #ifdef _ITSU_DEBUG_
1513         AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1514         AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
1515 #endif
1516         return -1;
1517       }
1518       //
1519       double xClus = sens->GetXTF()+clus->GetX();
1520       if (!transportedToLayer) {
1521         if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1522 #ifdef _ITSU_DEBUG_
1523           AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1524           AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
1525 #endif
1526           return -1; // go to the entrance to the layer
1527         }
1528         lrStart = ilr;
1529         transportedToLayer = kTRUE;
1530       }
1531       //
1532 #ifdef _ITSU_DEBUG_
1533       if (AliDebugLevelClass()>1) {
1534         AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1535       }
1536 #endif
1537       //
1538       if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1539 #ifdef _ITSU_DEBUG_
1540         AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1541         AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
1542 #endif
1543         return -1;
1544       }
1545       //
1546       Double_t p[2]={clus->GetY(), clus->GetZ()};
1547       Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
1548       double chi2cl = tmpTr.GetPredictedChi2(p,cov);
1549       chi2 += chi2cl;
1550       //
1551 #ifdef  _ITSU_TUNING_MODE_
1552       TObjArray* dest = trc->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1553       if (dest && fTrackPhaseID>kClus2Tracks) {
1554         //
1555         double htrPt = tmpTr.Pt();
1556         double dy = p[0]-tmpTr.GetY();
1557         double dz = p[1]-tmpTr.GetZ();
1558         ((TH2*)dest->At( GetHistoID(ilrA,kHResY,0,fTrackPhaseID) ))->Fill(htrPt,dy);
1559         ((TH2*)dest->At( GetHistoID(ilrA,kHResZ,0,fTrackPhaseID) ))->Fill(htrPt,dz);
1560         double errY = tmpTr.GetSigmaY2();
1561         double errZ = tmpTr.GetSigmaZ2();
1562         if (errY>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResYP,0,fTrackPhaseID) ))->Fill(htrPt,dy/Sqrt(errY));
1563         if (errZ>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResZP,0,fTrackPhaseID) ))->Fill(htrPt,dz/Sqrt(errZ));
1564         ((TH2*)dest->At( GetHistoID(ilrA,kHChi2Cl,0,fTrackPhaseID) ))->Fill(htrPt,chi2cl);
1565       }
1566 #endif  
1567       //      
1568       if ( !tmpTr.Update(p,cov) ) {
1569 #ifdef _ITSU_DEBUG_
1570         AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1571         AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
1572 #endif
1573         return -1;
1574       }
1575 #ifdef _ITSU_DEBUG_
1576       if (AliDebugLevelClass()>1) {
1577         printf("AfterRefit: "); tmpTr.Print();
1578       }
1579 #endif
1580       if (++nclFit==nCl && stopCond<0) {
1581         *trc = (AliKalmanTrack&)tmpTr;
1582         return chi2; // it was requested to not propagate after last update
1583       }
1584     }
1585     //
1586   }
1587   // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1588   // Still, try to go as close as possible to rDest.
1589   //
1590   if (lrStart!=lrStop) {
1591     if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1592 #ifdef _ITSU_DEBUG_
1593       AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1594 #endif
1595       return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1596     }    
1597     if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1598 #ifdef _ITSU_DEBUG_
1599       AliDebug(1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1600 #endif
1601       return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1602     }
1603   }
1604   // go to the destination radius. Note that here we don't select direction to avoid precision problems
1605   if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
1606     return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1607   }
1608   *trc = (AliKalmanTrack&)tmpTr;
1609 #ifdef _ITSU_DEBUG_
1610   if (AliDebugLevelClass()>2) {
1611     printf("After refit (now at lr %d): ",lrStop); trc->Print();
1612   }
1613 #endif
1614   return chi2;
1615 }
1616
1617 //__________________________________________________________________
1618 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp) 
1619 {
1620   // build MC label
1621   //
1622   const int kMaxLbPerCl = 3;
1623   int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1624   Int_t lr,nLab=0,nCl=0;
1625   AliITSUSeed *seed = hyp->GetWinner();
1626   while(seed) {
1627     int clID = seed->GetLrCluster(lr);
1628     if (clID>=0) {
1629       AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1630       nCl++;
1631       for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1632         int trLb = cl->GetLabel(imc);
1633         if (trLb<0) break;
1634         // search this mc track in already accounted ones
1635         int iLab;
1636         for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1637         if (iLab<nLab) lbStat[iLab]++;
1638         else {
1639           lbID[nLab] = trLb;
1640           lbStat[nLab++] = 1;
1641         }
1642       } // loop over given cluster's labels
1643     }
1644     seed = (AliITSUSeed*)seed->GetParent();
1645   } // loop over clusters
1646   // 
1647   AliESDtrack* esdTr = hyp->GetESDTrack();
1648   int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1649   if (nCl && nLab) {
1650     int maxLab=0,nTPCok=0;
1651     for (int ilb=nLab;ilb--;) {
1652       int st = lbStat[ilb];
1653       if (lbStat[maxLab]<st) maxLab=ilb;
1654       if (tpcLab==lbID[ilb]) nTPCok += st;
1655     }
1656     hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1657     hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1658     hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
1659     return;
1660   }
1661   //
1662   hyp->SetFakeRatio(-1.);
1663   hyp->SetLabel( -tpcLab );
1664   hyp->SetITSLabel( kDummyLabel );
1665 }
1666
1667 //__________________________________________________________________
1668 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1669 {
1670   // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1671   if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
1672     AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
1673     fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
1674     AliITSUSeed** tmpArr = fLayerCandidates;
1675     fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
1676     memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
1677     delete[] tmpArr; // delete only array, not objects
1678   }
1679   AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1680   int slot=fNBranchesAdded++;
1681   for (int slotF=slot;slotF--;) { // slotF is always slot-1
1682     AliITSUSeed* si = branches[slotF];
1683     if (si->Compare(seed)<0) break; // found the last seed with better quality
1684     // otherwise, shift the worse seed to the next slot
1685     branches[slot] = si;
1686     slot = slotF; // slot should be slotF+1
1687   }
1688   // if needed, move worse seeds
1689   branches[slot] = seed;
1690   return kTRUE;
1691   //
1692 }
1693
1694 //__________________________________________________________________
1695 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1696 {
1697   // keep allowed number of branches for current seed and disable extras
1698   AliCodeTimerAuto("",0); 
1699   int nb = Min(fNBranchesAdded,acceptMax);
1700   //  if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1701   // disable unused branches
1702   AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1703   //
1704 #ifdef  _ITSU_TUNING_MODE_
1705   for (int ib=0;ib<fNBranchesAdded;ib++) branches[ib]->SetOrdBranch(ib);
1706 #endif
1707   //
1708   for (int ib=nb;ib<fNBranchesAdded;ib++) MarkSeedFree(branches[ib]);
1709   fNCandidatesAdded += nb; // update total candidates counter
1710   fNBranchesAdded = 0; // reset branches counter
1711   //
1712 }
1713
1714 //__________________________________________________________________
1715 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1716 {
1717   // transfer allowed number of branches to hypothesis container
1718   //
1719   AliCodeTimerAuto("",0); 
1720   // sort candidates in increasing order of chi2
1721   static int lastSize = 0;
1722   static int *index = 0;
1723   static float *chi2 = 0;
1724   if (fLayerMaxCandidates>lastSize) {
1725     lastSize = fLayerMaxCandidates;
1726     delete[] index;
1727     delete[] chi2;
1728     index = new int[lastSize];
1729     chi2  = new float[lastSize];
1730   }
1731   for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1732   Sort(fNCandidatesAdded,chi2,index,kFALSE);
1733   //
1734   int nacc=0,nb=0;
1735   if (ilr>fCurrTrackCond->GetActiveLrInner()) { // just take 1st acceptMax candidates
1736     nb = Min(fNCandidatesAdded,acceptMax);
1737     for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1738   }
1739   else { // for innermost layer test ITS SA backward chi2 and matching to TPC
1740     AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
1741     for (nb=0;nb<fNCandidatesAdded;nb++) {
1742       AliITSUSeed* sd = fLayerCandidates[index[nb]];    
1743       fCurrHyp->SetWinner(sd);
1744       if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
1745       else {
1746         AddProlongationHypothesis(sd,ilr);
1747         if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
1748       }
1749     }
1750     fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
1751   }
1752   // disable unused candidates
1753   for (int ib=nb;ib<fNCandidatesAdded;ib++) MarkSeedFree(fLayerCandidates[index[ib]]);    
1754   fNCandidatesAdded = 0; // reset candidates counter
1755   //
1756 }
1757
1758 //__________________________________________________________________
1759 Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
1760 {
1761   // check seed backward propagation chi2 and matching to TPC 
1762   double bz0 = GetBz();
1763   double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
1764   AliITSUTrackHyp trback;
1765   trback.AliExternalTrackParam::operator=(*seed);
1766   trback.ResetCovariance(10000);
1767   int nclFit = 0;
1768   double chi2sa = RefitTrack(&trback,rDest,nclFit,1);
1769   if (chi2sa<0) return kFALSE;
1770   int ndf = nclFit*2-5;
1771   if (ndf>0) chi2sa /= ndf;
1772   if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2(nclFit)) return kFALSE;
1773   //
1774   // relate to TPC track at outer layer
1775   AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
1776   if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
1777 #ifdef _ITSU_DEBUG_
1778     if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
1779       AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
1780       seed->Print("etp");
1781       trback.Print();
1782     }
1783 #endif
1784     return kFALSE; 
1785   }
1786   double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
1787   if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
1788   //
1789   seed->SetChi2ITSSA(chi2sa);
1790   seed->SetChi2ITSTPC(chi2Match);
1791   return kTRUE;
1792 }
1793
1794 //__________________________________________________________________
1795 void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
1796 {
1797   // remove those hypothesis seeds which dont lead to candidates at final layer
1798   //
1799   // 1st, flag the seeds to save
1800   int lr0 = fCurrTrackCond->GetActiveLrInner();
1801   for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
1802     AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
1803     if (!seed) continue;
1804     seed->FlagTree(AliITSUSeed::kSave);
1805     dest->AddSeed(seed,lr0);
1806   }
1807   for (int ilr=0;ilr<fNLrActive;ilr++) {
1808     int nsd = fCurrHyp->GetNSeeds(ilr);
1809     for (int isd=0;isd<nsd;isd++) {
1810       AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
1811       if (!seed) continue; // already discarded or saved
1812       if (seed->IsSaved()) dest->AddSeed(seed,ilr);
1813       else MarkSeedFree(seed);
1814     }
1815   }
1816   //
1817   //  AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
1818 }
1819
1820 //__________________________________________________________________
1821 void AliITSUTrackerGlo::CleanHypothesis(AliITSUTrackHyp* hyp)
1822 {
1823   // clean hypothesis
1824   hyp->SetWinner(0);
1825   for (int ilr=0;ilr<fNLrActive;ilr++) {
1826     int nsd = hyp->GetNSeeds(ilr);
1827     for (int isd=0;isd<nsd;isd++) {
1828       AliITSUSeed* seed = hyp->RemoveSeed(ilr,isd);
1829       if (!seed) continue; // already discarded or saved
1830       MarkSeedFree(seed);
1831     }
1832   }
1833 }
1834
1835 //__________________________________________________________________
1836 void AliITSUTrackerGlo::FlagSplitClusters()
1837 {
1838   // set special bit on split clusters using MC info
1839   for (int ilr=fNLrActive;ilr--;) {
1840     int nsplit=0;
1841     AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1842     for (int isn=lr->GetNSensors();isn--;) {
1843       AliITSURecoSens* sens = lr->GetSensor(isn);
1844       int nCl = sens->GetNClusters();
1845       if (!nCl) continue;
1846       int cl0 = sens->GetFirstClusterId();
1847       for (int icl=nCl;icl--;) {
1848         AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1849         for (int icl1=icl;icl1--;) {
1850           AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1851           if (cl->HasCommonTrack(cl1)) {
1852             if (!cl->IsSplit())  nsplit++;
1853             if (!cl1->IsSplit()) nsplit++;
1854             cl->SetSplit();
1855             cl1->SetSplit();
1856           }
1857         }
1858       }
1859     }
1860     AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1861   }
1862   //
1863 }
1864
1865 //__________________________________________________________________
1866 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1867 {
1868   // check if the seed contains split cluster with size < maxSize
1869   int lrID,clID;
1870   if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1871     AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1872     if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1873   }
1874   return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1875   //
1876 }
1877
1878 //__________________________________________________________________
1879 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1880 {
1881   // print seeds clusters
1882   int lrID,clID;
1883   if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1884     AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1885     cl->Print(option);
1886   }
1887   if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1888   //
1889 }
1890
1891 //__________________________________________________________________
1892 void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg, UShort_t ref)
1893 {
1894   // mark used clusters
1895   int lrID,clID;
1896   while (seed) {
1897     if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1898       AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1899       if (ref) { // do we need to set or delete cluster->track ref?
1900         if (flg) {
1901           if (!cl->GetRecoInfo()) cl->SetRecoInfo(ref);  // set ref only if cluster already does not refer to other track, inc.counter
1902         }
1903         else {
1904           if (cl->GetRecoInfo()==ref) cl->SetRecoInfo(0); // unset reference only if it refers to ref, decrease counter     
1905         }
1906       }
1907     }
1908     seed = (AliITSUSeed*)seed->GetParent();
1909   }
1910   //
1911 }
1912
1913 //__________________________________________________________________
1914 void AliITSUTrackerGlo::CheckClusterUsage()
1915 {
1916   // check cluster usage info
1917   printf("ClusterUsage at pass %d\n",fCurrPassID);
1918   for (int ilr=0;ilr<fNLrActive;ilr++) {
1919     AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1920     int ncl = lr->GetNClusters();
1921     int nclUsed = 0;
1922     int nclUs0CntShr1 = 0;
1923     int nclUs1Cnt0 = 0;
1924     int nclUs1Shr1 = 0;
1925     int nclUseL = 0;
1926     for (int icl=0;icl<ncl;icl++) {
1927       AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(icl);
1928       int usc = cl->GetClUsage();
1929       if (usc<1) {
1930         if (cl->IsClusterUsed() || cl->IsClusterShared()) {
1931           nclUs0CntShr1++;
1932           printf("Lr%d Cl%4d: UseCnt=0 but UseFlg=%d UseShr=%d\n",ilr,icl,cl->IsClusterUsed(),cl->IsClusterShared());
1933         }
1934         continue;
1935       }
1936       nclUsed++;
1937       if (!cl->IsClusterUsed()) {
1938         nclUs1Cnt0++;
1939         printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1940       }
1941       if (cl->IsClusterShared()) {
1942         nclUs1Shr1++;
1943         printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1944       }
1945       if (usc>1) {
1946         nclUseL++;
1947         printf("Lr%d Cl%4d: UseCnt=%d!! but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1948       }
1949     }
1950     printf("ClUsage Lr%d: Ncl=%5d Nusd=%5d (%5.2f) Use0CS1:%4d Use1C0:%4d Use1S1:%4d UseL:%d\n",
1951            ilr,ncl,nclUsed,ncl? float(nclUsed)/ncl : 0, nclUs0CntShr1,nclUs1Cnt0,nclUs1Shr1,nclUseL);
1952   }
1953 }
1954
1955
1956 #ifdef  _ITSU_TUNING_MODE_
1957 //__________________________________________________________________
1958 TObjArray* AliITSUTrackerGlo::BookControlHistos(const char* pref)
1959 {
1960   // book special control histos
1961   TString prefS = pref;
1962   prefS.ToLower();
1963   TObjArray* dest = new TObjArray;
1964   dest->SetOwner(kTRUE);
1965   //
1966   const int kNResDef=7;
1967   const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1968   const double ptMax=10;
1969   const double plMax=10;
1970   const double chiMax=100;
1971   const int nptbins=50;
1972   const int nresbins=400;
1973   const int nplbins=50;
1974   const int nchbins=200;
1975   const int maxBr  = 15;
1976   const int maxCand = 200;
1977   TString ttl;
1978   for (int phase=0;phase<kNTrackingPhases;phase++) {
1979     for (int pass=0;pass<AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();pass++) {
1980       //
1981       for (int ilr=0;ilr<fNLrActive;ilr++) {
1982         //
1983         // ----------------- These are histos to be filled in Cluster2Tracks of each pass. 
1984         // PropagateBack and RefitInward will be stored among the histos of 1st pass
1985         if (pass>0 && phase!=kClus2Tracks) continue;
1986         //
1987         double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1988         ttl = Form("Pass%d_S%d_residY%d_%s",pass,phase,ilr,pref);
1989         TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1990         dest->AddAtAndExpand(hdy,GetHistoID(ilr,kHResY,pass,phase));
1991         hdy->SetDirectory(0);
1992         //
1993         ttl = Form("Pass%d_S%d_residYPull%d_%s",pass,phase,ilr,pref);   
1994         TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1995         dest->AddAtAndExpand(hdyp,GetHistoID(ilr,kHResYP,pass,phase));
1996         hdyp->SetDirectory(0);
1997         //
1998         ttl = Form("Pass%d_S%d_residZ%d_%s",pass,phase,ilr,pref);       
1999         TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
2000         dest->AddAtAndExpand(hdz,GetHistoID(ilr,kHResZ,pass,phase));
2001         hdz->SetDirectory(0);
2002         //
2003         ttl = Form("Pass%d_S%d_residZPull%d_%s",pass,phase,ilr,pref);           
2004         TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
2005         hdzp->SetDirectory(0);
2006         dest->AddAtAndExpand(hdzp,GetHistoID(ilr,kHResZP,pass,phase));
2007         //
2008         ttl = Form("Pass%d_S%d_chi2Cl%d_%s",pass,phase,ilr,pref);               
2009         TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2010         hchi->SetDirectory(0);
2011         dest->AddAtAndExpand(hchi,GetHistoID(ilr,kHChi2Cl,pass,phase));
2012         //
2013         // ------------------- These histos are filled for Clusters2Tracks only
2014         if (phase!=kClus2Tracks) continue;
2015         //
2016         ttl = Form("Pass%d_S%d_chi2Nrm%d_%s",pass,phase,ilr,pref);              
2017         TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2018         hchiN->SetDirectory(0);
2019         dest->AddAtAndExpand(hchiN,GetHistoID(ilr,kHChi2Nrm,pass,phase));
2020         //
2021         ttl = Form("Pass%d_S%d_bestInBranch%d_%s",pass,phase,ilr,pref);         
2022         TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
2023         hnbr->SetDirectory(0);
2024         dest->AddAtAndExpand(hnbr,GetHistoID(ilr,kHBestInBranch,pass,phase));
2025         //
2026         ttl = Form("Pass%d_S%d_bestInCands%d_%s",pass,phase,ilr,pref);          
2027         TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
2028         hncn->SetDirectory(0);
2029         dest->AddAtAndExpand(hncn,GetHistoID(ilr,kHBestInCand,pass,phase));
2030         //
2031       } // loop over layers
2032       //
2033       // these are special histos, filled not per layer but in the end of track fit in Clusters2Tracks in EVERY pass
2034       //  
2035       if (phase==kClus2Tracks) {
2036         TH2* hchiMatch = 0; 
2037         ttl = Form("Pass%d_Chi2Match_%s",pass,pref);
2038         hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2039         hchiMatch->SetDirectory(0);
2040         dest->AddAtAndExpand(hchiMatch,GetHistoID(-1,kHChiMatch,pass,phase));
2041         // 
2042         TH2* hchiSA = 0; 
2043         ttl = Form("Pass%d_Chi2ITSSA_%s",pass,pref);
2044         hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
2045         hchiSA->SetDirectory(0);
2046         dest->AddAtAndExpand(hchiSA,GetHistoID(-1,kHChiITSSA,pass,phase));
2047         //      
2048       }
2049     } // loop over tracking passes
2050   }// loop over tracking phases
2051   //
2052   return dest;
2053 }
2054
2055 //__________________________________________________________________
2056 Int_t AliITSUTrackerGlo::GetHistoID(Int_t lr, Int_t hid, Int_t pass, Int_t phase)
2057 {
2058   // get id for the requested histo
2059   if (lr<0) lr=-1;
2060   lr++;
2061   return pass*kHistosPass + phase*kHistosPhase + lr*kMaxHID + hid;
2062   //
2063 }
2064
2065 #endif