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