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