1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
22 //-------------------------------------------------------------------------
23 // Implementation of the ITS Upgrade tracker mother class.
24 //-------------------------------------------------------------------------
26 #include <Riostream.h>
30 #include "AliITSUTrackerGlo.h"
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
33 #include "AliITSURecoDet.h"
34 #include "AliITSURecoSens.h"
35 #include "AliITSUReconstructor.h"
36 #include "AliITSReconstructor.h"
37 #include "AliITSUSeed.h"
38 #include "AliITSUClusterPix.h"
39 #include "AliITSUGeomTGeo.h"
40 #include "AliCodeTimer.h"
41 #include "AliRefArray.h"
42 using namespace AliITSUAux;
43 using namespace TMath;
45 //----------------- tmp stuff -----------------
47 ClassImp(AliITSUTrackerGlo)
49 const Double_t AliITSUTrackerGlo::fgkToler = 1e-6;// tolerance for layer on-surface check
52 //_________________________________________________________________________
53 AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
57 ,fCurrESDtrMClb(kDummyLabel)
59 ,fCountProlongationTrials(0)
65 ,fLayerMaxCandidates(1000)
71 ,fSeedsPool("AliITSUSeed",0)
83 #ifdef _ITSU_TUNING_MODE_
88 // Default constructor
92 //_________________________________________________________________________
93 AliITSUTrackerGlo::~AliITSUTrackerGlo()
97 delete[] fLayerCandidates;
98 if (fWorkHyp) fWorkHyp->SetTPCSeed(0); // this hypothesis does not own the seed
101 #ifdef _ITSU_TUNING_MODE_
102 if (fCHistoArrCorr || fCHistoArrFake) {
103 TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
105 AliInfo("Storing control histos");
106 // ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
107 if (fCHistoArrCorr) {
108 fCHistoArrCorr->Write();
109 delete fCHistoArrCorr;
111 if (fCHistoArrFake) {
112 fCHistoArrFake->Write();
113 delete fCHistoArrFake;
124 //_________________________________________________________________________
125 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
127 // init with external reconstructor
129 fITS = rec->GetITSInterface();
130 fNLrActive = fITS->GetNLayersActive();
131 fWorkHyp = new AliITSUTrackHyp(fNLrActive);
133 if (fLayerMaxCandidates<1) fLayerMaxCandidates = 1000;
134 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
135 fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
136 fFreeSeedsID.Set(1000);
142 //_________________________________________________________________________
143 void AliITSUTrackerGlo::CreateDefaultTrackCond()
145 // creates default tracking conditions to be used when recoparam does not provide them
147 AliITSUTrackCond* cond = new AliITSUTrackCond();
149 cond->SetNLayers(fNLrActive);
150 cond->AddNewCondition(fNLrActive);
151 cond->AddGroupPattern( 0xffff, fNLrActive ); // require all layers hit
154 fDefTrackConds.AddLast(cond);
156 AliInfo("Created conditions: ");
157 for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
162 //_________________________________________________________________________
163 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
166 AliCodeTimerAuto("",0);
167 SetTrackingPhase(kClus2Tracks);
169 #ifdef _ITSU_TUNING_MODE_
170 if (!fCHistoArrCorr) fCHistoArrCorr = BookControlHistos("Corr");
171 if (!fCHistoArrFake) fCHistoArrFake = BookControlHistos("Fake");
174 static TArrayF esdTrPt(fESDIndex.GetSize());
176 TObjArray *trackConds = 0;
178 fCountProlongationTrials = 0;
184 fNTracksESD = esdEv->GetNumberOfTracks();
185 AliInfo(Form("Will try to find prolongations for %d tracks",fNTracksESD));
186 int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
188 if (!fDefTrackConds.GetEntriesFast()) {
189 AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
190 CreateDefaultTrackCond();
192 trackConds = &fDefTrackConds;
193 nTrackCond = trackConds->GetEntriesFast();
195 else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
197 static Bool_t first = kTRUE;
199 AliITSUReconstructor::GetRecoParam()->Print();
203 if (fHypStore.GetSize()<fNTracksESD) fHypStore.Expand(fNTracksESD+100);
205 fITS->ProcessClusters();
207 #ifdef _ITSU_TUNING_MODE_
208 FlagSplitClusters(); // tmp RS
211 // the tracks will be reconstructed in decreasing pt order, sort them
212 if (fESDIndex.GetSize()<fNTracksESD) {
213 fESDIndex.Set(fNTracksESD+200);
214 esdTrPt.Set(fNTracksESD+200);
216 int* esdTrackIndex = fESDIndex.GetArray();
217 float* trPt = esdTrPt.GetArray();
218 for (int itr=fNTracksESD;itr--;) {
219 AliESDtrack* esdTr = esdEv->GetTrack(itr);
220 esdTr->SetStatus(AliESDtrack::kITSupg); // Flag all tracks as participating in the ITSup reco
221 trPt[itr] = esdTr->Pt();
223 Sort(fNTracksESD,trPt,esdTrackIndex,kTRUE);
225 for (int icnd=0;icnd<nTrackCond;icnd++) {
227 fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
228 if (!fCurrTrackCond->IsInitDone()) fCurrTrackCond->Init();
229 // select ESD tracks to propagate
230 for (int itr=0;itr<fNTracksESD;itr++) {
231 int trID = esdTrackIndex[itr];
232 fCurrESDtrack = esdEv->GetTrack(trID);
233 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
236 if (!NeedToProlong(fCurrESDtrack,trID)) continue; // are we interested in this track?
238 // if specific tracks should be checked, create a global array int wtc[] = {evITS*10000+trID, ...};
240 int nwtc = sizeof(wtc)/sizeof(int);
242 for (int k=0;k<nwtc;k++) {
243 if (wtc[k]==evID*10000+trID) {
245 printf("At esdTr: %d MC: %d\n",wtc[k],fCurrESDtrMClb);
249 AliLog::SetClassDebugLevel("AliITSUTrackerGlo",dbg ? 10:0);
251 AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMass(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
252 FindTrack(fCurrESDtrack, trID);
255 if (AliDebugLevelClass()>+2) {
256 AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
259 FinalizeHypotheses();
260 #ifdef _ITSU_TUNING_MODE_
261 CheckClusterUsage(); //!!RS
264 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0); // in case wtc array was defined, uncomment this
267 AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,fNTracksESD));
272 //_________________________________________________________________________
273 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
276 // Do outward fits in ITS
277 AliCodeTimerAuto("",0);
279 SetTrackingPhase(kPropBack);
280 fNTracksESD = esdEv->GetNumberOfTracks();
281 AliDebug(1,Form("Will propagate back %d tracks",fNTracksESD));
283 double bz0 = GetBz();
284 Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
285 AliITSUTrackHyp dummyTr;
286 const double kWatchStep=10.; // for larger steps watch arc vs segment difference
287 Double_t times[AliPID::kSPECIES];
289 for (int itr=0;itr<fNTracksESD;itr++) {
290 fCurrESDtrack = esdEv->GetTrack(itr);
291 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
292 // Start time integral and add distance from current position to vertex
293 if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
295 fCurrESDtrack->GetXYZ(xyzTrk);
296 Double_t dst = 0.; // set initial track lenght, tof
298 double dxs = xyzTrk[0] - xyzVtx[0];
299 double dys = xyzTrk[1] - xyzVtx[1];
300 double dzs = xyzTrk[2] - xyzVtx[2];
301 // RS: for large segment steps d use approximation of cicrular arc s by
302 // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
303 // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
304 // Hence s^2/d^2 = (1+1/6 p^2)^2
305 dst = dxs*dxs + dys*dys;
306 if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
307 double crv = Abs(fCurrESDtrack->GetC(bz0));
308 double fcarc = 1.+crv*crv*dst/6.;
315 fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
317 if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
318 dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
319 dummyTr.StartTimeIntegral();
320 dummyTr.AddTimeStep(dst);
321 dummyTr.GetIntegratedTimes(times);
322 fCurrESDtrack->SetIntegratedTimes(times);
323 fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
327 fCurrHyp = GetTrackHyp(itr);
328 fCurrMass = fCurrHyp->GetMass();
329 fCurrHyp->StartTimeIntegral();
330 fCurrHyp->AddTimeStep(dst);
331 fCurrHyp->ResetCovariance(10000);
333 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax(),nclFit);
334 if (chi2>0) { // propagate to exit from the ITS/TPC screen
335 int ndf = nclFit*2-5;
336 if (ndf>0) chi2 /= ndf;
337 fCurrHyp->SetChi2(chi2);
338 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSout);
341 AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
342 //fCurrHyp->AliExternalTrackParam::Print();
343 //fCurrHyp->GetWinner()->Print();
348 AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
349 fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
354 //_________________________________________________________________________
355 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
358 // refit the tracks inward, using current cov.matrix
359 AliCodeTimerAuto("",0);
361 SetTrackingPhase(kRefitInw);
362 fNTracksESD = esdEv->GetNumberOfTracks();
363 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
365 AliDebug(1,Form("Will refit inward %d tracks",fNTracksESD));
367 for (int itr=0;itr<fNTracksESD;itr++) {
368 fCurrESDtrack = esdEv->GetTrack(itr);
369 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
370 // Start time integral and add distance from current position to vertex
371 UInt_t trStat = fCurrESDtrack->GetStatus();
372 if ( !(trStat & AliESDtrack::kITSout) ) continue;
373 if ( trStat & AliESDtrack::kITSrefit ) continue; // already done
374 if ( (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
376 fCurrHyp = GetTrackHyp(itr);
377 fCurrHyp->AliExternalTrackParam::operator=(*fCurrESDtrack); // fetch current ESDtrack kinematics
378 fCurrMass = fCurrHyp->GetMass();
381 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin(),nclFit);
382 if (chi2>0) { // propagate up to inside radius of the beam pipe
383 int ndf = nclFit*2-5;
384 if (ndf>0) chi2 /= ndf;
385 fCurrHyp->SetChi2(chi2);
386 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
389 AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
393 AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
394 fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
396 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
400 //_________________________________________________________________________
401 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
403 // read from tree (if pointer provided) or directly from the ITS reco interface
405 return fReconstructor->LoadClusters(treeRP);
408 //_________________________________________________________________________
409 void AliITSUTrackerGlo::UnloadClusters()
415 Info("UnloadClusters","To be implemented");
417 //_________________________________________________________________________
418 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
424 Info("GetCluster","To be implemented");
428 //_________________________________________________________________________
429 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr, Int_t esdID)
431 // do we need to match this track to ITS?
433 static double bz = GetBz();
435 // is it already prolonged
436 if (esdTr->IsOn(AliESDtrack::kITSin)) return kFALSE;
437 AliITSUTrackHyp* hyp = GetTrackHyp(esdID); // in case there is unfinished hypothesis from prev.pass, clean it
439 CleanHypothesis(hyp);
440 if (hyp->GetSkip()) return kFALSE; // need to skip this
444 if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
445 esdTr->IsOn(AliESDtrack::kTPCout) ||
446 esdTr->IsOn(AliESDtrack::kITSin) ||
447 esdTr->GetKinkIndex(0)>0) return kFALSE;
449 if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
452 esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz);
453 // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
454 if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation())
455 && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
456 Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
458 // if (esdTr->Pt()<3) return kFALSE;//RS
462 //_________________________________________________________________________
463 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
465 // find prolongaion candidates finding for single seed
467 AliITSUSeed seedUC; // copy of the seed from the upper layer
468 AliITSUSeed seedT; // transient seed between the seedUC and new prolongation hypothesis
470 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
471 AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID); // initialize prolongations hypotheses tree
472 if (!hypTr || hypTr->GetSkip()) return;
476 fCurrHyp->InitFrom(hypTr);
478 AliITSURecoSens *hitSens[AliITSURecoSens::kNNeighbors+1];
480 int ilaUp = fNLrActive; // previous active layer really checked (some may be excluded!)
482 for (int ila=fNLrActive;ila--;ilaUp=fCurrActLrID) {
483 if (fCurrTrackCond->IsLayerExcluded(ila)) continue;
485 fCurrLayer = fITS->GetLayerActive(ila);
486 Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
488 // for the outermost layer the seed is created from the ESD track
489 int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
490 int maxNBranches = fCurrTrackCond->GetMaxBranches(ila);
491 int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
493 for (int isd=0;isd<nSeedsUp;isd++) {
495 if (ilaUp==fNLrActive) {
497 seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
500 seedU = fCurrHyp->GetSeed(ilaUp,isd); // seed on prev.active layer to prolong
501 if (seedU->IsKilled()) continue;
502 seedUC = *seedU; // its copy will be prolonged
503 seedUC.SetParent(seedU);
505 seedUC.ResetFMatrix(); // reset the matrix for propagation to next layer
506 // go till next active layer
507 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()));
508 if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
510 AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
511 // Check if the seed satisfies to track definition
512 if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
513 continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
515 if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
516 if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
522 if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
523 mcquest = fCurrESDtrMClb;
525 printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
528 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest); // find detectors which may be hit by the track
530 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens); // find detectors which may be hit by the track
531 AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
535 for (int isn=nsens;isn--;) {
536 AliITSURecoSens* sens = hitSens[isn];
537 int ncl = sens->GetNClusters();
541 // We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
542 // since the transport matrix should be defined in this frame.
543 double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
544 if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
545 if (AliDebugLevelClass()>2) {
546 printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
551 if (xs<seedT.GetX()) {
552 if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
554 else { // some low precision tracks may hit the sensor plane outside of the layer radius
555 if (AliDebugLevelClass()>2) {
556 if (!seedT.ContainsFake()) {
557 printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
562 if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
564 // if (!seedT.PropagateToX(xs,bz)) continue;
565 // if (!seedT.Rotate(sens->GetPhiTF())) continue;
566 if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
568 int clID0 = sens->GetFirstClusterId();
569 for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
570 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
571 int clID = clID0 + icl;
572 if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
573 int res = CheckCluster(&seedT,ila,clID0+icl);
574 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
576 if (res==kStopSearchOnSensor) break; // stop looking on this sensor
577 if (res==kClusterNotMatching) continue; // cluster does not match
578 // cluster is matching and it was added to the hypotheses tree
581 // cluster search is done. Do we need to have a version of this seed skipping current layer
582 if (!NeedToKill(&seedUC,kMissingCluster)) {
583 AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
584 if (nsens) { // to do: make penalty to account for probability to miss the cluster for good reason
585 double penalty = -fCurrTrackCond->GetMissPenalty(ila);
586 seedSkp->SetChi2Cl(penalty);
588 if (seedSkp->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(ila)) {
589 MarkSeedFree(seedSkp);
591 else AddSeedBranch(seedSkp);
593 // transfer the new branches of the seed to the hypothesis container
594 if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
597 if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
598 // ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
599 if (AliDebugLevelClass()>2) { //RS
600 printf(">>> All hypotheses on lr %d: \n",ila);
601 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
605 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
606 if (ila!=0) continue;
607 double vecL[5] = {0};
608 double matL[15] = {0};
609 AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
610 while(sp->GetParent()) {
611 sp->Smooth(vecL,matL);
612 if (sp->GetLayerID()>=fNLrActive-1) break;
613 sp = (AliITSUSeed*)sp->GetParent();
618 SaveReducedHypothesesTree(hypTr);
619 if (AliDebugLevelClass()>1) {
620 printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
624 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
628 //_________________________________________________________________________
629 AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
631 // init prolongaion candidates finding for single seed
632 AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
635 fCountProlongationTrials++;
637 fCurrMass = esdTr->GetMass();
638 if (fCurrMass<kPionMass*0.9) fCurrMass = kPionMass; // don't trust to mu, e identification from TPCin
640 hyp = new AliITSUTrackHyp(fNLrActive);
641 hyp->SetESDTrack(esdTr);
642 hyp->SetUniqueID(esdID);
643 hyp->SetMass(fCurrMass);
644 hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
645 SetTrackHyp(hyp,esdID);
647 if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
653 //_________________________________________________________________________
654 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
656 // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
659 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
661 int dir = lTo > lFrom ? 1:-1;
662 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
663 Bool_t checkFirst = kTRUE;
664 Bool_t limReached = kFALSE;
667 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
670 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
671 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
673 // go the entrance of the layer, assuming no materials in between
674 double xToGo = lrTo->GetR(-dir);
677 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
680 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
683 // double xts = xToGo;
684 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
685 //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
686 //seed->Print("etp");
690 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
691 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
692 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
693 //seed->Print("etp");
697 if (limReached) break;
703 //_________________________________________________________________________
704 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
706 // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
708 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
710 int dir = lTo > lFrom ? 1:-1;
711 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
712 Bool_t checkFirst = kTRUE;
713 Bool_t limReached = kFALSE;
716 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
719 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
720 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
722 // go the entrance of the layer, assuming no materials in between
723 double xToGo = lrTo->GetR(-dir);
726 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
729 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
732 // double xts = xToGo;
733 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
734 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
735 // seed->Print("etp");
738 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
739 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
740 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
741 //seed->Print("etp");
745 if (limReached) break;
751 //_________________________________________________________________________
752 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
754 // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
756 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
758 int dir = lTo > lFrom ? 1:-1;
759 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
760 Bool_t checkFirst = kTRUE;
763 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
766 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
767 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
769 // go the entrance of the layer, assuming no materials in between
770 double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
772 // double xts = xToGo;
773 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
774 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
775 // seed->Print("etp");
778 if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
780 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
781 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
782 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
783 //seed->Print("etp");
792 //_________________________________________________________________________
793 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
795 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
796 // If check is requested, do this only provided the track has not exited the layer already
797 double xToGo = lr->GetR(dir);
798 if (check) { // do we need to track till the surface of the current layer ?
799 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
800 AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
801 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
802 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
804 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
805 // go via layer to its boundary, applying material correction.
806 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
812 //_________________________________________________________________________
813 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
815 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
816 // If check is requested, do this only provided the track has not exited the layer already
817 double xToGo = lr->GetR(dir);
818 if (check) { // do we need to track till the surface of the current layer ?
819 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
820 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
821 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
823 AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
824 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
825 // go via layer to its boundary, applying material correction.
826 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
833 //_________________________________________________________________________
834 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
836 // go to the entrance of lr in direction dir, w/o applying material corrections.
837 // If check is requested, do this only provided the track did not reach the layer already
838 double xToGo = lr->GetR(-dir);
839 if (check) { // do we need to track till the surface of the current layer ?
840 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
841 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
842 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
844 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
845 // go via layer to its boundary, applying material correction.
846 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
851 //_________________________________________________________________________
852 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
854 // go to the entrance of lr in direction dir, w/o applying material corrections.
855 // If check is requested, do this only provided the track did not reach the layer already
856 double xToGo = lr->GetR(-dir);
857 if (check) { // do we need to track till the surface of the current layer ?
858 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
859 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
860 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
862 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
863 // go via layer to its boundary, applying material correction.
864 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
869 //_________________________________________________________________________
870 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
872 // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
873 // as well as some aux info
875 AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
876 seed->GetXYZ(&fTrImpData[kTrXIn]); // lab position at the entrance from above
877 static AliExternalTrackParam sc; // seed copy for manipulations
880 fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
881 if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
882 double dr = lrA->GetDR(); // approximate X dist at the inner radius
883 if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
884 // special case: track does not reach inner radius, might be tangential
885 double r = sc.GetD(0,0,bz);
887 if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
888 dr = Abs(sc.GetX() - x);
889 if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
892 fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
893 double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
894 double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
895 sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
896 sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
897 double dphi0 = 0.5*Abs(fTrImpData[kTrPhiOut]-fTrImpData[kTrPhiIn]);
898 double phi0 = 0.5*(fTrImpData[kTrPhiOut]+fTrImpData[kTrPhiIn]);
899 if ( dphi0>(0.5*Pi()) ) {
900 // special case of angles around pi
905 fTrImpData[kTrPhi0] = phi0;
906 fTrImpData[kTrZ0] = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
907 dphi0 += sgy/lrA->GetR();
908 fTrImpData[kTrDPhi] = dphi0<PiOver2() ? dphi0 : PiOver2();
909 fTrImpData[kTrDZ] = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn]) + sgz;
914 //________________________________________
915 void AliITSUTrackerGlo::ResetSeedsPool()
917 // mark all seeds in the pool as unused
918 AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
920 fSeedsPool.Clear(); // seeds don't allocate memory
924 //________________________________________
925 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd)
927 // account that this seed is "deleted"
928 int id = sd->GetPoolID();
930 AliError(Form("Freeing of seed %p NOT from the pool is requested",sd));
933 // AliInfo(Form("%d %p",id, seed));
934 fSeedsPool.RemoveAt(id);
935 if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
936 fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
939 //_________________________________________________________________________
940 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID)
942 // Check if the cluster (in tracking frame!) is matching to track.
943 // The track must be already propagated to sensor tracking frame.
944 // Returns: kStopSearchOnSensor if the search on given sensor should be stopped,
945 // kClusterMatching if the cluster is matching
946 // kClusterMatching otherwise
948 const double kTolerX = 5e-4;
949 // AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
950 AliCluster *cl = fCurrLayer->GetCluster(clID);
952 Bool_t goodCl = kFALSE;
953 int currLabel = Abs(fCurrESDtrMClb);
955 if (cl->GetLabel(0)>=0) {
956 for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
959 if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
960 if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
962 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
963 AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
967 return kStopSearchOnSensor; // propagation failed, seedT is intact
970 double dy = cl->GetY()-track->GetY();
971 double dz = cl->GetZ()-track->GetZ();
974 double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
975 fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
976 if (dy2>tol2) { // the clusters are sorted in Z(col) then in Y(row).
978 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
979 AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
982 PrintSeedClusters(track);
984 // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
985 // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
986 if (dy<0) return kStopSearchOnSensor; // No chance that other cluster of this sensor will match (all Y's will be even larger)
987 else return kClusterNotMatching; // Other clusters may match
990 tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
991 fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
993 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
994 AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
997 PrintSeedClusters(track);
999 return kClusterNotMatching; // Other clusters may match
1003 Double_t p[2]={cl->GetY(), cl->GetZ()};
1004 Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
1005 double chi2 = track->GetPredictedChi2(p,cov);
1007 if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
1008 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1009 AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1010 lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1011 track->Print("etp");
1013 PrintSeedClusters(track);
1015 return kClusterNotMatching;
1018 track = NewSeedFromPool(track); // input track will be reused, use its clone for updates
1019 if (!track->Update()) {
1020 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1021 AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
1022 track->Print("etp");
1024 PrintSeedClusters(track);
1026 MarkSeedFree(track);
1027 return kClusterNotMatching;
1029 track->SetChi2Cl(chi2);
1030 track->SetLrClusterID(lr,clID);
1032 if (track->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(lr)) {
1033 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1034 AliDebug(2,Form("Lost good cl on L:%d : Chi2Glo=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1035 lr,track->GetChi2GloNrm(),fCurrTrackCond->GetMaxChi2GloNrm(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1036 track->Print("etp");
1038 PrintSeedClusters(track);
1040 MarkSeedFree(track);
1041 return kClusterNotMatching;
1043 // cl->IncreaseClusterUsage(); // do this only for winners
1045 track->SetFake(!goodCl);
1047 AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
1048 goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
1050 AddSeedBranch(track);
1052 return kClusterMatching;
1055 //_________________________________________________________________________
1056 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
1058 // check if the seed should not be discarded
1059 const UShort_t kMask = 0xffff;
1060 if (flag==kMissingCluster) {
1061 int lastChecked = seed->GetLayerID();
1062 UShort_t patt = seed->GetHitsPattern();
1063 if (lastChecked) patt |= (~(kMask<<lastChecked))&fCurrTrackCond->GetAllowedLayers(); // not all layers were checked, complete unchecked ones by potential hits
1064 Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
1070 //______________________________________________________________________________
1071 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1073 // propagate seed to given x applying material correction if requested
1074 const Double_t kEpsilon = 1e-5;
1075 Double_t xpos = seed->GetX();
1076 Int_t dir = (xpos<xToGo) ? 1:-1;
1077 Double_t xyz0[3],xyz1[3],param[7];
1079 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1080 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1081 while ( (xToGo-xpos)*dir > kEpsilon){
1082 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1083 Double_t x = xpos+step;
1084 Double_t bz=GetBz(); // getting the local Bz
1085 if (!seed->PropagateToX(x,bz)) return kFALSE;
1087 if (matCorr || updTime) {
1088 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1091 seed->GetXYZ(xyz1); // // global pos at the end of step
1093 MeanMaterialBudget(xyz0,xyz1,param);
1094 Double_t xrho=param[0]*param[4], xx0=param[1];
1095 if (dir>0) xrho = -xrho; // outward should be negative
1096 if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
1099 else { // matCorr is not requested but time integral is
1100 double d0 = xyz1[0]-xyz0[0];
1101 double d1 = xyz1[1]-xyz0[1];
1102 double d2 = xyz1[2]-xyz0[2];
1103 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1106 if (updTime) seed->AddTimeStep(ds);
1107 xpos = seed->GetX();
1112 //______________________________________________________________________________
1113 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1115 // propagate seed to given x applying material correction if requested
1116 const Double_t kEpsilon = 1e-5;
1117 Double_t xpos = seed->GetX();
1118 Int_t dir = (xpos<xToGo) ? 1:-1;
1119 Double_t xyz0[3],xyz1[3],param[7];
1121 if (AliDebugLevelClass()>1) {
1122 AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
1123 seed->AliExternalTrackParam::Print();
1125 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1126 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1127 while ( (xToGo-xpos)*dir > kEpsilon){
1128 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1129 Double_t x = xpos+step;
1130 Double_t bz=GetBz(); // getting the local Bz
1131 if (!seed->PropagateTo(x,bz)) return kFALSE;
1133 if (matCorr || updTime) {
1134 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1137 seed->GetXYZ(xyz1); // // global pos at the end of step
1140 MeanMaterialBudget(xyz0,xyz1,param);
1141 Double_t xrho=param[0]*param[4], xx0=param[1];
1142 if (dir>0) xrho = -xrho; // outward should be negative
1143 if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1146 else { // matCorr is not requested but time integral is
1147 double d0 = xyz1[0]-xyz0[0];
1148 double d1 = xyz1[1]-xyz0[1];
1149 double d2 = xyz1[2]-xyz0[2];
1150 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1153 if (updTime) seed->AddTimeStep(ds);
1155 xpos = seed->GetX();
1157 if (AliDebugLevelClass()>1) {
1158 AliDebug(2,Form("After Propagate to X=%f",xToGo));
1159 seed->AliExternalTrackParam::Print();
1164 //______________________________________________________________________________
1165 Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
1167 // finalize single TPC track hypothesis
1168 if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
1169 AliITSUSeed* winner = 0;
1171 fCurrMass = hyp->GetMass();
1172 if (!(winner=hyp->DefineWinner(fCurrTrackCond->GetActiveLrInner()))) return kFALSE;
1175 #ifdef _ITSU_TUNING_MODE_ // fill tuning histos
1176 TObjArray* dest = hyp->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1178 AliITSUSeed* sd = winner;
1179 double htrPt = hyp->Pt();
1182 if ( (clID=sd->GetLrCluster(lrID))<0 ) continue;
1183 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Nrm ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2GloNrm());
1184 ((TH2*)dest->At( GetHistoID(lrID,kHBestInBranch,fCurrPassID) ))->Fill(htrPt,sd->GetOrdBranch());
1185 ((TH2*)dest->At( GetHistoID(lrID,kHBestInCand ,fCurrPassID) ))->Fill(htrPt,sd->GetOrdCand());
1187 if (dest==fCHistoArrFake && !sd->IsFake()) continue; // for the fake seeds fill only fake clusters part
1189 ((TH2*)dest->At( GetHistoID(lrID,kHResY ,fCurrPassID) ))->Fill(htrPt,sd->GetResidY());
1190 ((TH2*)dest->At( GetHistoID(lrID,kHResZ ,fCurrPassID) ))->Fill(htrPt,sd->GetResidZ());
1191 ((TH2*)dest->At( GetHistoID(lrID,kHResYP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullY());
1192 ((TH2*)dest->At( GetHistoID(lrID,kHResZP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullZ());
1193 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Cl ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2Cl());
1195 } while((sd=(AliITSUSeed*)sd->GetParent()));
1197 ((TH2*)dest->At( GetHistoID(-1,kHChiMatch,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSTPC());
1198 ((TH2*)dest->At( GetHistoID(-1,kHChiITSSA,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSSA());
1203 CheckClusterSharingConflicts(hyp);
1204 return hyp->GetWinner(); // winner might change of disappear after resolving conflicts
1208 //______________________________________________________________________________
1209 void AliITSUTrackerGlo::FinalizeHypotheses()
1211 // select winner for each hypothesis, remove cl. sharing conflicts
1212 AliCodeTimerAuto("",0);
1214 int* esdTrackIndex = fESDIndex.GetArray();
1215 for (int ih=0;ih<fNTracksESD;ih++) FinalizeHypothesis(GetTrackHyp(esdTrackIndex[ih])); // finlize in decreasing pt order
1217 for (int ih=fNTracksESD;ih--;) UpdateESDTrack(GetTrackHyp(esdTrackIndex[ih]),AliESDtrack::kITSin);
1221 //______________________________________________________________________________
1222 void AliITSUTrackerGlo::CheckClusterSharingConflicts(AliITSUTrackHyp* hyp)
1224 // remove eventual cluster sharing conflicts
1225 AliITSUSeed* winner = 0;
1226 if (!(winner=hyp->GetWinner())) return;
1227 UShort_t idH = (UShort_t)hyp->GetUniqueID()+1;
1229 int inLr = fCurrTrackCond->GetActiveLrInner();
1230 AliITSUSeed* winSD=winner;
1232 if ( (clID=winSD->GetLrCluster(lrID))<0 ) continue;
1233 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1234 Int_t refID = cl->GetRecoInfo(); // was it already referred by some track (with refID-1 if >0)
1235 if (!refID) {cl->SetRecoInfo(idH); continue;} // if not, refer from track to cluster
1236 if (refID==idH) continue; // ignore reference to itself
1237 // the cluster is already used by other track, need to resolve the conflict
1238 AliITSUTrackHyp* hypC = GetTrackHyp(refID-1); // competitor
1239 // printf("Ref to %d (%p) from %d Cl %d %d\n",refID-1,hypC, idH-1,clID,lrID);
1240 AliITSUSeed* winnerC = hypC->GetWinner();
1242 printf("Missing winner, ERROR\n");
1246 if (AliDebugLevelClass()>1) { //------------------------ DEBUG -------------------
1247 AliInfo(Form("Shared cl#%4d on Lr:%d: Hyp%3d/Q:%6.3f vs Hyp%3d/Q:%6.3f",
1248 clID,lrID,idH-1,winner->GetQualityVar(),refID-1,winnerC->GetQualityVar()));
1249 // dump winner of hyp
1250 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1251 idH-1,winner->Pt(),winner->GetChi2GloNrm(),winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),
1252 winner->GetNLayersHit(),hyp->GetITSLabel(),hyp->GetESDTrack()->GetTPCLabel());
1254 AliITSUSeed* sd = winner;
1257 int clIDt = sd->GetLrCluster(lrs);
1258 if (clIDt<0) continue;
1259 while( lrs>++prevL ) printf("%4s ","----");
1260 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1261 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1263 // dump winner of hypC
1264 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1265 refID-1,winnerC->Pt(),winnerC->GetChi2GloNrm(),winnerC->GetChi2ITSTPC(),winnerC->GetChi2ITSSA(),
1266 winnerC->GetNLayersHit(),hypC->GetITSLabel(),hypC->GetESDTrack()->GetTPCLabel());
1271 int clIDt = sd->GetLrCluster(lrs);
1272 if (clIDt<0) continue;
1273 while( lrs>++prevL ) printf("%4s ","----");
1274 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1275 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1277 } //-------------------------------------------------- DEBUG -------------------
1279 if (winner->GetQualityVar()<winnerC->GetQualityVar()) { // current tracks is better then competitor track
1280 FlagSeedClusters(winnerC,kFALSE,refID); // unflag cluster usage by loser
1281 cl->SetRecoInfo(idH);
1283 //if ( hypC->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",refID-1));
1285 if (hypC->DefineWinner(inLr)) { // find new winner instead of suppressed candidate
1287 CheckClusterSharingConflicts(hypC); // and check its sharing conflicts again
1288 if (winner->IsKilled()) break; // the current winner might have been killed during check of new winner of hypC!
1291 else { // competitor hypC is better than the hyp
1292 FlagSeedClusters(winner,kFALSE,idH); // unflag cluster usage by loser
1294 //if ( hyp->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",idH-1));
1295 if (hyp->DefineWinner(inLr)) {
1297 CheckClusterSharingConflicts(hyp);
1302 } while ((winSD=(AliITSUSeed*)winSD->GetParent()));
1307 //______________________________________________________________________________
1308 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1310 // update ESD track with current best hypothesis
1312 AliESDtrack* esdTr = hyp->GetESDTrack();
1313 if (!esdTr || esdTr->IsOn(flag)) return;
1314 AliITSUSeed* win = hyp->GetWinner();
1315 if (!win || win->IsKilled()) return;
1319 case AliESDtrack::kITSin:
1320 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1322 // set fakes cluster info
1324 UShort_t clfake = 0;
1325 AliITSUSeed* sd = win;
1327 if (sd->IsFake()) clfake |= 0x1<<sd->GetLayerID();
1328 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1330 // RS: Temporary set special flag for tracks from the afterburner
1331 if (fCurrPassID>0) clfake |= 0x1<<7;
1333 esdTr->SetITSSharedMap(clfake);
1335 // TEMPORARY: store iteration id
1336 esdTr->SetITSModuleIndex(11,fCurrPassID);
1339 case AliESDtrack::kITSout:
1340 // here the stored chi2 will correspond to backward ITS-SA fit
1341 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1343 // TODO: avoid setting friend
1346 case AliESDtrack::kITSrefit:
1348 // at the moment fill as a chi2 the TPC-ITS matching chi2
1349 chiSav = hyp->GetChi2();
1350 hyp->SetChi2(win->GetChi2ITSTPC());
1351 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1352 hyp->SetChi2(chiSav);
1354 // TODO: avoid setting cluster info
1357 AliFatal(Form("Unknown flag %d",flag));
1360 esdTr->SetITSLabel(hyp->GetITSLabel());
1361 // transfer module indices
1365 //______________________________________________________________________________
1366 Double_t AliITSUTrackerGlo::RefitTrack(AliExternalTrackParam* trc, Double_t rDest, Int_t &nclFit, Int_t stopCond)
1368 // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
1369 // if stopCond<0 : propagate till last cluster then stop
1370 // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
1371 // if stopCond>0 : rDest must be reached
1373 double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1374 int dir,lrStart,lrStop;
1376 dir = rCurr<rDest ? 1 : -1;
1377 lrStart = fITS->FindFirstLayerID(rCurr,dir);
1378 lrStop = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
1380 if (AliDebugLevelClass()>2) {
1381 printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1382 printf("Before refit: "); trc->Print();
1384 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));
1386 Int_t clInfo[2*AliITSUAux::kMaxLayers];
1387 Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
1388 fCurrMass = fCurrHyp->GetMass();
1389 AliExternalTrackParam tmpTr(*trc);
1394 int lrStop1 = lrStop+1;
1395 for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
1396 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1397 if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1398 int ilrA2,ilrA = lr->GetActiveID();
1399 // passive layer or active w/o hits will be traversed on the way to next cluster
1400 if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue;
1402 // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1404 if (dir>0) { // clusters are stored in increasing radius order
1405 iclLr[nclLr++]=clInfo[ilrA2++];
1406 if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
1409 if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
1410 iclLr[nclLr++]=clInfo[ilrA2];
1413 Bool_t transportedToLayer = kFALSE;
1414 for (int icl=0;icl<nclLr;icl++) {
1415 AliITSUClusterPix* clus = (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1416 AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1417 if (!tmpTr.Rotate(sens->GetPhiTF())) {
1418 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1419 AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
1423 double xClus = sens->GetXTF()+clus->GetX();
1424 if (!transportedToLayer) {
1425 if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1426 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1427 AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
1428 return -1; // go to the entrance to the layer
1431 transportedToLayer = kTRUE;
1434 if (AliDebugLevelClass()>1) {
1435 AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1438 if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1439 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1440 AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
1444 Double_t p[2]={clus->GetY(), clus->GetZ()};
1445 Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
1446 double chi2cl = tmpTr.GetPredictedChi2(p,cov);
1449 #ifdef _ITSU_TUNING_MODE_
1450 TObjArray* dest = trc->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1451 if (dest && fTrackPhaseID>kClus2Tracks) {
1453 double htrPt = tmpTr.Pt();
1454 double dy = p[0]-tmpTr.GetY();
1455 double dz = p[1]-tmpTr.GetZ();
1456 ((TH2*)dest->At( GetHistoID(ilrA,kHResY,0,fTrackPhaseID) ))->Fill(htrPt,dy);
1457 ((TH2*)dest->At( GetHistoID(ilrA,kHResZ,0,fTrackPhaseID) ))->Fill(htrPt,dz);
1458 double errY = tmpTr.GetSigmaY2();
1459 double errZ = tmpTr.GetSigmaZ2();
1460 if (errY>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResYP,0,fTrackPhaseID) ))->Fill(htrPt,dy/Sqrt(errY));
1461 if (errZ>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResZP,0,fTrackPhaseID) ))->Fill(htrPt,dz/Sqrt(errZ));
1462 ((TH2*)dest->At( GetHistoID(ilrA,kHChi2Cl,0,fTrackPhaseID) ))->Fill(htrPt,chi2cl);
1466 if ( !tmpTr.Update(p,cov) ) {
1467 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1468 AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
1471 if (AliDebugLevelClass()>1) {
1472 printf("AfterRefit: "); tmpTr.Print();
1474 if (++nclFit==nCl && stopCond<0) {
1476 return chi2; // it was requested to not propagate after last update
1481 // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1482 // Still, try to go as close as possible to rDest.
1484 if (lrStart!=lrStop) {
1485 if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1486 AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1487 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1489 if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1490 AliDebug(1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1491 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1494 // go to the destination radius. Note that here we don't select direction to avoid precision problems
1495 if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
1496 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1499 if (AliDebugLevelClass()>2) {
1500 printf("After refit (now at lr %d): ",lrStop); trc->Print();
1505 //__________________________________________________________________
1506 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp)
1510 const int kMaxLbPerCl = 3;
1511 int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1512 Int_t lr,nLab=0,nCl=0;
1513 AliITSUSeed *seed = hyp->GetWinner();
1515 int clID = seed->GetLrCluster(lr);
1517 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1519 for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1520 int trLb = cl->GetLabel(imc);
1522 // search this mc track in already accounted ones
1524 for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1525 if (iLab<nLab) lbStat[iLab]++;
1530 } // loop over given cluster's labels
1532 seed = (AliITSUSeed*)seed->GetParent();
1533 } // loop over clusters
1535 AliESDtrack* esdTr = hyp->GetESDTrack();
1536 int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1538 int maxLab=0,nTPCok=0;
1539 for (int ilb=nLab;ilb--;) {
1540 int st = lbStat[ilb];
1541 if (lbStat[maxLab]<st) maxLab=ilb;
1542 if (tpcLab==lbID[ilb]) nTPCok += st;
1544 hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1545 hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1546 hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
1550 hyp->SetFakeRatio(-1.);
1551 hyp->SetLabel( -tpcLab );
1552 hyp->SetITSLabel( kDummyLabel );
1555 //__________________________________________________________________
1556 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1558 // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1559 if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
1560 AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
1561 fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
1562 AliITSUSeed** tmpArr = fLayerCandidates;
1563 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
1564 memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
1565 delete tmpArr; // delete only array, not objects
1567 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1568 int slot=fNBranchesAdded++;
1569 for (int slotF=slot;slotF--;) { // slotF is always slot-1
1570 AliITSUSeed* si = branches[slotF];
1571 if (si->Compare(seed)<0) break; // found the last seed with better quality
1572 // otherwise, shift the worse seed to the next slot
1573 branches[slot] = si;
1574 slot = slotF; // slot should be slotF+1
1576 // if needed, move worse seeds
1577 branches[slot] = seed;
1582 //__________________________________________________________________
1583 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1585 // keep allowed number of branches for current seed and disable extras
1586 AliCodeTimerAuto("",0);
1587 int nb = Min(fNBranchesAdded,acceptMax);
1588 // if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1589 // disable unused branches
1590 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1592 #ifdef _ITSU_TUNING_MODE_
1593 for (int ib=0;ib<fNBranchesAdded;ib++) branches[ib]->SetOrdBranch(ib);
1596 for (int ib=nb;ib<fNBranchesAdded;ib++) MarkSeedFree(branches[ib]);
1597 fNCandidatesAdded += nb; // update total candidates counter
1598 fNBranchesAdded = 0; // reset branches counter
1602 //__________________________________________________________________
1603 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1605 // transfer allowed number of branches to hypothesis container
1607 AliCodeTimerAuto("",0);
1608 // sort candidates in increasing order of chi2
1609 static int lastSize = 0;
1610 static int *index = 0;
1611 static float *chi2 = 0;
1612 if (fLayerMaxCandidates>lastSize) {
1613 lastSize = fLayerMaxCandidates;
1616 index = new int[lastSize];
1617 chi2 = new float[lastSize];
1619 for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1620 Sort(fNCandidatesAdded,chi2,index,kFALSE);
1623 if (ilr>fCurrTrackCond->GetActiveLrInner()) { // just take 1st acceptMax candidates
1624 nb = Min(fNCandidatesAdded,acceptMax);
1625 for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1627 else { // for innermost layer test ITS SA backward chi2 and matching to TPC
1628 AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
1629 for (nb=0;nb<fNCandidatesAdded;nb++) {
1630 AliITSUSeed* sd = fLayerCandidates[index[nb]];
1631 fCurrHyp->SetWinner(sd);
1632 if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
1634 AddProlongationHypothesis(sd,ilr);
1635 if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
1638 fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
1640 // disable unused candidates
1641 for (int ib=nb;ib<fNCandidatesAdded;ib++) MarkSeedFree(fLayerCandidates[index[ib]]);
1642 fNCandidatesAdded = 0; // reset candidates counter
1646 //__________________________________________________________________
1647 Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
1649 // check seed backward propagation chi2 and matching to TPC
1650 double bz0 = GetBz();
1651 double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
1652 static AliExternalTrackParam trback;
1654 trback.ResetCovariance(10000);
1656 double chi2sa = RefitTrack(&trback,rDest,nclFit,1);
1657 if (chi2sa<0) return kFALSE;
1658 int ndf = nclFit*2-5;
1659 if (ndf>0) chi2sa /= ndf;
1660 if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2(nclFit)) return kFALSE;
1662 // relate to TPC track at outer layer
1663 AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
1664 if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
1665 if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
1666 AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
1672 double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
1673 if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
1675 seed->SetChi2ITSSA(chi2sa);
1676 seed->SetChi2ITSTPC(chi2Match);
1680 //__________________________________________________________________
1681 void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
1683 // remove those hypothesis seeds which dont lead to candidates at final layer
1685 // 1st, flag the seeds to save
1686 int lr0 = fCurrTrackCond->GetActiveLrInner();
1687 for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
1688 AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
1689 if (!seed) continue;
1690 seed->FlagTree(AliITSUSeed::kSave);
1691 dest->AddSeed(seed,lr0);
1693 for (int ilr=0;ilr<fNLrActive;ilr++) {
1694 int nsd = fCurrHyp->GetNSeeds(ilr);
1695 for (int isd=0;isd<nsd;isd++) {
1696 AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
1697 if (!seed) continue; // already discarded or saved
1698 if (seed->IsSaved()) dest->AddSeed(seed,ilr);
1699 else MarkSeedFree(seed);
1703 // AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
1706 //__________________________________________________________________
1707 void AliITSUTrackerGlo::CleanHypothesis(AliITSUTrackHyp* hyp)
1711 for (int ilr=0;ilr<fNLrActive;ilr++) {
1712 int nsd = hyp->GetNSeeds(ilr);
1713 for (int isd=0;isd<nsd;isd++) {
1714 AliITSUSeed* seed = hyp->RemoveSeed(ilr,isd);
1715 if (!seed) continue; // already discarded or saved
1721 //__________________________________________________________________
1722 void AliITSUTrackerGlo::FlagSplitClusters()
1724 // set special bit on split clusters using MC info
1725 for (int ilr=fNLrActive;ilr--;) {
1727 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1728 for (int isn=lr->GetNSensors();isn--;) {
1729 AliITSURecoSens* sens = lr->GetSensor(isn);
1730 int nCl = sens->GetNClusters();
1732 int cl0 = sens->GetFirstClusterId();
1733 for (int icl=nCl;icl--;) {
1734 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1735 for (int icl1=icl;icl1--;) {
1736 AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1737 if (cl->HasCommonTrack(cl1)) {
1738 if (!cl->IsSplit()) nsplit++;
1739 if (!cl1->IsSplit()) nsplit++;
1746 AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1751 //__________________________________________________________________
1752 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1754 // check if the seed contains split cluster with size < maxSize
1756 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1757 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1758 if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1760 return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1764 //__________________________________________________________________
1765 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1767 // print seeds clusters
1769 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1770 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1773 if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1777 //__________________________________________________________________
1778 void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg, UShort_t ref)
1780 // mark used clusters
1783 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1784 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1785 if (ref) { // do we need to set or delete cluster->track ref?
1787 if (!cl->GetRecoInfo()) cl->SetRecoInfo(ref); // set ref only if cluster already does not refer to other track, inc.counter
1790 if (cl->GetRecoInfo()==ref) cl->SetRecoInfo(0); // unset reference only if it refers to ref, decrease counter
1794 seed = (AliITSUSeed*)seed->GetParent();
1799 //__________________________________________________________________
1800 void AliITSUTrackerGlo::CheckClusterUsage()
1802 // check cluster usage info
1803 for (int ilr=0;ilr<fNLrActive;ilr++) {
1804 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1805 int ncl = lr->GetNClusters();
1807 int nclUs0CntShr1 = 0;
1811 for (int icl=0;icl<ncl;icl++) {
1812 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(icl);
1813 int usc = cl->GetClUsage();
1815 if (cl->IsClusterUsed() || cl->IsClusterShared()) {
1817 printf("Lr%d Cl%4d: UseCnt=0 but UseFlg=%d UseShr=%d\n",ilr,icl,cl->IsClusterUsed(),cl->IsClusterShared());
1822 if (!cl->IsClusterUsed()) {
1824 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1826 if (cl->IsClusterShared()) {
1828 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1832 printf("Lr%d Cl%4d: UseCnt=%d!! but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1835 printf("ClUsage Lr%d: Ncl=%5d Nusd=%5d (%5.2f) Use0CS1:%4d Use1C0:%4d Use1S1:%4d UseL:%d\n",
1836 ilr,ncl,nclUsed,ncl? float(nclUsed)/ncl : 0, nclUs0CntShr1,nclUs1Cnt0,nclUs1Shr1,nclUseL);
1841 #ifdef _ITSU_TUNING_MODE_
1842 //__________________________________________________________________
1843 TObjArray* AliITSUTrackerGlo::BookControlHistos(const char* pref)
1845 // book special control histos
1846 TString prefS = pref;
1848 TObjArray* dest = new TObjArray;
1849 dest->SetOwner(kTRUE);
1851 const int kNResDef=7;
1852 const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1853 const double ptMax=10;
1854 const double plMax=10;
1855 const double chiMax=100;
1856 const int nptbins=50;
1857 const int nresbins=400;
1858 const int nplbins=50;
1859 const int nchbins=200;
1860 const int maxBr = 15;
1861 const int maxCand = 200;
1863 for (int phase=0;phase<kNTrackingPhases;phase++) {
1864 for (int pass=0;pass<AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();pass++) {
1866 for (int ilr=0;ilr<fNLrActive;ilr++) {
1868 // ----------------- These are histos to be filled in Cluster2Tracks of each pass.
1869 // PropagateBack and RefitInward will be stored among the histos of 1st pass
1870 if (pass>0 && phase!=kClus2Tracks) continue;
1872 double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1873 ttl = Form("Pass%d_S%d_residY%d_%s",pass,phase,ilr,pref);
1874 TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1875 dest->AddAtAndExpand(hdy,GetHistoID(ilr,kHResY,pass,phase));
1876 hdy->SetDirectory(0);
1878 ttl = Form("Pass%d_S%d_residYPull%d_%s",pass,phase,ilr,pref);
1879 TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1880 dest->AddAtAndExpand(hdyp,GetHistoID(ilr,kHResYP,pass,phase));
1881 hdyp->SetDirectory(0);
1883 ttl = Form("Pass%d_S%d_residZ%d_%s",pass,phase,ilr,pref);
1884 TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1885 dest->AddAtAndExpand(hdz,GetHistoID(ilr,kHResZ,pass,phase));
1886 hdz->SetDirectory(0);
1888 ttl = Form("Pass%d_S%d_residZPull%d_%s",pass,phase,ilr,pref);
1889 TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1890 hdzp->SetDirectory(0);
1891 dest->AddAtAndExpand(hdzp,GetHistoID(ilr,kHResZP,pass,phase));
1893 ttl = Form("Pass%d_S%d_chi2Cl%d_%s",pass,phase,ilr,pref);
1894 TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1895 hchi->SetDirectory(0);
1896 dest->AddAtAndExpand(hchi,GetHistoID(ilr,kHChi2Cl,pass,phase));
1898 // ------------------- These histos are filled for Clusters2Tracks only
1899 if (phase!=kClus2Tracks) continue;
1901 ttl = Form("Pass%d_S%d_chi2Nrm%d_%s",pass,phase,ilr,pref);
1902 TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1903 hchiN->SetDirectory(0);
1904 dest->AddAtAndExpand(hchiN,GetHistoID(ilr,kHChi2Nrm,pass,phase));
1906 ttl = Form("Pass%d_S%d_bestInBranch%d_%s",pass,phase,ilr,pref);
1907 TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
1908 hnbr->SetDirectory(0);
1909 dest->AddAtAndExpand(hnbr,GetHistoID(ilr,kHBestInBranch,pass,phase));
1911 ttl = Form("Pass%d_S%d_bestInCands%d_%s",pass,phase,ilr,pref);
1912 TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
1913 hncn->SetDirectory(0);
1914 dest->AddAtAndExpand(hncn,GetHistoID(ilr,kHBestInCand,pass,phase));
1916 } // loop over layers
1918 // these are special histos, filled not per layer but in the end of track fit in Clusters2Tracks in EVERY pass
1920 if (phase==kClus2Tracks) {
1922 ttl = Form("Pass%d_Chi2Match_%s",pass,pref);
1923 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1924 hchiMatch->SetDirectory(0);
1925 dest->AddAtAndExpand(hchiMatch,GetHistoID(-1,kHChiMatch,pass,phase));
1928 ttl = Form("Pass%d_Chi2ITSSA_%s",pass,pref);
1929 hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
1930 hchiSA->SetDirectory(0);
1931 dest->AddAtAndExpand(hchiSA,GetHistoID(-1,kHChiITSSA,pass,phase));
1934 } // loop over tracking passes
1935 }// loop over tracking phases
1940 //__________________________________________________________________
1941 Int_t AliITSUTrackerGlo::GetHistoID(Int_t lr, Int_t hid, Int_t pass, Int_t phase)
1943 // get id for the requested histo
1946 return pass*kHistosPass + phase*kHistosPhase + lr*kMaxHID + hid;