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