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)
58 ,fCurrESDtrMClb(kDummyLabel)
60 ,fCountProlongationTrials(0)
66 ,fLayerMaxCandidates(1000)
72 ,fSeedsPool("AliITSUSeed",0)
85 #ifdef _ITSU_TUNING_MODE_
90 // Default constructor
94 //_________________________________________________________________________
95 AliITSUTrackerGlo::~AliITSUTrackerGlo()
99 delete[] fLayerCandidates;
100 if (fWorkHyp) fWorkHyp->SetTPCSeed(0); // this hypothesis does not own the seed
104 #ifdef _ITSU_TUNING_MODE_
105 if (fCHistoArrCorr || fCHistoArrFake) {
106 TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
108 AliInfo("Storing control histos");
109 // ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
110 if (fCHistoArrCorr) {
111 fCHistoArrCorr->Write();
112 delete fCHistoArrCorr;
114 if (fCHistoArrFake) {
115 fCHistoArrFake->Write();
116 delete fCHistoArrFake;
127 //_________________________________________________________________________
128 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
130 // init with external reconstructor
132 fITS = rec->GetITSInterface();
133 fNLrActive = fITS->GetNLayersActive();
134 fWorkHyp = new AliITSUTrackHyp(fNLrActive);
135 // create material lookup table
136 const int kNTest = 1000;
137 const double kStepsPerCM=5;
138 fMatLUT = new AliITSUMatLUT(fITS->GetRMin(),fITS->GetRMax(),Nint(kStepsPerCM*(fITS->GetRMax()-fITS->GetRMin())));
140 for (int ilr=fITS->GetNLayers();ilr--;) {
141 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
142 if (zmn>Abs(lr->GetZMin())) zmn = Abs(lr->GetZMin());
143 if (zmn>Abs(lr->GetZMax())) zmn = Abs(lr->GetZMax());
145 fMatLUT->FillData(kNTest,-zmn,zmn);
147 if (fLayerMaxCandidates<1) fLayerMaxCandidates = 1000;
148 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
149 fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
150 fFreeSeedsID.Set(1000);
155 //_________________________________________________________________________
156 void AliITSUTrackerGlo::CreateDefaultTrackCond()
158 // creates default tracking conditions to be used when recoparam does not provide them
160 AliITSUTrackCond* cond = new AliITSUTrackCond();
162 cond->SetNLayers(fNLrActive);
163 cond->AddNewCondition(fNLrActive);
164 cond->AddGroupPattern( 0xffff, fNLrActive ); // require all layers hit
167 fDefTrackConds.AddLast(cond);
169 AliInfo("Created conditions: ");
170 for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
175 //_________________________________________________________________________
176 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
179 AliCodeTimerAuto("",0);
180 SetTrackingPhase(kClus2Tracks);
182 #ifdef _ITSU_TUNING_MODE_
183 if (!fCHistoArrCorr) fCHistoArrCorr = BookControlHistos("Corr");
184 if (!fCHistoArrFake) fCHistoArrFake = BookControlHistos("Fake");
187 static TArrayF esdTrPt(fESDIndex.GetSize());
189 TObjArray *trackConds = 0;
191 fCountProlongationTrials = 0;
197 fNTracksESD = esdEv->GetNumberOfTracks();
198 AliInfo(Form("Will try to find prolongations for %d tracks",fNTracksESD));
199 int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
200 fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
202 if (!fDefTrackConds.GetEntriesFast()) {
203 AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
204 CreateDefaultTrackCond();
206 trackConds = &fDefTrackConds;
207 nTrackCond = trackConds->GetEntriesFast();
209 else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
211 static Bool_t first = kTRUE;
213 AliITSUReconstructor::GetRecoParam()->Print();
217 if (fHypStore.GetSize()<fNTracksESD) fHypStore.Expand(fNTracksESD+100);
219 fITS->ProcessClusters();
221 #ifdef _ITSU_TUNING_MODE_
222 FlagSplitClusters(); // tmp RS
225 // the tracks will be reconstructed in decreasing pt order, sort them
226 if (fESDIndex.GetSize()<fNTracksESD) {
227 fESDIndex.Set(fNTracksESD+200);
228 esdTrPt.Set(fNTracksESD+200);
230 int* esdTrackIndex = fESDIndex.GetArray();
231 float* trPt = esdTrPt.GetArray();
232 for (int itr=fNTracksESD;itr--;) {
233 AliESDtrack* esdTr = esdEv->GetTrack(itr);
234 esdTr->SetStatus(AliESDtrack::kITSupg); // Flag all tracks as participating in the ITSup reco
235 trPt[itr] = esdTr->Pt();
237 Sort(fNTracksESD,trPt,esdTrackIndex,kTRUE);
239 for (int icnd=0;icnd<nTrackCond;icnd++) {
241 fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
242 if (!fCurrTrackCond->IsInitDone()) fCurrTrackCond->Init();
243 // select ESD tracks to propagate
244 for (int itr=0;itr<fNTracksESD;itr++) {
245 int trID = esdTrackIndex[itr];
246 fCurrESDtrack = esdEv->GetTrack(trID);
247 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
250 if (!NeedToProlong(fCurrESDtrack,trID)) continue; // are we interested in this track?
252 // if specific tracks should be checked, create a global array int wtc[] = {evITS*10000+trID, ...};
254 int nwtc = sizeof(wtc)/sizeof(int);
256 for (int k=0;k<nwtc;k++) {
257 if (wtc[k]==evID*10000+trID) {
259 printf("At esdTr: %d MC: %d\n",wtc[k],fCurrESDtrMClb);
263 AliLog::SetClassDebugLevel("AliITSUTrackerGlo",dbg ? 10:0);
266 AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMassForTracking(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
268 FindTrack(fCurrESDtrack, trID);
272 if (AliDebugLevelClass()>+2) {
273 AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
277 FinalizeHypotheses();
278 #ifdef _ITSU_TUNING_MODE_
279 CheckClusterUsage(); //!!RS
282 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0); // in case wtc array was defined, uncomment this
285 AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,fNTracksESD));
290 //_________________________________________________________________________
291 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
294 // Do outward fits in ITS
295 AliCodeTimerAuto("",0);
297 SetTrackingPhase(kPropBack);
298 fNTracksESD = esdEv->GetNumberOfTracks();
299 fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
301 AliDebug(1,Form("Will propagate back %d tracks",fNTracksESD));
304 double bz0 = GetBz();
305 Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
306 AliITSUTrackHyp dummyTr;
307 const double kWatchStep=10.; // for larger steps watch arc vs segment difference
308 Double_t times[AliPID::kSPECIESC];
310 for (int itr=0;itr<fNTracksESD;itr++) {
311 fCurrESDtrack = esdEv->GetTrack(itr);
312 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
313 // Start time integral and add distance from current position to vertex
314 if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
315 if (!fCurrESDtrack->IsOn(AliESDtrack::kTPCin)) continue; // skip ITS s.a.
317 fCurrESDtrack->GetXYZ(xyzTrk);
318 Double_t dst = 0.; // set initial track lenght, tof
320 double dxs = xyzTrk[0] - xyzVtx[0];
321 double dys = xyzTrk[1] - xyzVtx[1];
322 // double dzs = xyzTrk[2] - xyzVtx[2];
323 // RS: for large segment steps d use approximation of cicrular arc s by
324 // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
325 // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
326 // Hence s^2/d^2 = (1+1/6 p^2)^2
327 dst = dxs*dxs + dys*dys;
328 if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
329 double crv = Abs(fCurrESDtrack->GetC(bz0));
330 double fcarc = 1.+crv*crv*dst/6.;
333 // RS: temporary hack since we don't have SPD vertex:
335 dst *= 1+fCurrESDtrack->GetTgl()*fCurrESDtrack->GetTgl();
339 fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
341 if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
342 dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
343 dummyTr.StartTimeIntegral();
344 dummyTr.AddTimeStep(dst);
345 dummyTr.GetIntegratedTimes(times,AliPID::kSPECIESC);
346 fCurrESDtrack->SetIntegratedTimes(times);
347 fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
351 fCurrHyp = GetTrackHyp(itr);
352 fCurrMass = fCurrHyp->GetMass();
353 fCurrHyp->StartTimeIntegral();
354 fCurrHyp->AddTimeStep(dst);
355 fCurrHyp->ResetCovariance(10000);
357 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax(),nclFit);
358 if (chi2>0) { // propagate to exit from the ITS/TPC screen
359 int ndf = nclFit*2-5;
360 if (ndf>0) chi2 /= ndf;
361 fCurrHyp->SetChi2(chi2);
362 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSout);
366 AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
368 //fCurrHyp->AliExternalTrackParam::Print();
369 //fCurrHyp->GetWinner()->Print();
374 AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
375 fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
380 //_________________________________________________________________________
381 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
384 // refit the tracks inward, using current cov.matrix
385 AliCodeTimerAuto("",0);
387 SetTrackingPhase(kRefitInw);
388 fUseMatLUT = AliITSUReconstructor::GetRecoParam()->GetUseMatLUT(fTrackPhaseID);
389 fNTracksESD = esdEv->GetNumberOfTracks();
390 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
392 AliDebug(1,Form("Will refit inward %d tracks",fNTracksESD));
393 // Bool_t uselogMS = AliExternalTrackParam::GetUseLogTermMS();
394 // AliExternalTrackParam::SetUseLogTermMS(kTRUE);
396 for (int itr=0;itr<fNTracksESD;itr++) {
397 fCurrESDtrack = esdEv->GetTrack(itr);
398 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
399 // Start time integral and add distance from current position to vertex
400 UInt_t trStat = fCurrESDtrack->GetStatus();
401 if ( !(trStat & AliESDtrack::kITSout) ) continue;
402 if ( trStat & AliESDtrack::kITSrefit ) continue; // already done
403 if ( !(trStat & AliESDtrack::kTPCin) ) continue; // skip ITS s.a.
404 if ( (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
406 fCurrHyp = GetTrackHyp(itr);
407 fCurrMass = fCurrHyp->GetMass();
408 *fCurrHyp = *fCurrESDtrack; // fetch current ESDtrack kinematics
411 // fCurrHyp->ResetCovariance(1000); // in case we want to do SA refit
412 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin(),nclFit);
413 if (chi2>0) { // propagate up to inside radius of the beam pipe
414 int ndf = nclFit*2-5;
415 if (ndf>0) chi2 /= ndf;
416 fCurrHyp->SetChi2(chi2);
417 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
421 AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
425 // AliExternalTrackParam::SetUseLogTermMS(uselogMS);
427 AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
428 fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,fNTracksESD));
430 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
434 //_________________________________________________________________________
435 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
437 // read from tree (if pointer provided) or directly from the ITS reco interface
439 return fReconstructor->LoadClusters(treeRP);
442 //_________________________________________________________________________
443 void AliITSUTrackerGlo::UnloadClusters()
446 // Remove clusters from the memory
448 AliITSURecoDet *det=fReconstructor->GetITSInterface();
449 Int_t nlayers=det->GetNLayersActive();
450 for (Int_t i=0; i<nlayers; i++) {
451 TClonesArray *clusters=*(det->GetLayerActive(i)->GetClustersAddress());
456 //_________________________________________________________________________
457 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
463 Info("GetCluster","To be implemented");
467 //_________________________________________________________________________
468 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr, Int_t esdID)
470 // do we need to match this track to ITS?
472 static double bz = GetBz();
474 // is it already prolonged
475 if (esdTr->IsOn(AliESDtrack::kITSin)) return kFALSE;
476 AliITSUTrackHyp* hyp = GetTrackHyp(esdID); // in case there is unfinished hypothesis from prev.pass, clean it
478 CleanHypothesis(hyp);
479 if (hyp->GetSkip()) return kFALSE; // need to skip this
483 if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
484 esdTr->IsOn(AliESDtrack::kTPCout) ||
485 esdTr->IsOn(AliESDtrack::kITSin) ||
486 esdTr->GetKinkIndex(0)>0) return kFALSE;
488 if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
491 esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz);
492 // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
493 if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation())
494 && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
495 Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
497 // if (esdTr->Pt()<3) return kFALSE;//RS
501 //_________________________________________________________________________
502 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
504 // find prolongaion candidates finding for single seed
506 AliITSUSeed seedUC; // copy of the seed from the upper layer
507 AliITSUSeed seedT; // transient seed between the seedUC and new prolongation hypothesis
509 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
510 AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID); // initialize prolongations hypotheses tree
511 if (!hypTr || hypTr->GetSkip()) return;
515 fCurrHyp->InitFrom(hypTr);
517 AliITSURecoSens *hitSens[AliITSURecoLayer::kMaxSensMatching];
519 int ilaUp = fNLrActive; // previous active layer really checked (some may be excluded!)
521 for (int ila=fNLrActive;ila--;ilaUp=fCurrActLrID) {
522 if (fCurrTrackCond->IsLayerExcluded(ila)) continue;
524 fCurrLayer = fITS->GetLayerActive(ila);
525 Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
527 // for the outermost layer the seed is created from the ESD track
528 int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
529 int maxNBranches = fCurrTrackCond->GetMaxBranches(ila);
530 int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
532 for (int isd=0;isd<nSeedsUp;isd++) {
534 if (ilaUp==fNLrActive) {
536 seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
539 seedU = fCurrHyp->GetSeed(ilaUp,isd); // seed on prev.active layer to prolong
540 if (seedU->IsKilled()) continue;
541 seedUC = *seedU; // its copy will be prolonged
542 seedUC.SetParent(seedU);
544 seedUC.ResetFMatrix(); // reset the matrix for propagation to next layer
545 // go till next active layer
547 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()));
549 if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
552 AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
554 // Check if the seed satisfies to track definition
555 if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
556 continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
558 if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
559 if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
565 if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
566 mcquest = fCurrESDtrMClb;
568 printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
571 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest); // find detectors which may be hit by the track
573 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens); // find detectors which may be hit by the track
575 AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
580 for (int isn=nsens;isn--;) {
581 AliITSURecoSens* sens = hitSens[isn];
582 int ncl = sens->GetNClusters();
586 // We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
587 // since the transport matrix should be defined in this frame.
588 double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
589 if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
591 if (AliDebugLevelClass()>2) {
592 printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
598 if (xs<seedT.GetX()) {
599 if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
601 else { // some low precision tracks may hit the sensor plane outside of the layer radius
603 if (AliDebugLevelClass()>2) {
604 if (!seedT.ContainsFake()) {
605 printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
611 if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
613 // if (!seedT.PropagateToX(xs,bz)) continue;
614 // if (!seedT.Rotate(sens->GetPhiTF())) continue;
615 if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
617 int clID0 = sens->GetFirstClusterId();
618 for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
619 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
620 int clID = clID0 + icl;
621 if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
622 int res = CheckCluster(&seedT,ila,clID0+icl);
623 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
625 if (res==kStopSearchOnSensor) break; // stop looking on this sensor
626 if (res==kClusterNotMatching) continue; // cluster does not match
627 // cluster is matching and it was added to the hypotheses tree
630 // cluster search is done. Do we need to have a version of this seed skipping current layer
631 if (!NeedToKill(&seedUC,kMissingCluster)) {
632 AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
633 if (nsens) { // to do: make penalty to account for probability to miss the cluster for good reason
634 double penalty = -fCurrTrackCond->GetMissPenalty(ila);
635 seedSkp->SetChi2Cl(penalty);
637 if (seedSkp->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(ila)) {
638 MarkSeedFree(seedSkp);
640 else AddSeedBranch(seedSkp);
642 // transfer the new branches of the seed to the hypothesis container
643 if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
646 if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
647 // ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
649 if (AliDebugLevelClass()>2) { //RS
650 printf(">>> All hypotheses on lr %d: \n",ila);
651 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
656 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
657 if (ila!=0) continue;
658 double vecL[5] = {0};
659 double matL[15] = {0};
660 AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
661 while(sp->GetParent()) {
662 sp->Smooth(vecL,matL);
663 if (sp->GetLayerID()>=fNLrActive-1) break;
664 sp = (AliITSUSeed*)sp->GetParent();
669 SaveReducedHypothesesTree(hypTr);
671 if (AliDebugLevelClass()>1) {
672 printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
677 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
681 //_________________________________________________________________________
682 AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
684 // init prolongaion candidates finding for single seed
685 AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
688 fCountProlongationTrials++;
690 fCurrMass = esdTr->GetMassForTracking();
692 hyp = new AliITSUTrackHyp(fNLrActive);
693 hyp->SetESDTrack(esdTr);
694 hyp->SetUniqueID(esdID);
695 hyp->SetMass(fCurrMass);
696 hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
697 SetTrackHyp(hyp,esdID);
698 Bool_t lut = fUseMatLUT; // for propagation from TPC use TGeo material budget always
700 if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
707 //_________________________________________________________________________
708 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
710 // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
713 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
715 int dir = lTo > lFrom ? 1:-1;
716 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
717 Bool_t checkFirst = kTRUE;
718 Bool_t limReached = kFALSE;
721 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
724 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
725 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
727 // go the entrance of the layer, assuming no materials in between
728 double xToGo = lrTo->GetR(-dir);
731 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
734 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
737 // double xts = xToGo;
738 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
739 //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
740 //seed->Print("etp");
745 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
747 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
748 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
749 //seed->Print("etp");
753 if (limReached) break;
759 //_________________________________________________________________________
760 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
762 // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
764 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
766 int dir = lTo > lFrom ? 1:-1;
767 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
768 Bool_t checkFirst = kTRUE;
769 Bool_t limReached = kFALSE;
772 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
775 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
776 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
778 // go the entrance of the layer, assuming no materials in between
779 double xToGo = lrTo->GetR(-dir);
782 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
785 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
788 // double xts = xToGo;
789 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
790 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
791 // seed->Print("etp");
795 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
797 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
798 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
799 //seed->Print("etp");
803 if (limReached) break;
809 //_________________________________________________________________________
810 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
812 // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
814 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
816 int dir = lTo > lFrom ? 1:-1;
817 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
818 Bool_t checkFirst = kTRUE;
821 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
824 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
825 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
827 // go the entrance of the layer, assuming no materials in between
828 double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
830 // double xts = xToGo;
831 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
832 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
833 // seed->Print("etp");
836 if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
839 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
841 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
842 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
843 //seed->Print("etp");
852 //_________________________________________________________________________
853 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
855 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
856 // If check is requested, do this only provided the track has not exited the layer already
857 double xToGo = lr->GetR(dir);
858 if (check) { // do we need to track till the surface of the current layer ?
859 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
861 AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
863 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
864 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
866 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
867 // go via layer to its boundary, applying material correction.
868 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
874 //_________________________________________________________________________
875 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
877 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
878 // If check is requested, do this only provided the track has not exited the layer already
879 double xToGo = lr->GetR(dir);
880 if (check) { // do we need to track till the surface of the current layer ?
881 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
882 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
883 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
886 AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
888 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
889 // go via layer to its boundary, applying material correction.
890 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
897 //_________________________________________________________________________
898 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
900 // go to the entrance of lr in direction dir, w/o applying material corrections.
901 // If check is requested, do this only provided the track did not reach the layer already
902 double xToGo = lr->GetR(-dir);
903 if (check) { // do we need to track till the surface of the current layer ?
904 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
905 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
906 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
908 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
909 // go via layer to its boundary, applying material correction.
910 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
915 //_________________________________________________________________________
916 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
918 // go to the entrance of lr in direction dir, w/o applying material corrections.
919 // If check is requested, do this only provided the track did not reach the layer already
920 double xToGo = lr->GetR(-dir);
921 if (check) { // do we need to track till the surface of the current layer ?
922 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
923 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
924 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
926 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
927 // go via layer to its boundary, applying material correction.
928 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
933 //_________________________________________________________________________
934 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
936 // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
937 // as well as some aux info
939 AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
940 seed->GetXYZ(&fTrImpData[kTrXIn]); // lab position at the entrance from above
941 static AliExternalTrackParam sc; // seed copy for manipulations
944 fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
945 if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
946 BringTo02Pi(fTrImpData[kTrPhiIn]);
947 double dr = lrA->GetDR(); // approximate X dist at the inner radius
948 if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
949 // special case: track does not reach inner radius, might be tangential
950 double r = sc.GetD(0,0,bz);
952 if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
953 dr = Abs(sc.GetX() - x);
954 if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
957 fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
958 BringTo02Pi(fTrImpData[kTrPhiOut]);
959 double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
960 double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
961 sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
962 sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
963 double phi0 = MeanPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
964 double dphi0 = DeltaPhiSmall(fTrImpData[kTrPhiOut],fTrImpData[kTrPhiIn]);
966 fTrImpData[kTrPhi0] = phi0;
967 fTrImpData[kTrZ0] = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
968 dphi0 += sgy/lrA->GetR();
969 fTrImpData[kTrDPhi] = dphi0<PiOver2() ? dphi0 : PiOver2();
970 fTrImpData[kTrDZ] = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn]) + sgz;
975 //________________________________________
976 void AliITSUTrackerGlo::ResetSeedsPool()
978 // mark all seeds in the pool as unused
979 AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
981 fSeedsPool.Clear(); // seeds don't allocate memory
985 //________________________________________
986 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd)
988 // account that this seed is "deleted"
989 int id = sd->GetPoolID();
991 AliError(Form("Freeing of seed %p NOT from the pool is requested",sd));
994 // AliInfo(Form("%d %p",id, seed));
995 fSeedsPool.RemoveAt(id);
996 if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
997 fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
1000 //_________________________________________________________________________
1001 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID)
1003 // Check if the cluster (in tracking frame!) is matching to track.
1004 // The track must be already propagated to sensor tracking frame.
1005 // Returns: kStopSearchOnSensor if the search on given sensor should be stopped,
1006 // kClusterMatching if the cluster is matching
1007 // kClusterMatching otherwise
1009 const double kTolerX = 5e-4;
1010 // AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1011 AliCluster *cl = fCurrLayer->GetCluster(clID);
1013 Bool_t goodCl = kFALSE;
1014 int currLabel = Abs(fCurrESDtrMClb);
1016 if (cl->GetLabel(0)>=0) {
1017 for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
1020 if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
1021 if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
1024 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1025 AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1026 track->Print("etp");
1030 return kStopSearchOnSensor; // propagation failed, seedT is intact
1033 double dy = cl->GetY()-track->GetY();
1034 double dz = cl->GetZ()-track->GetZ();
1037 double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
1038 fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
1039 if (dy2>tol2) { // the clusters are sorted in Z(col) then in Y(row).
1042 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1043 AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1044 track->Print("etp");
1046 PrintSeedClusters(track);
1049 // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
1050 // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
1051 if (dy<0) return kStopSearchOnSensor; // No chance that other cluster of this sensor will match (all Y's will be even larger)
1052 else return kClusterNotMatching; // Other clusters may match
1055 tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
1056 fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
1059 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1060 AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1061 track->Print("etp");
1063 PrintSeedClusters(track);
1066 return kClusterNotMatching; // Other clusters may match
1070 Double_t p[2]={cl->GetY(), cl->GetZ()};
1071 Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
1072 double chi2 = track->GetPredictedChi2(p,cov);
1074 if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
1076 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1077 AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1078 lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1079 track->Print("etp");
1081 PrintSeedClusters(track);
1084 return kClusterNotMatching;
1087 track = NewSeedFromPool(track); // input track will be reused, use its clone for updates
1088 if (!track->Update()) {
1090 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1091 AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
1092 track->Print("etp");
1094 PrintSeedClusters(track);
1097 MarkSeedFree(track);
1098 return kClusterNotMatching;
1100 track->SetChi2Cl(chi2);
1101 track->SetLrClusterID(lr,clID);
1103 if (track->GetChi2GloNrm()>fCurrTrackCond->GetMaxChi2GloNrm(lr)) {
1105 if (AliDebugLevelClass()>2 && goodCl && !track->ContainsFake()) {
1106 AliDebug(2,Form("Lost good cl on L:%d : Chi2Glo=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
1107 lr,track->GetChi2GloNrm(),fCurrTrackCond->GetMaxChi2GloNrm(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
1108 track->Print("etp");
1110 PrintSeedClusters(track);
1113 MarkSeedFree(track);
1114 return kClusterNotMatching;
1116 // cl->IncreaseClusterUsage(); // do this only for winners
1118 track->SetFake(!goodCl);
1121 AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
1122 goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
1125 AddSeedBranch(track);
1127 return kClusterMatching;
1130 //_________________________________________________________________________
1131 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
1133 // check if the seed should not be discarded
1134 const UShort_t kMask = 0xffff;
1135 if (flag==kMissingCluster) {
1136 int lastChecked = seed->GetLayerID();
1137 UShort_t patt = seed->GetHitsPattern();
1138 if (lastChecked) patt |= (~(kMask<<lastChecked))&fCurrTrackCond->GetAllowedLayers(); // not all layers were checked, complete unchecked ones by potential hits
1139 Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
1145 //______________________________________________________________________________
1146 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1148 // propagate seed to given x applying material correction if requested
1149 const Double_t kEpsilon = 1e-5;
1150 Double_t xpos = seed->GetX();
1151 Int_t dir = (xpos<xToGo) ? 1:-1;
1152 Double_t xyz0[3],xyz1[3];
1154 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1155 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1156 while ( (xToGo-xpos)*dir > kEpsilon){
1157 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1158 Double_t x = xpos+step;
1159 Double_t bz=GetBz(); // getting the local Bz
1160 if (!seed->PropagateToX(x,bz)) return kFALSE;
1162 if (matCorr || updTime) {
1163 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1166 seed->GetXYZ(xyz1); // // global pos at the end of step
1169 ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);
1170 if (dir>0) xrho = -xrho; // outward should be negative
1171 if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
1173 else { // matCorr is not requested but time integral is
1174 double d0 = xyz1[0]-xyz0[0];
1175 double d1 = xyz1[1]-xyz0[1];
1176 double d2 = xyz1[2]-xyz0[2];
1177 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1180 if (updTime) seed->AddTimeStep(ds);
1181 xpos = seed->GetX();
1186 //______________________________________________________________________________
1187 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1189 // propagate seed to given x applying material correction if requested
1190 const Double_t kEpsilon = 1e-5;
1191 Double_t xpos = seed->GetX();
1192 Int_t dir = (xpos<xToGo) ? 1:-1;
1193 Double_t xyz0[3],xyz1[3];
1196 if (AliDebugLevelClass()>1) {
1197 AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
1198 seed->AliExternalTrackParam::Print();
1201 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1202 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1203 while ( (xToGo-xpos)*dir > kEpsilon){
1204 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1205 Double_t x = xpos+step;
1206 Double_t bz=GetBz(); // getting the local Bz
1207 if (!seed->PropagateTo(x,bz)) return kFALSE;
1209 if (matCorr || updTime) {
1210 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1213 seed->GetXYZ(xyz1); // // global pos at the end of step
1217 ds = GetMaterialBudget(xyz0,xyz1,xx0,xrho);
1218 if (dir>0) xrho = -xrho; // outward should be negative
1219 if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1221 else { // matCorr is not requested but time integral is
1222 double d0 = xyz1[0]-xyz0[0];
1223 double d1 = xyz1[1]-xyz0[1];
1224 double d2 = xyz1[2]-xyz0[2];
1225 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1228 if (updTime) seed->AddTimeStep(ds);
1230 xpos = seed->GetX();
1233 if (AliDebugLevelClass()>1) {
1234 AliDebug(2,Form("After Propagate to X=%f",xToGo));
1235 seed->AliExternalTrackParam::Print();
1241 //______________________________________________________________________________
1242 Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
1244 // finalize single TPC track hypothesis
1245 if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
1246 AliITSUSeed* winner = 0;
1248 fCurrMass = hyp->GetMass();
1249 if (!(winner=hyp->DefineWinner(fCurrTrackCond->GetActiveLrInner()))) return kFALSE;
1252 #ifdef _ITSU_TUNING_MODE_ // fill tuning histos
1253 TObjArray* dest = hyp->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1255 AliITSUSeed* sd = winner;
1256 double htrPt = hyp->Pt();
1259 if ( (clID=sd->GetLrCluster(lrID))<0 ) continue;
1260 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Nrm ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2GloNrm());
1261 ((TH2*)dest->At( GetHistoID(lrID,kHBestInBranch,fCurrPassID) ))->Fill(htrPt,sd->GetOrdBranch());
1262 ((TH2*)dest->At( GetHistoID(lrID,kHBestInCand ,fCurrPassID) ))->Fill(htrPt,sd->GetOrdCand());
1264 if (dest==fCHistoArrFake && !sd->IsFake()) continue; // for the fake seeds fill only fake clusters part
1266 ((TH2*)dest->At( GetHistoID(lrID,kHResY ,fCurrPassID) ))->Fill(htrPt,sd->GetResidY());
1267 ((TH2*)dest->At( GetHistoID(lrID,kHResZ ,fCurrPassID) ))->Fill(htrPt,sd->GetResidZ());
1268 ((TH2*)dest->At( GetHistoID(lrID,kHResYP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullY());
1269 ((TH2*)dest->At( GetHistoID(lrID,kHResZP ,fCurrPassID) ))->Fill(htrPt,sd->GetPullZ());
1270 ((TH2*)dest->At( GetHistoID(lrID,kHChi2Cl ,fCurrPassID) ))->Fill(htrPt,sd->GetChi2Cl());
1272 } while((sd=(AliITSUSeed*)sd->GetParent()));
1274 ((TH2*)dest->At( GetHistoID(-1,kHChiMatch,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSTPC());
1275 ((TH2*)dest->At( GetHistoID(-1,kHChiITSSA,fCurrPassID) ))->Fill(htrPt,winner->GetChi2ITSSA());
1280 CheckClusterSharingConflicts(hyp);
1281 return hyp->GetWinner(); // winner might change of disappear after resolving conflicts
1285 //______________________________________________________________________________
1286 void AliITSUTrackerGlo::FinalizeHypotheses()
1288 // select winner for each hypothesis, remove cl. sharing conflicts
1289 AliCodeTimerAuto("",0);
1291 int* esdTrackIndex = fESDIndex.GetArray();
1292 for (int ih=0;ih<fNTracksESD;ih++) FinalizeHypothesis(GetTrackHyp(esdTrackIndex[ih])); // finlize in decreasing pt order
1294 for (int ih=fNTracksESD;ih--;) UpdateESDTrack(GetTrackHyp(esdTrackIndex[ih]),AliESDtrack::kITSin);
1298 //______________________________________________________________________________
1299 void AliITSUTrackerGlo::CheckClusterSharingConflicts(AliITSUTrackHyp* hyp)
1301 // remove eventual cluster sharing conflicts
1302 AliITSUSeed* winner = 0;
1303 if (!(winner=hyp->GetWinner())) return;
1304 UShort_t idH = (UShort_t)hyp->GetUniqueID()+1;
1306 int inLr = fCurrTrackCond->GetActiveLrInner();
1307 AliITSUSeed* winSD=winner;
1309 if ( (clID=winSD->GetLrCluster(lrID))<0 ) continue;
1310 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1311 Int_t refID = cl->GetRecoInfo(); // was it already referred by some track (with refID-1 if >0)
1312 if (!refID) {cl->SetRecoInfo(idH); continue;} // if not, refer from track to cluster
1313 if (refID==idH) continue; // ignore reference to itself
1314 // the cluster is already used by other track, need to resolve the conflict
1315 AliITSUTrackHyp* hypC = GetTrackHyp(refID-1); // competitor
1316 // printf("Ref to %d (%p) from %d Cl %d %d\n",refID-1,hypC, idH-1,clID,lrID);
1317 AliITSUSeed* winnerC = hypC->GetWinner();
1319 printf("Missing winner, ERROR\n");
1324 if (AliDebugLevelClass()>1) { //------------------------ DEBUG -------------------
1325 AliInfo(Form("Shared cl#%4d on Lr:%d: Hyp%3d/Q:%6.3f vs Hyp%3d/Q:%6.3f",
1326 clID,lrID,idH-1,winner->GetQualityVar(),refID-1,winnerC->GetQualityVar()));
1327 // dump winner of hyp
1328 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1329 idH-1,winner->Pt(),winner->GetChi2GloNrm(),winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),
1330 winner->GetNLayersHit(),hyp->GetITSLabel(),hyp->GetESDTrack()->GetTPCLabel());
1332 AliITSUSeed* sd = winner;
1335 int clIDt = sd->GetLrCluster(lrs);
1336 if (clIDt<0) continue;
1337 while( lrs>++prevL ) printf("%4s ","----");
1338 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1339 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1341 // dump winner of hypC
1342 printf("#%4d Pt:%.3f Chi:%6.3f/%6.3f/%6.3f Ncl:%d MCits%+5d MCtpc:%+5d |",
1343 refID-1,winnerC->Pt(),winnerC->GetChi2GloNrm(),winnerC->GetChi2ITSTPC(),winnerC->GetChi2ITSSA(),
1344 winnerC->GetNLayersHit(),hypC->GetITSLabel(),hypC->GetESDTrack()->GetTPCLabel());
1349 int clIDt = sd->GetLrCluster(lrs);
1350 if (clIDt<0) continue;
1351 while( lrs>++prevL ) printf("%4s ","----");
1352 printf("%4d (%5.1f)",clIDt,sd->GetChi2Cl());
1353 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1355 } //-------------------------------------------------- DEBUG -------------------
1358 if (winner->GetQualityVar()<winnerC->GetQualityVar()) { // current tracks is better then competitor track
1359 FlagSeedClusters(winnerC,kFALSE,refID); // unflag cluster usage by loser
1360 cl->SetRecoInfo(idH);
1362 //if ( hypC->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",refID-1));
1364 if (hypC->DefineWinner(inLr)) { // find new winner instead of suppressed candidate
1366 CheckClusterSharingConflicts(hypC); // and check its sharing conflicts again
1367 if (winner->IsKilled()) break; // the current winner might have been killed during check of new winner of hypC!
1370 else { // competitor hypC is better than the hyp
1371 FlagSeedClusters(winner,kFALSE,idH); // unflag cluster usage by loser
1373 //if ( hyp->GetLabel()>=0) AliInfo(Form("KILLING CORRECT: %d",idH-1));
1374 if (hyp->DefineWinner(inLr)) {
1376 CheckClusterSharingConflicts(hyp);
1381 } while ((winSD=(AliITSUSeed*)winSD->GetParent()));
1386 //______________________________________________________________________________
1387 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1389 // update ESD track with current best hypothesis
1391 AliESDtrack* esdTr = hyp->GetESDTrack();
1392 if (!esdTr || esdTr->IsOn(flag)) return;
1393 AliITSUSeed* win = hyp->GetWinner();
1394 if (!win || win->IsKilled()) return;
1398 case AliESDtrack::kITSin:
1399 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1401 // set fakes cluster info
1403 UShort_t clfake = 0;
1405 AliITSUSeed* sd = win;
1408 int lr, clID = sd->GetLrCluster(lr);
1409 if (sd->IsFake()) clfake |= 0x1<<lr;
1411 esdTr->SetITSModuleIndex(ip++, sd->GetLrClusterID());
1412 AliITSUClusterPix *cl = (AliITSUClusterPix*)fITS->GetLayerActive(lr)->GetCluster(clID);
1413 #ifdef _ITSU_TUNING_MODE_
1414 if (cl->IsSplit()) clSplit |= 0x1<<lr;
1417 } while ((sd=(AliITSUSeed*)sd->GetParent()));
1419 // RS: Temporary set special flag for tracks from the afterburner
1420 if (fCurrPassID>0) clfake |= 0x1<<7;
1422 esdTr->SetITSSharedMap(clfake);
1423 esdTr->SetITSModuleIndex(10,clSplit);
1425 // TEMPORARY: store iteration id
1426 esdTr->SetITSModuleIndex(11,fCurrPassID);
1429 case AliESDtrack::kITSout:
1430 // here the stored chi2 will correspond to backward ITS-SA fit
1431 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1433 // TODO: avoid setting friend
1436 case AliESDtrack::kITSrefit:
1438 // at the moment fill as a chi2 the TPC-ITS matching chi2
1439 chiSav = hyp->GetChi2();
1440 hyp->SetChi2(win->GetChi2ITSTPC());
1441 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1442 hyp->SetChi2(chiSav);
1444 // TODO: avoid setting cluster info
1447 AliFatal(Form("Unknown flag %d",flag));
1450 esdTr->SetITSLabel(hyp->GetITSLabel());
1451 // transfer chip indices
1455 //______________________________________________________________________________
1456 Double_t AliITSUTrackerGlo::RefitTrack(AliITSUTrackHyp* trc, Double_t rDest, Int_t &nclFit, Int_t stopCond)
1458 // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
1459 // if stopCond<0 : propagate till last cluster then stop
1460 // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
1461 // if stopCond>0 : rDest must be reached
1463 double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1464 int dir,lrStart,lrStop;
1466 dir = rCurr<rDest ? 1 : -1;
1467 lrStart = fITS->FindFirstLayerID(rCurr,dir);
1468 lrStop = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
1471 if (AliDebugLevelClass()>2) {
1472 printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1473 printf("Before refit: "); trc->Print();
1476 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));
1478 Int_t clInfo[2*AliITSUAux::kMaxLayers];
1479 Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
1480 fCurrMass = fCurrHyp->GetMass();
1481 AliITSUTrackHyp tmpTr(*(AliKalmanTrack*)trc);
1486 int lrStop1 = lrStop+dir;
1487 for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
1488 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1489 if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1490 int ilrA2,ilrA = lr->GetActiveID();
1491 // passive layer or active w/o hits will be traversed on the way to next cluster
1492 if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue;
1494 // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1496 if (dir>0) { // clusters are stored in increasing radius order
1497 iclLr[nclLr++]=clInfo[ilrA2++];
1498 if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
1501 if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
1502 iclLr[nclLr++]=clInfo[ilrA2];
1505 Bool_t transportedToLayer = kFALSE;
1506 for (int icl=0;icl<nclLr;icl++) {
1507 AliITSUClusterPix* clus = (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1508 AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1509 if (!tmpTr.Rotate(sens->GetPhiTF())) {
1511 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1512 AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
1517 double xClus = sens->GetXTF()+clus->GetX();
1518 if (!transportedToLayer) {
1519 if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1521 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1522 AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
1524 return -1; // go to the entrance to the layer
1527 transportedToLayer = kTRUE;
1531 if (AliDebugLevelClass()>1) {
1532 AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1536 if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1538 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1539 AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
1544 Double_t p[2]={clus->GetY(), clus->GetZ()};
1545 Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
1546 double chi2cl = tmpTr.GetPredictedChi2(p,cov);
1549 #ifdef _ITSU_TUNING_MODE_
1550 TObjArray* dest = trc->GetLabel()>=0 ? fCHistoArrCorr : fCHistoArrFake;
1551 if (dest && fTrackPhaseID>kClus2Tracks) {
1553 double htrPt = tmpTr.Pt();
1554 double dy = p[0]-tmpTr.GetY();
1555 double dz = p[1]-tmpTr.GetZ();
1556 ((TH2*)dest->At( GetHistoID(ilrA,kHResY,0,fTrackPhaseID) ))->Fill(htrPt,dy);
1557 ((TH2*)dest->At( GetHistoID(ilrA,kHResZ,0,fTrackPhaseID) ))->Fill(htrPt,dz);
1558 double errY = tmpTr.GetSigmaY2();
1559 double errZ = tmpTr.GetSigmaZ2();
1560 if (errY>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResYP,0,fTrackPhaseID) ))->Fill(htrPt,dy/Sqrt(errY));
1561 if (errZ>0) ((TH2*)dest->At( GetHistoID(ilrA,kHResZP,0,fTrackPhaseID) ))->Fill(htrPt,dz/Sqrt(errZ));
1562 ((TH2*)dest->At( GetHistoID(ilrA,kHChi2Cl,0,fTrackPhaseID) ))->Fill(htrPt,chi2cl);
1566 if ( !tmpTr.Update(p,cov) ) {
1568 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1569 AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
1574 if (AliDebugLevelClass()>1) {
1575 printf("AfterRefit: "); tmpTr.Print();
1578 if (++nclFit==nCl && stopCond<0) {
1579 *trc = (AliKalmanTrack&)tmpTr;
1580 return chi2; // it was requested to not propagate after last update
1585 // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1586 // Still, try to go as close as possible to rDest.
1588 if (lrStart!=lrStop) {
1589 if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1591 AliDebug(1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1593 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1595 if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1597 AliDebug(1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1599 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1602 // go to the destination radius. Note that here we don't select direction to avoid precision problems
1603 if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
1604 return (stopCond>0) ? -chi2 : chi2; // rDest was obligatory
1606 *trc = (AliKalmanTrack&)tmpTr;
1608 if (AliDebugLevelClass()>2) {
1609 printf("After refit (now at lr %d): ",lrStop); trc->Print();
1615 //__________________________________________________________________
1616 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp)
1620 const int kMaxLbPerCl = 3;
1621 int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1622 Int_t lr,nLab=0,nCl=0;
1623 AliITSUSeed *seed = hyp->GetWinner();
1625 int clID = seed->GetLrCluster(lr);
1627 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1629 for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1630 int trLb = cl->GetLabel(imc);
1632 // search this mc track in already accounted ones
1634 for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1635 if (iLab<nLab) lbStat[iLab]++;
1640 } // loop over given cluster's labels
1642 seed = (AliITSUSeed*)seed->GetParent();
1643 } // loop over clusters
1645 AliESDtrack* esdTr = hyp->GetESDTrack();
1646 int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1648 int maxLab=0,nTPCok=0;
1649 for (int ilb=nLab;ilb--;) {
1650 int st = lbStat[ilb];
1651 if (lbStat[maxLab]<st) maxLab=ilb;
1652 if (tpcLab==lbID[ilb]) nTPCok += st;
1654 hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1655 hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1656 hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
1660 hyp->SetFakeRatio(-1.);
1661 hyp->SetLabel( -tpcLab );
1662 hyp->SetITSLabel( kDummyLabel );
1665 //__________________________________________________________________
1666 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1668 // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1669 if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
1670 AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
1671 fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
1672 AliITSUSeed** tmpArr = fLayerCandidates;
1673 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
1674 memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
1675 delete[] tmpArr; // delete only array, not objects
1677 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1678 int slot=fNBranchesAdded++;
1679 for (int slotF=slot;slotF--;) { // slotF is always slot-1
1680 AliITSUSeed* si = branches[slotF];
1681 if (si->Compare(seed)<0) break; // found the last seed with better quality
1682 // otherwise, shift the worse seed to the next slot
1683 branches[slot] = si;
1684 slot = slotF; // slot should be slotF+1
1686 // if needed, move worse seeds
1687 branches[slot] = seed;
1692 //__________________________________________________________________
1693 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1695 // keep allowed number of branches for current seed and disable extras
1696 AliCodeTimerAuto("",0);
1697 int nb = Min(fNBranchesAdded,acceptMax);
1698 // if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1699 // disable unused branches
1700 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1702 #ifdef _ITSU_TUNING_MODE_
1703 for (int ib=0;ib<fNBranchesAdded;ib++) branches[ib]->SetOrdBranch(ib);
1706 for (int ib=nb;ib<fNBranchesAdded;ib++) MarkSeedFree(branches[ib]);
1707 fNCandidatesAdded += nb; // update total candidates counter
1708 fNBranchesAdded = 0; // reset branches counter
1712 //__________________________________________________________________
1713 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1715 // transfer allowed number of branches to hypothesis container
1717 AliCodeTimerAuto("",0);
1718 // sort candidates in increasing order of chi2
1719 static int lastSize = 0;
1720 static int *index = 0;
1721 static float *chi2 = 0;
1722 if (fLayerMaxCandidates>lastSize) {
1723 lastSize = fLayerMaxCandidates;
1726 index = new int[lastSize];
1727 chi2 = new float[lastSize];
1729 for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1730 Sort(fNCandidatesAdded,chi2,index,kFALSE);
1733 if (ilr>fCurrTrackCond->GetActiveLrInner()) { // just take 1st acceptMax candidates
1734 nb = Min(fNCandidatesAdded,acceptMax);
1735 for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1737 else { // for innermost layer test ITS SA backward chi2 and matching to TPC
1738 AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
1739 for (nb=0;nb<fNCandidatesAdded;nb++) {
1740 AliITSUSeed* sd = fLayerCandidates[index[nb]];
1741 fCurrHyp->SetWinner(sd);
1742 if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
1744 AddProlongationHypothesis(sd,ilr);
1745 if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
1748 fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
1750 // disable unused candidates
1751 for (int ib=nb;ib<fNCandidatesAdded;ib++) MarkSeedFree(fLayerCandidates[index[ib]]);
1752 fNCandidatesAdded = 0; // reset candidates counter
1756 //__________________________________________________________________
1757 Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
1759 // check seed backward propagation chi2 and matching to TPC
1760 double bz0 = GetBz();
1761 double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
1762 AliITSUTrackHyp trback;
1763 trback.AliExternalTrackParam::operator=(*seed);
1764 trback.ResetCovariance(10000);
1766 double chi2sa = RefitTrack(&trback,rDest,nclFit,1);
1767 if (chi2sa<0) return kFALSE;
1768 int ndf = nclFit*2-5;
1769 if (ndf>0) chi2sa /= ndf;
1770 if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2(nclFit)) return kFALSE;
1772 // relate to TPC track at outer layer
1773 AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
1774 if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
1776 if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
1777 AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
1784 double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
1785 if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
1787 seed->SetChi2ITSSA(chi2sa);
1788 seed->SetChi2ITSTPC(chi2Match);
1792 //__________________________________________________________________
1793 void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
1795 // remove those hypothesis seeds which dont lead to candidates at final layer
1797 // 1st, flag the seeds to save
1798 int lr0 = fCurrTrackCond->GetActiveLrInner();
1799 for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
1800 AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
1801 if (!seed) continue;
1802 seed->FlagTree(AliITSUSeed::kSave);
1803 dest->AddSeed(seed,lr0);
1805 for (int ilr=0;ilr<fNLrActive;ilr++) {
1806 int nsd = fCurrHyp->GetNSeeds(ilr);
1807 for (int isd=0;isd<nsd;isd++) {
1808 AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
1809 if (!seed) continue; // already discarded or saved
1810 if (seed->IsSaved()) dest->AddSeed(seed,ilr);
1811 else MarkSeedFree(seed);
1815 // AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
1818 //__________________________________________________________________
1819 void AliITSUTrackerGlo::CleanHypothesis(AliITSUTrackHyp* hyp)
1823 for (int ilr=0;ilr<fNLrActive;ilr++) {
1824 int nsd = hyp->GetNSeeds(ilr);
1825 for (int isd=0;isd<nsd;isd++) {
1826 AliITSUSeed* seed = hyp->RemoveSeed(ilr,isd);
1827 if (!seed) continue; // already discarded or saved
1833 //__________________________________________________________________
1834 void AliITSUTrackerGlo::FlagSplitClusters()
1836 // set special bit on split clusters using MC info
1837 for (int ilr=fNLrActive;ilr--;) {
1839 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1840 for (int isn=lr->GetNSensors();isn--;) {
1841 AliITSURecoSens* sens = lr->GetSensor(isn);
1842 int nCl = sens->GetNClusters();
1844 int cl0 = sens->GetFirstClusterId();
1845 for (int icl=nCl;icl--;) {
1846 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1847 for (int icl1=icl;icl1--;) {
1848 AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1849 if (cl->HasCommonTrack(cl1)) {
1850 if (!cl->IsSplit()) nsplit++;
1851 if (!cl1->IsSplit()) nsplit++;
1858 AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1863 //__________________________________________________________________
1864 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1866 // check if the seed contains split cluster with size < maxSize
1868 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1869 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1870 if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1872 return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1876 //__________________________________________________________________
1877 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1879 // print seeds clusters
1881 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1882 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1885 if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1889 //__________________________________________________________________
1890 void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg, UShort_t ref)
1892 // mark used clusters
1895 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1896 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1897 if (ref) { // do we need to set or delete cluster->track ref?
1899 if (!cl->GetRecoInfo()) cl->SetRecoInfo(ref); // set ref only if cluster already does not refer to other track, inc.counter
1902 if (cl->GetRecoInfo()==ref) cl->SetRecoInfo(0); // unset reference only if it refers to ref, decrease counter
1906 seed = (AliITSUSeed*)seed->GetParent();
1911 //__________________________________________________________________
1912 void AliITSUTrackerGlo::CheckClusterUsage()
1914 // check cluster usage info
1915 printf("ClusterUsage at pass %d\n",fCurrPassID);
1916 for (int ilr=0;ilr<fNLrActive;ilr++) {
1917 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1918 int ncl = lr->GetNClusters();
1920 int nclUs0CntShr1 = 0;
1924 for (int icl=0;icl<ncl;icl++) {
1925 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(icl);
1926 int usc = cl->GetClUsage();
1928 if (cl->IsClusterUsed() || cl->IsClusterShared()) {
1930 printf("Lr%d Cl%4d: UseCnt=0 but UseFlg=%d UseShr=%d\n",ilr,icl,cl->IsClusterUsed(),cl->IsClusterShared());
1935 if (!cl->IsClusterUsed()) {
1937 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1939 if (cl->IsClusterShared()) {
1941 printf("Lr%d Cl%4d: UseCnt=%d but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1945 printf("Lr%d Cl%4d: UseCnt=%d!! but UseFlg=%d UseShr=%d\n",ilr,icl,usc, cl->IsClusterUsed(),cl->IsClusterShared());
1948 printf("ClUsage Lr%d: Ncl=%5d Nusd=%5d (%5.2f) Use0CS1:%4d Use1C0:%4d Use1S1:%4d UseL:%d\n",
1949 ilr,ncl,nclUsed,ncl? float(nclUsed)/ncl : 0, nclUs0CntShr1,nclUs1Cnt0,nclUs1Shr1,nclUseL);
1954 #ifdef _ITSU_TUNING_MODE_
1955 //__________________________________________________________________
1956 TObjArray* AliITSUTrackerGlo::BookControlHistos(const char* pref)
1958 // book special control histos
1959 TString prefS = pref;
1961 TObjArray* dest = new TObjArray;
1962 dest->SetOwner(kTRUE);
1964 const int kNResDef=7;
1965 const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1966 const double ptMax=10;
1967 const double plMax=10;
1968 const double chiMax=100;
1969 const int nptbins=50;
1970 const int nresbins=400;
1971 const int nplbins=50;
1972 const int nchbins=200;
1973 const int maxBr = 15;
1974 const int maxCand = 200;
1976 for (int phase=0;phase<kNTrackingPhases;phase++) {
1977 for (int pass=0;pass<AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();pass++) {
1979 for (int ilr=0;ilr<fNLrActive;ilr++) {
1981 // ----------------- These are histos to be filled in Cluster2Tracks of each pass.
1982 // PropagateBack and RefitInward will be stored among the histos of 1st pass
1983 if (pass>0 && phase!=kClus2Tracks) continue;
1985 double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1986 ttl = Form("Pass%d_S%d_residY%d_%s",pass,phase,ilr,pref);
1987 TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1988 dest->AddAtAndExpand(hdy,GetHistoID(ilr,kHResY,pass,phase));
1989 hdy->SetDirectory(0);
1991 ttl = Form("Pass%d_S%d_residYPull%d_%s",pass,phase,ilr,pref);
1992 TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1993 dest->AddAtAndExpand(hdyp,GetHistoID(ilr,kHResYP,pass,phase));
1994 hdyp->SetDirectory(0);
1996 ttl = Form("Pass%d_S%d_residZ%d_%s",pass,phase,ilr,pref);
1997 TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1998 dest->AddAtAndExpand(hdz,GetHistoID(ilr,kHResZ,pass,phase));
1999 hdz->SetDirectory(0);
2001 ttl = Form("Pass%d_S%d_residZPull%d_%s",pass,phase,ilr,pref);
2002 TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
2003 hdzp->SetDirectory(0);
2004 dest->AddAtAndExpand(hdzp,GetHistoID(ilr,kHResZP,pass,phase));
2006 ttl = Form("Pass%d_S%d_chi2Cl%d_%s",pass,phase,ilr,pref);
2007 TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2008 hchi->SetDirectory(0);
2009 dest->AddAtAndExpand(hchi,GetHistoID(ilr,kHChi2Cl,pass,phase));
2011 // ------------------- These histos are filled for Clusters2Tracks only
2012 if (phase!=kClus2Tracks) continue;
2014 ttl = Form("Pass%d_S%d_chi2Nrm%d_%s",pass,phase,ilr,pref);
2015 TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2016 hchiN->SetDirectory(0);
2017 dest->AddAtAndExpand(hchiN,GetHistoID(ilr,kHChi2Nrm,pass,phase));
2019 ttl = Form("Pass%d_S%d_bestInBranch%d_%s",pass,phase,ilr,pref);
2020 TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
2021 hnbr->SetDirectory(0);
2022 dest->AddAtAndExpand(hnbr,GetHistoID(ilr,kHBestInBranch,pass,phase));
2024 ttl = Form("Pass%d_S%d_bestInCands%d_%s",pass,phase,ilr,pref);
2025 TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
2026 hncn->SetDirectory(0);
2027 dest->AddAtAndExpand(hncn,GetHistoID(ilr,kHBestInCand,pass,phase));
2029 } // loop over layers
2031 // these are special histos, filled not per layer but in the end of track fit in Clusters2Tracks in EVERY pass
2033 if (phase==kClus2Tracks) {
2035 ttl = Form("Pass%d_Chi2Match_%s",pass,pref);
2036 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
2037 hchiMatch->SetDirectory(0);
2038 dest->AddAtAndExpand(hchiMatch,GetHistoID(-1,kHChiMatch,pass,phase));
2041 ttl = Form("Pass%d_Chi2ITSSA_%s",pass,pref);
2042 hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
2043 hchiSA->SetDirectory(0);
2044 dest->AddAtAndExpand(hchiSA,GetHistoID(-1,kHChiITSSA,pass,phase));
2047 } // loop over tracking passes
2048 }// loop over tracking phases
2053 //__________________________________________________________________
2054 Int_t AliITSUTrackerGlo::GetHistoID(Int_t lr, Int_t hid, Int_t pass, Int_t phase)
2056 // get id for the requested histo
2059 return pass*kHistosPass + phase*kHistosPhase + lr*kMaxHID + hid;