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);
252 AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMass(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
253 FindTrack(fCurrESDtrack, trID);
256 if (AliDebugLevelClass()>+2) {
257 AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
260 FinalizeHypotheses();
261 #ifdef _ITSU_TUNING_MODE_
262 CheckClusterUsage(); //!!RS
265 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0); // in case wtc array was defined, uncomment this
268 AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,fNTracksESD));
273 //_________________________________________________________________________
274 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
277 // Do outward fits in ITS
278 AliCodeTimerAuto("",0);
280 SetTrackingPhase(kPropBack);
281 fNTracksESD = esdEv->GetNumberOfTracks();
282 AliDebug(1,Form("Will propagate back %d tracks",fNTracksESD));
284 double bz0 = GetBz();
285 Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
286 AliITSUTrackHyp dummyTr;
287 const double kWatchStep=10.; // for larger steps watch arc vs segment difference
288 Double_t times[AliPID::kSPECIES];
290 for (int itr=0;itr<fNTracksESD;itr++) {
291 fCurrESDtrack = esdEv->GetTrack(itr);
292 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
293 // Start time integral and add distance from current position to vertex
294 if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
296 fCurrESDtrack->GetXYZ(xyzTrk);
297 Double_t dst = 0.; // set initial track lenght, tof
299 double dxs = xyzTrk[0] - xyzVtx[0];
300 double dys = xyzTrk[1] - xyzVtx[1];
301 double dzs = xyzTrk[2] - xyzVtx[2];
302 // RS: for large segment steps d use approximation of cicrular arc s by
303 // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
304 // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
305 // Hence s^2/d^2 = (1+1/6 p^2)^2
306 dst = dxs*dxs + dys*dys;
307 if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
308 double crv = Abs(fCurrESDtrack->GetC(bz0));
309 double fcarc = 1.+crv*crv*dst/6.;
316 fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
318 if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
319 dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
320 dummyTr.StartTimeIntegral();
321 dummyTr.AddTimeStep(dst);
322 dummyTr.GetIntegratedTimes(times);
323 fCurrESDtrack->SetIntegratedTimes(times);
324 fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
328 fCurrHyp = GetTrackHyp(itr);
329 fCurrMass = fCurrHyp->GetMass();
330 fCurrHyp->StartTimeIntegral();
331 fCurrHyp->AddTimeStep(dst);
332 fCurrHyp->ResetCovariance(10000);
333 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax());
334 if (chi2>0) { // propagate to exit from the ITS/TPC screen
335 int ndf = fCurrHyp->GetWinner()->GetNLayersHit()*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();
380 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin());
381 if (chi2>0) { // propagate up to inside radius of the beam pipe
382 fCurrHyp->SetChi2(chi2);
383 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
386 AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
390 AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
391 fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
393 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
397 //_________________________________________________________________________
398 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
400 // read from tree (if pointer provided) or directly from the ITS reco interface
402 return fReconstructor->LoadClusters(treeRP);
405 //_________________________________________________________________________
406 void AliITSUTrackerGlo::UnloadClusters()
412 Info("UnloadClusters","To be implemented");
414 //_________________________________________________________________________
415 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
421 Info("GetCluster","To be implemented");
425 //_________________________________________________________________________
426 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr, Int_t esdID)
428 // do we need to match this track to ITS?
430 static double bz = GetBz();
432 // is it already prolonged
433 if (esdTr->IsOn(AliESDtrack::kITSin)) return kFALSE;
434 AliITSUTrackHyp* hyp = GetTrackHyp(esdID); // in case there is unfinished hypothesis from prev.pass, clean it
436 CleanHypothesis(hyp);
437 if (hyp->GetSkip()) return kFALSE; // need to skip this
441 if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
442 esdTr->IsOn(AliESDtrack::kTPCout) ||
443 esdTr->IsOn(AliESDtrack::kITSin) ||
444 esdTr->GetKinkIndex(0)>0) return kFALSE;
446 if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
449 esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz);
450 // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
451 if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation())
452 && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
453 Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
455 // if (esdTr->Pt()<3) return kFALSE;//RS
459 //_________________________________________________________________________
460 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
462 // find prolongaion candidates finding for single seed
464 AliITSUSeed seedUC; // copy of the seed from the upper layer
465 AliITSUSeed seedT; // transient seed between the seedUC and new prolongation hypothesis
467 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
468 AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID); // initialize prolongations hypotheses tree
469 if (!hypTr || hypTr->GetSkip()) return;
473 fCurrHyp->InitFrom(hypTr);
475 AliITSURecoSens *hitSens[AliITSURecoSens::kNNeighbors+1];
477 int ilaUp = fNLrActive; // previous active layer really checked (some may be excluded!)
479 for (int ila=fNLrActive;ila--;ilaUp=fCurrActLrID) {
480 if (fCurrTrackCond->IsLayerExcluded(ila)) continue;
482 fCurrLayer = fITS->GetLayerActive(ila);
483 Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
485 // for the outermost layer the seed is created from the ESD track
486 int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
487 int maxNBranches = fCurrTrackCond->GetMaxBranches(ila);
488 int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
490 for (int isd=0;isd<nSeedsUp;isd++) {
492 if (ilaUp==fNLrActive) {
494 seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
497 seedU = fCurrHyp->GetSeed(ilaUp,isd); // seed on prev.active layer to prolong
498 if (seedU->IsKilled()) continue;
499 seedUC = *seedU; // its copy will be prolonged
500 seedUC.SetParent(seedU);
502 seedUC.ResetFMatrix(); // reset the matrix for propagation to next layer
503 // go till next active layer
504 AliDebug(2,Form("working on Lr:%d Seed:%d of %d for esdID=%d (MClb:%d) | pT=%.3f",ila,isd,nSeedsUp,esdID,fCurrESDtrMClb,seedUC.Pt()));
505 if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
507 AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
508 // Check if the seed satisfies to track definition
509 if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
510 continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
512 if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
513 if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
519 if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
520 mcquest = fCurrESDtrMClb;
522 printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
525 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest); // find detectors which may be hit by the track
527 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens); // find detectors which may be hit by the track
528 AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
532 for (int isn=nsens;isn--;) {
533 AliITSURecoSens* sens = hitSens[isn];
534 int ncl = sens->GetNClusters();
538 // We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
539 // since the transport matrix should be defined in this frame.
540 double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
541 if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
542 if (AliDebugLevelClass()>2) {
543 printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
548 if (xs<seedT.GetX()) {
549 if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
551 else { // some low precision tracks may hit the sensor plane outside of the layer radius
552 if (AliDebugLevelClass()>2) {
553 if (!seedT.ContainsFake()) {
554 printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
559 if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
561 // if (!seedT.PropagateToX(xs,bz)) continue;
562 // if (!seedT.Rotate(sens->GetPhiTF())) continue;
563 if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
565 int clID0 = sens->GetFirstClusterId();
566 for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
567 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
568 int clID = clID0 + icl;
569 if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
570 int res = CheckCluster(&seedT,ila,clID0+icl);
571 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
573 if (res==kStopSearchOnSensor) break; // stop looking on this sensor
574 if (res==kClusterNotMatching) continue; // cluster does not match
575 // cluster is matching and it was added to the hypotheses tree
578 // cluster search is done. Do we need to have a version of this seed skipping current layer
579 if (!NeedToKill(&seedUC,kMissingCluster)) {
580 AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
581 double penalty = -fCurrTrackCond->GetMissPenalty(ila);
582 // to do: make penalty to account for probability to miss the cluster for good reason
583 seedSkp->SetChi2Cl(penalty);
584 if (seedSkp->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(ila)) {
585 MarkSeedFree(seedSkp);
587 else AddSeedBranch(seedSkp);
589 // transfer the new branches of the seed to the hypothesis container
590 if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
593 if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
594 // ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
595 if (AliDebugLevelClass()>2) { //RS
596 printf(">>> All hypotheses on lr %d: \n",ila);
597 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
601 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
602 if (ila!=0) continue;
603 double vecL[5] = {0};
604 double matL[15] = {0};
605 AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
606 while(sp->GetParent()) {
607 sp->Smooth(vecL,matL);
608 if (sp->GetLayerID()>=fNLrActive-1) break;
609 sp = (AliITSUSeed*)sp->GetParent();
614 SaveReducedHypothesesTree(hypTr);
615 if (AliDebugLevelClass()>1) {
616 printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
620 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
624 //_________________________________________________________________________
625 AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
627 // init prolongaion candidates finding for single seed
628 AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
631 fCountProlongationTrials++;
633 fCurrMass = esdTr->GetMass();
634 if (fCurrMass<kPionMass*0.9) fCurrMass = kPionMass; // don't trust to mu, e identification from TPCin
636 hyp = new AliITSUTrackHyp(fNLrActive);
637 hyp->SetESDTrack(esdTr);
638 hyp->SetUniqueID(esdID);
639 hyp->SetMass(fCurrMass);
640 hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
641 SetTrackHyp(hyp,esdID);
643 if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
649 //_________________________________________________________________________
650 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
652 // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
655 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
657 int dir = lTo > lFrom ? 1:-1;
658 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
659 Bool_t checkFirst = kTRUE;
660 Bool_t limReached = kFALSE;
663 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
666 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
667 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
669 // go the entrance of the layer, assuming no materials in between
670 double xToGo = lrTo->GetR(-dir);
673 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
676 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
679 // double xts = xToGo;
680 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
681 //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
682 //seed->Print("etp");
686 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
687 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
688 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
689 //seed->Print("etp");
693 if (limReached) break;
699 //_________________________________________________________________________
700 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
702 // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
704 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
706 int dir = lTo > lFrom ? 1:-1;
707 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
708 Bool_t checkFirst = kTRUE;
709 Bool_t limReached = kFALSE;
712 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
715 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
716 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
718 // go the entrance of the layer, assuming no materials in between
719 double xToGo = lrTo->GetR(-dir);
722 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
725 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
728 // double xts = xToGo;
729 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
730 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
731 // seed->Print("etp");
734 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
735 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
736 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
737 //seed->Print("etp");
741 if (limReached) break;
747 //_________________________________________________________________________
748 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
750 // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
752 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
754 int dir = lTo > lFrom ? 1:-1;
755 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
756 Bool_t checkFirst = kTRUE;
759 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
762 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
763 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
765 // go the entrance of the layer, assuming no materials in between
766 double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
768 // double xts = xToGo;
769 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
770 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
771 // seed->Print("etp");
774 if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
776 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
777 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
778 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
779 //seed->Print("etp");
788 //_________________________________________________________________________
789 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
791 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
792 // If check is requested, do this only provided the track has not exited the layer already
793 double xToGo = lr->GetR(dir);
794 if (check) { // do we need to track till the surface of the current layer ?
795 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
796 AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
797 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
798 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
800 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
801 // go via layer to its boundary, applying material correction.
802 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
808 //_________________________________________________________________________
809 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
811 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
812 // If check is requested, do this only provided the track has not exited the layer already
813 double xToGo = lr->GetR(dir);
814 if (check) { // do we need to track till the surface of the current layer ?
815 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
816 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
817 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
819 AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
820 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
821 // go via layer to its boundary, applying material correction.
822 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
829 //_________________________________________________________________________
830 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
832 // go to the entrance of lr in direction dir, w/o applying material corrections.
833 // If check is requested, do this only provided the track did not reach the layer already
834 double xToGo = lr->GetR(-dir);
835 if (check) { // do we need to track till the surface of the current layer ?
836 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
837 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
838 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
840 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
841 // go via layer to its boundary, applying material correction.
842 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
847 //_________________________________________________________________________
848 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
850 // go to the entrance of lr in direction dir, w/o applying material corrections.
851 // If check is requested, do this only provided the track did not reach the layer already
852 double xToGo = lr->GetR(-dir);
853 if (check) { // do we need to track till the surface of the current layer ?
854 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
855 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
856 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
858 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
859 // go via layer to its boundary, applying material correction.
860 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
865 //_________________________________________________________________________
866 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
868 // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
869 // as well as some aux info
871 AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
872 seed->GetXYZ(&fTrImpData[kTrXIn]); // lab position at the entrance from above
873 static AliExternalTrackParam sc; // seed copy for manipulations
876 fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
877 if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
878 double dr = lrA->GetDR(); // approximate X dist at the inner radius
879 if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
880 // special case: track does not reach inner radius, might be tangential
881 double r = sc.GetD(0,0,bz);
883 if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
884 dr = Abs(sc.GetX() - x);
885 if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
888 fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
889 double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
890 double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
891 sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
892 sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
893 double dphi0 = 0.5*Abs(fTrImpData[kTrPhiOut]-fTrImpData[kTrPhiIn]);
894 double phi0 = 0.5*(fTrImpData[kTrPhiOut]+fTrImpData[kTrPhiIn]);
895 if ( dphi0>(0.5*Pi()) ) {
896 // special case of angles around pi
901 fTrImpData[kTrPhi0] = phi0;
902 fTrImpData[kTrZ0] = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
903 dphi0 += sgy/lrA->GetR();
904 fTrImpData[kTrDPhi] = dphi0<PiOver2() ? dphi0 : PiOver2();
905 fTrImpData[kTrDZ] = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn]) + sgz;
910 //________________________________________
911 void AliITSUTrackerGlo::ResetSeedsPool()
913 // mark all seeds in the pool as unused
914 AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
916 fSeedsPool.Clear(); // seeds don't allocate memory
920 //________________________________________
921 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd)
923 // account that this seed is "deleted"
924 int id = sd->GetPoolID();
926 AliError(Form("Freeing of seed %p NOT from the pool is requested",sd));
929 // AliInfo(Form("%d %p",id, seed));
930 fSeedsPool.RemoveAt(id);
931 if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
932 fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
935 //_________________________________________________________________________
936 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID)
938 // Check if the cluster (in tracking frame!) is matching to track.
939 // The track must be already propagated to sensor tracking frame.
940 // Returns: kStopSearchOnSensor if the search on given sensor should be stopped,
941 // kClusterMatching if the cluster is matching
942 // kClusterMatching otherwise
944 const double kTolerX = 5e-4;
945 // AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
946 AliCluster *cl = fCurrLayer->GetCluster(clID);
948 Bool_t goodCl = kFALSE;
949 int currLabel = Abs(fCurrESDtrMClb);
951 if (cl->GetLabel(0)>=0) {
952 for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
955 if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
956 if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
958 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
959 AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
963 return kStopSearchOnSensor; // propagation failed, seedT is intact
966 double dy = cl->GetY()-track->GetY();
967 double dz = cl->GetZ()-track->GetZ();
970 double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
971 fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
972 if (dy2>tol2) { // the clusters are sorted in Z(col) then in Y(row).
974 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
975 AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
978 PrintSeedClusters(track);
980 // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
981 // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
982 if (dy<0) return kStopSearchOnSensor; // No chance that other cluster of this sensor will match (all Y's will be even larger)
983 else return kClusterNotMatching; // Other clusters may match
986 tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
987 fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
989 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
990 AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
993 PrintSeedClusters(track);
995 return kClusterNotMatching; // Other clusters may match
999 Double_t p[2]={cl->GetY(), cl->GetZ()};
1000 Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
1001 double chi2 = track->GetPredictedChi2(p,cov);
1003 if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
1004 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1005 AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1006 lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1007 track->Print("etp");
1009 PrintSeedClusters(track);
1011 return kClusterNotMatching;
1014 track = NewSeedFromPool(track); // input track will be reused, use its clone for updates
1015 if (!track->Update()) {
1016 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1017 AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
1018 track->Print("etp");
1020 PrintSeedClusters(track);
1022 MarkSeedFree(track);
1023 return kClusterNotMatching;
1025 track->SetChi2Cl(chi2);
1026 track->SetLrClusterID(lr,clID);
1028 if (track->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(lr)) {
1029 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1030 AliDebug(2,Form("Lost good cl on L:%d : Chi2Glo=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1031 lr,track->GetChi2GloNrm(),fCurrTrackCond->GetMaxChi2GloNrm(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1032 track->Print("etp");
1034 PrintSeedClusters(track);
1036 MarkSeedFree(track);
1037 return kClusterNotMatching;
1039 // cl->IncreaseClusterUsage(); // do this only for winners
1041 track->SetFake(!goodCl);
1043 AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
1044 goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
1046 AddSeedBranch(track);
1048 return kClusterMatching;
1051 //_________________________________________________________________________
1052 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
1054 // check if the seed should not be discarded
1055 const UShort_t kMask = 0xffff;
1056 if (flag==kMissingCluster) {
1057 int lastChecked = seed->GetLayerID();
1058 UShort_t patt = seed->GetHitsPattern();
1059 if (lastChecked) patt |= ~(kMask<<lastChecked); // not all layers were checked, complete unchecked ones by potential hits
1060 Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
1066 //______________________________________________________________________________
1067 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1069 // propagate seed to given x applying material correction if requested
1070 const Double_t kEpsilon = 1e-5;
1071 Double_t xpos = seed->GetX();
1072 Int_t dir = (xpos<xToGo) ? 1:-1;
1073 Double_t xyz0[3],xyz1[3],param[7];
1075 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1076 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1077 while ( (xToGo-xpos)*dir > kEpsilon){
1078 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1079 Double_t x = xpos+step;
1080 Double_t bz=GetBz(); // getting the local Bz
1081 if (!seed->PropagateToX(x,bz)) return kFALSE;
1083 if (matCorr || updTime) {
1084 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1087 seed->GetXYZ(xyz1); // // global pos at the end of step
1089 MeanMaterialBudget(xyz0,xyz1,param);
1090 Double_t xrho=param[0]*param[4], xx0=param[1];
1091 if (dir>0) xrho = -xrho; // outward should be negative
1092 if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
1095 else { // matCorr is not requested but time integral is
1096 double d0 = xyz1[0]-xyz0[0];
1097 double d1 = xyz1[1]-xyz0[1];
1098 double d2 = xyz1[2]-xyz0[2];
1099 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1102 if (updTime) seed->AddTimeStep(ds);
1103 xpos = seed->GetX();
1108 //______________________________________________________________________________
1109 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1111 // propagate seed to given x applying material correction if requested
1112 const Double_t kEpsilon = 1e-5;
1113 Double_t xpos = seed->GetX();
1114 Int_t dir = (xpos<xToGo) ? 1:-1;
1115 Double_t xyz0[3],xyz1[3],param[7];
1117 if (AliDebugLevelClass()>1) {
1118 AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
1119 seed->AliExternalTrackParam::Print();
1121 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1122 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1123 while ( (xToGo-xpos)*dir > kEpsilon){
1124 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1125 Double_t x = xpos+step;
1126 Double_t bz=GetBz(); // getting the local Bz
1127 if (!seed->PropagateTo(x,bz)) return kFALSE;
1129 if (matCorr || updTime) {
1130 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1133 seed->GetXYZ(xyz1); // // global pos at the end of step
1136 MeanMaterialBudget(xyz0,xyz1,param);
1137 Double_t xrho=param[0]*param[4], xx0=param[1];
1138 if (dir>0) xrho = -xrho; // outward should be negative
1139 if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1142 else { // matCorr is not requested but time integral is
1143 double d0 = xyz1[0]-xyz0[0];
1144 double d1 = xyz1[1]-xyz0[1];
1145 double d2 = xyz1[2]-xyz0[2];
1146 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1149 if (updTime) seed->AddTimeStep(ds);
1151 xpos = seed->GetX();
1153 if (AliDebugLevelClass()>1) {
1154 AliDebug(2,Form("After Propagate to X=%f",xToGo));
1155 seed->AliExternalTrackParam::Print();
1160 //______________________________________________________________________________
1161 Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
1163 // finalize single TPC track hypothesis
1164 if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
1165 AliITSUSeed* winner = 0;
1167 fCurrMass = hyp->GetMass();
1168 if (!(winner=hyp->DefineWinner(fCurrTrackCond->GetActiveLrInner()))) return kFALSE;
1171 #ifdef _ITSU_TUNING_MODE_ // fill tuning histos
1172 TObjArray* dest = hyp->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1174 AliITSUSeed* sd = winner;
1175 double htrPt = hyp->Pt();
1178 if ( (clID=sd->GetLrCluster(lrID))<0 ) continue;
1179 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Nrm ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2GloNrm());
1180 ((TH2*)dest->At( GetHistoID(lrID,kHBestInBranch,fCurrPassID) ))->Fill(htrPt,sd->GetOrdBranch());
1181 ((TH2*)dest->At( GetHistoID(lrID,kHBestInCand ,fCurrPassID) ))->Fill(htrPt,sd->GetOrdCand());
1183 if (dest==fCHistoArrFake && !sd->IsFake()) continue; // for the fake seeds fill only fake clusters part
1185 ((TH2*)dest->At( GetHistoID(lrID,kHResY ,fCurrPassID) ))->Fill(htrPt,sd->GetResidY());
1186 ((TH2*)dest->At( GetHistoID(lrID,kHResZ ,fCurrPassID) ))->Fill(htrPt,sd->GetResidZ());
1187 ((TH2*)dest->At( GetHistoID(lrID,kHResYP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullY());
1188 ((TH2*)dest->At( GetHistoID(lrID,kHResZP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullZ());
1189 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Cl ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2Cl());
1191 } while((sd=(AliITSUSeed*)sd->GetParent()));
1193 ((TH2*)dest->At( GetHistoID(-1,kHChiMatch,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSTPC());
1194 ((TH2*)dest->At( GetHistoID(-1,kHChiITSSA,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSSA());
1199 CheckClusterSharingConflicts(hyp);
1200 return hyp->GetWinner(); // winner might change of disappear after resolving conflicts
1204 //______________________________________________________________________________
1205 void AliITSUTrackerGlo::FinalizeHypotheses()
1207 // select winner for each hypothesis, remove cl. sharing conflicts
1208 AliCodeTimerAuto("",0);
1210 int* esdTrackIndex = fESDIndex.GetArray();
1211 for (int ih=0;ih<fNTracksESD;ih++) FinalizeHypothesis(GetTrackHyp(esdTrackIndex[ih])); // finlize in decreasing pt order
1213 for (int ih=fNTracksESD;ih--;) UpdateESDTrack(GetTrackHyp(esdTrackIndex[ih]),AliESDtrack::kITSin);
1217 //______________________________________________________________________________
1218 void AliITSUTrackerGlo::CheckClusterSharingConflicts(AliITSUTrackHyp* hyp)
1220 // remove eventual cluster sharing conflicts
1221 AliITSUSeed* winner = 0;
1222 if (!(winner=hyp->GetWinner())) return;
1223 UShort_t idH = (UShort_t)hyp->GetUniqueID()+1;
1225 int inLr = fCurrTrackCond->GetActiveLrInner();
1226 AliITSUSeed* winSD=winner;
1228 if ( (clID=winSD->GetLrCluster(lrID))<0 ) continue;
1229 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1230 Int_t refID = cl->GetRecoInfo(); // was it already referred by some track (with refID-1 if >0)
1231 if (!refID) {cl->SetRecoInfo(idH); continue;} // if not, refer from track to cluster
1232 if (refID==idH) continue; // ignore reference to itself
1233 // the cluster is already used by other track, need to resolve the conflict
1234 AliITSUTrackHyp* hypC = GetTrackHyp(refID-1); // competitor
1235 // printf("Ref to %d (%p) from %d Cl %d %d\n",refID-1,hypC, idH-1,clID,lrID);
1236 AliITSUSeed* winnerC = hypC->GetWinner();
1238 printf("Missing winner, ERROR\n");
1242 if (AliDebugLevelClass()>1) { //------------------------ DEBUG -------------------
1243 AliInfo(Form("Shared cl#%4d on Lr:%d: Hyp%3d/Q:%6.3f vs Hyp%3d/Q:%6.3f",
1244 clID,lrID,idH-1,winner->GetQualityVar(),refID-1,winnerC->GetQualityVar()));
1245 // dump winner of hyp
1246 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1247 idH-1,winner->Pt(),winner->GetChi2GloNrm(),winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),
1248 winner->GetNLayersHit(),hyp->GetITSLabel(),hyp->GetESDTrack()->GetTPCLabel());
1250 AliITSUSeed* sd = winner;
1253 int clIDt = sd->GetLrCluster(lrs);
1254 if (clIDt<0) continue;
1255 while( lrs>++prevL ) printf("%4s ","----");
1256 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1257 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1259 // dump winner of hypC
1260 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1261 refID-1,winnerC->Pt(),winnerC->GetChi2GloNrm(),winnerC->GetChi2ITSTPC(),winnerC->GetChi2ITSSA(),
1262 winnerC->GetNLayersHit(),hypC->GetITSLabel(),hypC->GetESDTrack()->GetTPCLabel());
1267 int clIDt = sd->GetLrCluster(lrs);
1268 if (clIDt<0) continue;
1269 while( lrs>++prevL ) printf("%4s ","----");
1270 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1271 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1273 } //-------------------------------------------------- DEBUG -------------------
1275 if (winner->GetQualityVar()<winnerC->GetQualityVar()) { // current tracks is better then competitor track
1276 FlagSeedClusters(winnerC,kFALSE,refID); // unflag cluster usage by loser
1277 cl->SetRecoInfo(idH);
1279 //if ( hypC->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",refID-1));
1281 if (hypC->DefineWinner(inLr)) { // find new winner instead of suppressed candidate
1283 CheckClusterSharingConflicts(hypC); // and check its sharing conflicts again
1284 if (winner->IsKilled()) break; // the current winner might have been killed during check of new winner of hypC!
1287 else { // competitor hypC is better than the hyp
1288 FlagSeedClusters(winner,kFALSE,idH); // unflag cluster usage by loser
1290 //if ( hyp->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",idH-1));
1291 if (hyp->DefineWinner(inLr)) {
1293 CheckClusterSharingConflicts(hyp);
1298 } while ((winSD=(AliITSUSeed*)winSD->GetParent()));
1303 //______________________________________________________________________________
1304 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1306 // update ESD track with current best hypothesis
1308 AliESDtrack* esdTr = hyp->GetESDTrack();
1309 if (!esdTr || esdTr->IsOn(flag)) return;
1310 AliITSUSeed* win = hyp->GetWinner();
1311 if (!win || win->IsKilled()) return;
1315 case AliESDtrack::kITSin:
1316 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1318 // set fakes cluster info
1320 UShort_t clfake = 0;
1321 AliITSUSeed* sd = win;
1323 if (sd->IsFake()) clfake |= 0x1<<sd->GetLayerID();
1324 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1326 // RS: Temporary set special flag for tracks from the afterburner
1327 if (fCurrPassID>0) clfake |= 0x1<<7;
1329 esdTr->SetITSSharedMap(clfake);
1333 case AliESDtrack::kITSout:
1334 // here the stored chi2 will correspond to backward ITS-SA fit
1335 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1337 // TODO: avoid setting friend
1340 case AliESDtrack::kITSrefit:
1342 // at the moment fill as a chi2 the TPC-ITS matching chi2
1343 chiSav = hyp->GetChi2();
1344 hyp->SetChi2(win->GetChi2ITSTPC());
1345 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1346 hyp->SetChi2(chiSav);
1348 // TODO: avoid setting cluster info
1351 AliFatal(Form("Unknown flag %d",flag));
1354 esdTr->SetITSLabel(hyp->GetITSLabel());
1355 // transfer module indices
1359 //______________________________________________________________________________
1360 Double_t AliITSUTrackerGlo::RefitTrack(AliExternalTrackParam* trc, Double_t rDest, Int_t stopCond)
1362 // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
1363 // if stopCond<0 : propagate till last cluster then stop
1364 // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
1365 // if stopCond>0 : rDest must be reached
1367 double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1368 int dir,lrStart,lrStop;
1370 dir = rCurr<rDest ? 1 : -1;
1371 lrStart = fITS->FindFirstLayerID(rCurr,dir);
1372 lrStop = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
1374 if (AliDebugLevelClass()>2) {
1375 printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1376 printf("Before refit: "); trc->Print();
1378 if (lrStop<0 || lrStart<0) AliFatal(Form("Failed to find start(%d) or last(%d) layers. Track from %.3f to %.3f",lrStart,lrStop,rCurr,rDest));
1380 Int_t clInfo[2*AliITSUAux::kMaxLayers];
1381 Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
1382 fCurrMass = fCurrHyp->GetMass();
1383 AliExternalTrackParam tmpTr(*trc);
1385 int iclLr[2],nclLr,clCount=0;
1387 int lrStop1 = lrStop+1;
1388 for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
1389 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1390 if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1391 int ilrA2,ilrA = lr->GetActiveID();
1392 // passive layer or active w/o hits will be traversed on the way to next cluster
1393 if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue;
1395 // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1397 if (dir>0) { // clusters are stored in increasing radius order
1398 iclLr[nclLr++]=clInfo[ilrA2++];
1399 if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
1402 if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
1403 iclLr[nclLr++]=clInfo[ilrA2];
1406 Bool_t transportedToLayer = kFALSE;
1407 for (int icl=0;icl<nclLr;icl++) {
1408 AliITSUClusterPix* clus = (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1409 AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1410 if (!tmpTr.Rotate(sens->GetPhiTF())) {
1411 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1412 AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
1416 double xClus = sens->GetXTF()+clus->GetX();
1417 if (!transportedToLayer) {
1418 if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1419 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1420 AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
1421 return -1; // go to the entrance to the layer
1424 transportedToLayer = kTRUE;
1427 if (AliDebugLevelClass()>1) {
1428 AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1431 if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1432 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1433 AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
1437 Double_t p[2]={clus->GetY(), clus->GetZ()};
1438 Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
1439 double chi2cl = tmpTr.GetPredictedChi2(p,cov);
1442 #ifdef _ITSU_TUNING_MODE_
1443 TObjArray* dest = trc->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1444 if (dest && fTrackPhaseID>kClus2Tracks) {
1446 double htrPt = tmpTr.Pt();
1447 double dy = p[0]-tmpTr.GetY();
1448 double dz = p[1]-tmpTr.GetZ();
1449 ((TH2*)dest->At( GetHistoID(ilrA,kHResY,0,fTrackPhaseID) ))->Fill(htrPt,dy);
1450 ((TH2*)dest->At( GetHistoID(ilrA,kHResZ,0,fTrackPhaseID) ))->Fill(htrPt,dz);
1451 double errY = tmpTr.GetSigmaY2();
1452 double errZ = tmpTr.GetSigmaZ2();
1453 if (errY>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResYP,0,fTrackPhaseID) ))->Fill(htrPt,dy/Sqrt(errY));
1454 if (errZ>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResZP,0,fTrackPhaseID) ))->Fill(htrPt,dz/Sqrt(errZ));
1455 ((TH2*)dest->At( GetHistoID(ilrA,kHChi2Cl,0,fTrackPhaseID) ))->Fill(htrPt,chi2cl);
1459 if ( !tmpTr.Update(p,cov) ) {
1460 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1461 AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
1464 if (AliDebugLevelClass()>1) {
1465 printf("AfterRefit: "); tmpTr.Print();
1467 if (++clCount==nCl && stopCond<0) {
1469 return chi2; // it was requested to not propagate after last update
1474 // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1475 // Still, try to go as close as possible to rDest.
1477 if (lrStart!=lrStop) {
1478 if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1479 AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1480 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1482 if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1483 AliDebug(1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1484 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1487 // go to the destination radius. Note that here we don't select direction to avoid precision problems
1488 if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
1489 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1492 if (AliDebugLevelClass()>2) {
1493 printf("After refit (now at lr %d): ",lrStop); trc->Print();
1498 //__________________________________________________________________
1499 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp)
1503 const int kMaxLbPerCl = 3;
1504 int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1505 Int_t lr,nLab=0,nCl=0;
1506 AliITSUSeed *seed = hyp->GetWinner();
1508 int clID = seed->GetLrCluster(lr);
1510 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1512 for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1513 int trLb = cl->GetLabel(imc);
1515 // search this mc track in already accounted ones
1517 for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1518 if (iLab<nLab) lbStat[iLab]++;
1523 } // loop over given cluster's labels
1525 seed = (AliITSUSeed*)seed->GetParent();
1526 } // loop over clusters
1528 AliESDtrack* esdTr = hyp->GetESDTrack();
1529 int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1531 int maxLab=0,nTPCok=0;
1532 for (int ilb=nLab;ilb--;) {
1533 int st = lbStat[ilb];
1534 if (lbStat[maxLab]<st) maxLab=ilb;
1535 if (tpcLab==lbID[ilb]) nTPCok += st;
1537 hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1538 hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1539 hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
1543 hyp->SetFakeRatio(-1.);
1544 hyp->SetLabel( -tpcLab );
1545 hyp->SetITSLabel( kDummyLabel );
1548 //__________________________________________________________________
1549 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1551 // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1552 if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
1553 AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
1554 fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
1555 AliITSUSeed** tmpArr = fLayerCandidates;
1556 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
1557 memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
1558 delete tmpArr; // delete only array, not objects
1560 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1561 int slot=fNBranchesAdded++;
1562 for (int slotF=slot;slotF--;) { // slotF is always slot-1
1563 AliITSUSeed* si = branches[slotF];
1564 if (si->Compare(seed)<0) break; // found the last seed with better quality
1565 // otherwise, shift the worse seed to the next slot
1566 branches[slot] = si;
1567 slot = slotF; // slot should be slotF+1
1569 // if needed, move worse seeds
1570 branches[slot] = seed;
1575 //__________________________________________________________________
1576 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1578 // keep allowed number of branches for current seed and disable extras
1579 AliCodeTimerAuto("",0);
1580 int nb = Min(fNBranchesAdded,acceptMax);
1581 // if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1582 // disable unused branches
1583 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1585 #ifdef _ITSU_TUNING_MODE_
1586 for (int ib=0;ib<fNBranchesAdded;ib++) branches[ib]->SetOrdBranch(ib);
1589 for (int ib=nb;ib<fNBranchesAdded;ib++) MarkSeedFree(branches[ib]);
1590 fNCandidatesAdded += nb; // update total candidates counter
1591 fNBranchesAdded = 0; // reset branches counter
1595 //__________________________________________________________________
1596 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1598 // transfer allowed number of branches to hypothesis container
1600 AliCodeTimerAuto("",0);
1601 // sort candidates in increasing order of chi2
1602 static int lastSize = 0;
1603 static int *index = 0;
1604 static float *chi2 = 0;
1605 if (fLayerMaxCandidates>lastSize) {
1606 lastSize = fLayerMaxCandidates;
1609 index = new int[lastSize];
1610 chi2 = new float[lastSize];
1612 for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1613 Sort(fNCandidatesAdded,chi2,index,kFALSE);
1616 if (ilr>0) { // just take 1st acceptMax candidates
1617 nb = Min(fNCandidatesAdded,acceptMax);
1618 for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1620 else { // for innermost layer test ITS SA backward chi2 and matching to TPC
1621 AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
1622 for (nb=0;nb<fNCandidatesAdded;nb++) {
1623 AliITSUSeed* sd = fLayerCandidates[index[nb]];
1624 fCurrHyp->SetWinner(sd);
1625 if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
1627 AddProlongationHypothesis(sd,ilr);
1628 if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
1631 fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
1633 // disable unused candidates
1634 for (int ib=nb;ib<fNCandidatesAdded;ib++) MarkSeedFree(fLayerCandidates[index[ib]]);
1635 fNCandidatesAdded = 0; // reset candidates counter
1639 //__________________________________________________________________
1640 Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
1642 // check seed backward propagation chi2 and matching to TPC
1643 double bz0 = GetBz();
1644 double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
1645 static AliExternalTrackParam trback;
1647 trback.ResetCovariance(10000);
1648 int ndf = seed->GetNLayersHit()*2-5;
1649 double chi2sa = RefitTrack(&trback,rDest,1);
1650 if (chi2sa<0) return kFALSE;
1651 if (ndf>0) chi2sa /= ndf;
1652 if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2()) return kFALSE;
1654 // relate to TPC track at outer layer
1655 AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
1656 if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
1657 if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
1658 AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
1664 double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
1665 if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
1667 seed->SetChi2ITSSA(chi2sa);
1668 seed->SetChi2ITSTPC(chi2Match);
1672 //__________________________________________________________________
1673 void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
1675 // remove those hypothesis seeds which dont lead to candidates at final layer
1677 // 1st, flag the seeds to save
1678 int lr0 = fCurrTrackCond->GetActiveLrInner();
1679 for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
1680 AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
1681 if (!seed) continue;
1682 seed->FlagTree(AliITSUSeed::kSave);
1683 dest->AddSeed(seed,lr0);
1685 for (int ilr=0;ilr<fNLrActive;ilr++) {
1686 int nsd = fCurrHyp->GetNSeeds(ilr);
1687 for (int isd=0;isd<nsd;isd++) {
1688 AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
1689 if (!seed) continue; // already discarded or saved
1690 if (seed->IsSaved()) dest->AddSeed(seed,ilr);
1691 else MarkSeedFree(seed);
1695 // AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
1698 //__________________________________________________________________
1699 void AliITSUTrackerGlo::CleanHypothesis(AliITSUTrackHyp* hyp)
1703 for (int ilr=0;ilr<fNLrActive;ilr++) {
1704 int nsd = hyp->GetNSeeds(ilr);
1705 for (int isd=0;isd<nsd;isd++) {
1706 AliITSUSeed* seed = hyp->RemoveSeed(ilr,isd);
1707 if (!seed) continue; // already discarded or saved
1713 //__________________________________________________________________
1714 void AliITSUTrackerGlo::FlagSplitClusters()
1716 // set special bit on split clusters using MC info
1717 for (int ilr=fNLrActive;ilr--;) {
1719 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1720 for (int isn=lr->GetNSensors();isn--;) {
1721 AliITSURecoSens* sens = lr->GetSensor(isn);
1722 int nCl = sens->GetNClusters();
1724 int cl0 = sens->GetFirstClusterId();
1725 for (int icl=nCl;icl--;) {
1726 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1727 for (int icl1=icl;icl1--;) {
1728 AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1729 if (cl->HasCommonTrack(cl1)) {
1730 if (!cl->IsSplit()) nsplit++;
1731 if (!cl1->IsSplit()) nsplit++;
1738 AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1743 //__________________________________________________________________
1744 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1746 // check if the seed contains split cluster with size < maxSize
1748 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1749 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1750 if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1752 return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1756 //__________________________________________________________________
1757 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1759 // print seeds clusters
1761 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1762 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1765 if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1769 //__________________________________________________________________
1770 void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg, UShort_t ref)
1772 // mark used clusters
1775 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1776 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1777 if (ref) { // do we need to set or delete cluster->track ref?
1779 if (!cl->GetRecoInfo()) cl->SetRecoInfo(ref); // set ref only if cluster already does not refer to other track, inc.counter
1782 if (cl->GetRecoInfo()==ref) cl->SetRecoInfo(0); // unset reference only if it refers to ref, decrease counter
1786 seed = (AliITSUSeed*)seed->GetParent();
1791 //__________________________________________________________________
1792 void AliITSUTrackerGlo::CheckClusterUsage()
1794 // check cluster usage info
1795 for (int ilr=0;ilr<fNLrActive;ilr++) {
1796 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1797 int ncl = lr->GetNClusters();
1799 int nclUs0CntShr1 = 0;
1803 for (int icl=0;icl<ncl;icl++) {
1804 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(icl);
1805 int usc = cl->GetClUsage();
1807 if (cl->IsClusterUsed() || cl->IsClusterShared()) {
1809 printf("Lr%d Cl%4d: UseCnt=0 but UseFlg=%d UseShr=%d\n",ilr,icl,cl->IsClusterUsed(),cl->IsClusterShared());
1814 if (!cl->IsClusterUsed()) {
1816 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1818 if (cl->IsClusterShared()) {
1820 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1824 printf("Lr%d Cl%4d: UseCnt=%d!! but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1827 printf("ClUsage Lr%d: Ncl=%5d Nusd=%5d (%5.2f) Use0CS1:%4d Use1C0:%4d Use1S1:%4d UseL:%d\n",
1828 ilr,ncl,nclUsed,ncl? float(nclUsed)/ncl : 0, nclUs0CntShr1,nclUs1Cnt0,nclUs1Shr1,nclUseL);
1833 #ifdef _ITSU_TUNING_MODE_
1834 //__________________________________________________________________
1835 TObjArray* AliITSUTrackerGlo::BookControlHistos(const char* pref)
1837 // book special control histos
1838 TString prefS = pref;
1840 TObjArray* dest = new TObjArray;
1841 dest->SetOwner(kTRUE);
1843 const int kNResDef=7;
1844 const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1845 const double ptMax=10;
1846 const double plMax=10;
1847 const double chiMax=100;
1848 const int nptbins=50;
1849 const int nresbins=400;
1850 const int nplbins=50;
1851 const int nchbins=200;
1852 const int maxBr = 15;
1853 const int maxCand = 200;
1855 for (int phase=0;phase<kNTrackingPhases;phase++) {
1856 for (int pass=0;pass<AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();pass++) {
1858 for (int ilr=0;ilr<fNLrActive;ilr++) {
1860 // ----------------- These are histos to be filled in Cluster2Tracks of each pass.
1861 // PropagateBack and RefitInward will be stored among the histos of 1st pass
1862 if (pass>0 && phase!=kClus2Tracks) continue;
1864 double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1865 ttl = Form("Pass%d_S%d_residY%d_%s",pass,phase,ilr,pref);
1866 TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1867 dest->AddAtAndExpand(hdy,GetHistoID(ilr,kHResY,pass,phase));
1868 hdy->SetDirectory(0);
1870 ttl = Form("Pass%d_S%d_residYPull%d_%s",pass,phase,ilr,pref);
1871 TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1872 dest->AddAtAndExpand(hdyp,GetHistoID(ilr,kHResYP,pass,phase));
1873 hdyp->SetDirectory(0);
1875 ttl = Form("Pass%d_S%d_residZ%d_%s",pass,phase,ilr,pref);
1876 TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1877 dest->AddAtAndExpand(hdz,GetHistoID(ilr,kHResZ,pass,phase));
1878 hdz->SetDirectory(0);
1880 ttl = Form("Pass%d_S%d_residZPull%d_%s",pass,phase,ilr,pref);
1881 TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1882 hdzp->SetDirectory(0);
1883 dest->AddAtAndExpand(hdzp,GetHistoID(ilr,kHResZP,pass,phase));
1885 ttl = Form("Pass%d_S%d_chi2Cl%d_%s",pass,phase,ilr,pref);
1886 TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1887 hchi->SetDirectory(0);
1888 dest->AddAtAndExpand(hchi,GetHistoID(ilr,kHChi2Cl,pass,phase));
1890 // ------------------- These histos are filled for Clusters2Tracks only
1891 if (phase!=kClus2Tracks) continue;
1893 ttl = Form("Pass%d_S%d_chi2Nrm%d_%s",pass,phase,ilr,pref);
1894 TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1895 hchiN->SetDirectory(0);
1896 dest->AddAtAndExpand(hchiN,GetHistoID(ilr,kHChi2Nrm,pass,phase));
1898 ttl = Form("Pass%d_S%d_bestInBranch%d_%s",pass,phase,ilr,pref);
1899 TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
1900 hnbr->SetDirectory(0);
1901 dest->AddAtAndExpand(hnbr,GetHistoID(ilr,kHBestInBranch,pass,phase));
1903 ttl = Form("Pass%d_S%d_bestInCands%d_%s",pass,phase,ilr,pref);
1904 TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
1905 hncn->SetDirectory(0);
1906 dest->AddAtAndExpand(hncn,GetHistoID(ilr,kHBestInCand,pass,phase));
1908 } // loop over layers
1910 // these are special histos, filled not per layer but in the end of track fit in Clusters2Tracks in EVERY pass
1912 if (phase==kClus2Tracks) {
1914 ttl = Form("Pass%d_Chi2Match_%s",pass,pref);
1915 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1916 hchiMatch->SetDirectory(0);
1917 dest->AddAtAndExpand(hchiMatch,GetHistoID(-1,kHChiMatch,pass,phase));
1920 ttl = Form("Pass%d_Chi2ITSSA_%s",pass,pref);
1921 hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
1922 hchiSA->SetDirectory(0);
1923 dest->AddAtAndExpand(hchiSA,GetHistoID(-1,kHChiITSSA,pass,phase));
1926 } // loop over tracking passes
1927 }// loop over tracking phases
1932 //__________________________________________________________________
1933 Int_t AliITSUTrackerGlo::GetHistoID(Int_t lr, Int_t hid, Int_t pass, Int_t phase)
1935 // get id for the requested histo
1938 return pass*kHistosPass + phase*kHistosPhase + lr*kMaxHID + hid;