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