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