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 "AliITSUAux.h"
39 #include "AliITSUClusterPix.h"
40 #include "AliITSUGeomTGeo.h"
41 #include "AliCodeTimer.h"
42 #include "AliRefArray.h"
43 using namespace AliITSUAux;
44 using namespace TMath;
46 //----------------- tmp stuff -----------------
48 ClassImp(AliITSUTrackerGlo)
50 const Double_t AliITSUTrackerGlo::fgkToler = 1e-6;// tolerance for layer on-surface check
53 //_________________________________________________________________________
54 AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
58 ,fCurrESDtrMClb(kDummyLabel)
60 ,fCountProlongationTrials(0)
65 ,fLayerMaxCandidates(1000)
71 ,fSeedsPool("AliITSUSeed",0)
82 #ifdef _FILL_CONTROL_HISTOS_
86 // Default constructor
90 //_________________________________________________________________________
91 AliITSUTrackerGlo::~AliITSUTrackerGlo()
95 delete[] fLayerCandidates;
96 if (fWorkHyp) fWorkHyp->SetTPCSeed(0); // this hypothesis does not own the seed
99 #ifdef _FILL_CONTROL_HISTOS_
101 TFile* ctrOut = TFile::Open("itsuControlHistos.root","recreate");
103 AliInfo("Storing control histos");
105 // ctrOut->WriteObject(fCHistoArr,"controlH","kSingleKey");
115 //_________________________________________________________________________
116 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
118 // init with external reconstructor
120 fITS = rec->GetITSInterface();
121 fNLrActive = fITS->GetNLayersActive();
122 fWorkHyp = new AliITSUTrackHyp(fNLrActive);
124 if (fLayerMaxCandidates<1) fLayerMaxCandidates = 1000;
125 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
126 fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
127 fFreeSeedsID.Set(1000);
133 //_________________________________________________________________________
134 void AliITSUTrackerGlo::CreateDefaultTrackCond()
136 // creates default tracking conditions to be used when recoparam does not provide them
138 AliITSUTrackCond* cond = new AliITSUTrackCond();
140 cond->SetNLayers(fNLrActive);
141 cond->AddNewCondition(fNLrActive);
142 cond->AddGroupPattern( 0xffff ); // require all layers hit
145 fDefTrackConds.AddLast(cond);
147 AliInfo("Created conditions: ");
148 for (int i=0;i<fDefTrackConds.GetEntriesFast();i++) fDefTrackConds[i]->Print();
153 //_________________________________________________________________________
154 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
157 AliCodeTimerAuto("",0);
158 SetTrackingPhase(kClus2Tracks);
160 #ifdef _FILL_CONTROL_HISTOS_
161 if (!fCHistoArr) BookControlHistos();
164 static TArrayF esdTrPt(fESDIndex.GetSize());
166 TObjArray *trackConds = 0;
168 fCountProlongationTrials = 0;
174 int nTrESD = esdEv->GetNumberOfTracks();
175 AliInfo(Form("Will try to find prolongations for %d tracks",nTrESD));
176 int nTrackCond = AliITSUReconstructor::GetRecoParam()->GetNTrackingConditions();
178 if (!fDefTrackConds.GetEntriesFast()) {
179 AliInfo("No tracking conditions found in recoparams, creating default one requesting all layers hit");
180 CreateDefaultTrackCond();
182 trackConds = &fDefTrackConds;
183 nTrackCond = trackConds->GetEntriesFast();
185 else trackConds = AliITSUReconstructor::GetRecoParam()->GetTrackingConditions();
187 static Bool_t first = kTRUE;
189 AliITSUReconstructor::GetRecoParam()->Print();
193 if (fHypStore.GetSize()<nTrESD) fHypStore.Expand(nTrESD+100);
195 fITS->ProcessClusters();
197 #ifdef _FILL_CONTROL_HISTOS_
198 FlagSplitClusters(); // tmp RS
201 // the tracks will be reconstructed in decreasing pt order, sort them
202 if (fESDIndex.GetSize()<nTrESD) {
203 fESDIndex.Set(nTrESD+200);
204 esdTrPt.Set(nTrESD+200);
206 int* esdTrackIndex = fESDIndex.GetArray();
207 float* trPt = esdTrPt.GetArray();
208 for (int itr=nTrESD;itr--;) trPt[itr] = esdEv->GetTrack(itr)->Pt();
209 Sort(nTrESD,trPt,esdTrackIndex,kTRUE);
211 for (int icnd=0;icnd<nTrackCond;icnd++) {
212 fCurrTrackCond = (AliITSUTrackCond*)trackConds->UncheckedAt(icnd);
213 if (!fCurrTrackCond->IsInitDone()) fCurrTrackCond->Init();
214 // select ESD tracks to propagate
215 for (int itr=0;itr<nTrESD;itr++) {
216 int trID = esdTrackIndex[itr];
217 fCurrESDtrack = esdEv->GetTrack(trID);
218 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
220 if (!NeedToProlong(fCurrESDtrack)) continue; // are we interested in this track?
222 // if specific tracks should be checked, create a global array int wtc[] = {evITS*10000+trID, ...};
224 int nwtc = sizeof(wtc)/sizeof(int);
226 for (int k=0;k<nwtc;k++) {
227 if (wtc[k]==evID*10000+trID) {
229 printf("\n\n\n\n\n\n\n");
230 printf("At esdTr: %d MC: %d\n",wtc[k],fCurrESDtrMClb);
231 printf("\n\n\n\n\n\n\n");
235 AliLog::SetClassDebugLevel("AliITSUTrackerGlo",dbg ? 10:0);
238 AliDebug(1,Form("Processing track %d(esd%d) | M=%.3f Pt=%.3f | MCLabel: %d",itr,trID,fCurrESDtrack->GetMass(kTRUE),fCurrESDtrack->Pt(),fCurrESDtrMClb));//RS
239 FindTrack(fCurrESDtrack, trID);
242 if (AliDebugLevelClass()>+2) {
243 AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
246 FinalizeHypotheses();
248 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0); // in case wtc array was defined, uncomment this
251 AliInfo(Form("%d ITSin for %d tried TPC seeds out of %d ESD tracks\n",fCountITSin,fCountProlongationTrials,nTrESD));
256 //_________________________________________________________________________
257 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
260 // Do outward fits in ITS
261 AliCodeTimerAuto("",0);
263 SetTrackingPhase(kPropBack);
264 int nTrESD = esdEv->GetNumberOfTracks();
265 AliDebug(1,Form("Will propagate back %d tracks",nTrESD));
267 double bz0 = GetBz();
268 Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
269 AliITSUTrackHyp dummyTr;
270 const double kWatchStep=10.; // for larger steps watch arc vs segment difference
271 Double_t times[AliPID::kSPECIES];
273 for (int itr=0;itr<nTrESD;itr++) {
274 fCurrESDtrack = esdEv->GetTrack(itr);
275 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
276 // Start time integral and add distance from current position to vertex
277 if (fCurrESDtrack->IsOn(AliESDtrack::kITSout)) continue; // already done
279 fCurrESDtrack->GetXYZ(xyzTrk);
280 Double_t dst = 0.; // set initial track lenght, tof
282 double dxs = xyzTrk[0] - xyzVtx[0];
283 double dys = xyzTrk[1] - xyzVtx[1];
284 double dzs = xyzTrk[2] - xyzVtx[2];
285 // RS: for large segment steps d use approximation of cicrular arc s by
286 // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
287 // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
288 // Hence s^2/d^2 = (1+1/6 p^2)^2
289 dst = dxs*dxs + dys*dys;
290 if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
291 double crv = Abs(fCurrESDtrack->GetC(bz0));
292 double fcarc = 1.+crv*crv*dst/6.;
299 fCurrESDtrack->SetStatus(AliESDtrack::kTIME);
301 if (!fCurrESDtrack->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
302 dummyTr.AliExternalTrackParam::operator=(*fCurrESDtrack);
303 dummyTr.StartTimeIntegral();
304 dummyTr.AddTimeStep(dst);
305 dummyTr.GetIntegratedTimes(times);
306 fCurrESDtrack->SetIntegratedTimes(times);
307 fCurrESDtrack->SetIntegratedLength(dummyTr.GetIntegratedLength());
311 fCurrHyp = GetTrackHyp(itr);
312 fCurrMass = fCurrHyp->GetMass();
313 fCurrHyp->StartTimeIntegral();
314 fCurrHyp->AddTimeStep(dst);
315 fCurrHyp->ResetCovariance(10000);
316 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMax());
317 if (chi2>0) { // propagate to exit from the ITS/TPC screen
318 int ndf = fCurrHyp->GetWinner()->GetNLayersHit()*2-5;
319 if (ndf>0) chi2 /= ndf;
320 fCurrHyp->SetChi2(chi2);
321 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSout);
325 AliDebug(2,Form("Refit Failed for track %d | ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
326 //fCurrHyp->AliExternalTrackParam::Print();
327 //fCurrHyp->GetWinner()->Print();
332 AliInfo(Form("%d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
333 fCountITSout,fCountITSin,fCountProlongationTrials,nTrESD));
338 //_________________________________________________________________________
339 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent *esdEv)
342 // refit the tracks inward, using current cov.matrix
343 AliCodeTimerAuto("",0);
345 SetTrackingPhase(kRefitInw);
346 Int_t nTrESD = esdEv->GetNumberOfTracks();
347 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
349 AliDebug(1,Form("Will refit inward %d tracks",nTrESD));
351 for (int itr=0;itr<nTrESD;itr++) {
352 fCurrESDtrack = esdEv->GetTrack(itr);
353 fCurrESDtrMClb = fCurrESDtrack->GetLabel();
354 // Start time integral and add distance from current position to vertex
355 UInt_t trStat = fCurrESDtrack->GetStatus();
356 if ( !(trStat & AliESDtrack::kITSout) ) continue;
357 if ( trStat & AliESDtrack::kITSrefit ) continue; // already done
358 if ( (trStat & AliESDtrack::kTPCout) && !(trStat & AliESDtrack::kTPCrefit) ) continue;
360 fCurrHyp = GetTrackHyp(itr);
361 fCurrHyp->AliExternalTrackParam::operator=(*fCurrESDtrack); // fetch current ESDtrack kinematics
362 fCurrMass = fCurrHyp->GetMass();
364 double chi2 = RefitTrack(fCurrHyp,fITS->GetRMin());
365 if (chi2>0) { // propagate up to inside radius of the beam pipe
366 fCurrHyp->SetChi2(chi2);
367 UpdateESDTrack(fCurrHyp,AliESDtrack::kITSrefit);
371 AliDebug(2,Form("Refit Failed for track %d |ESDtrack#%d (MClb:%d)",itr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
375 AliInfo(Form("%d ITSrefit in %d ITSout in %d ITSin tracks for %d tried TPC seeds out of %d ESD tracks\n",
376 fCountITSrefit,fCountITSout,fCountITSin,fCountProlongationTrials,nTrESD));
378 // AliLog::SetClassDebugLevel("AliITSUTrackerGlo",0);
382 //_________________________________________________________________________
383 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
385 // read from tree (if pointer provided) or directly from the ITS reco interface
387 return fReconstructor->LoadClusters(treeRP);
390 //_________________________________________________________________________
391 void AliITSUTrackerGlo::UnloadClusters()
397 Info("UnloadClusters","To be implemented");
399 //_________________________________________________________________________
400 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
406 Info("GetCluster","To be implemented");
410 //_________________________________________________________________________
411 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr)
413 // do we need to match this track to ITS?
415 static double bz = GetBz();
416 if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
417 esdTr->IsOn(AliESDtrack::kTPCout) ||
418 esdTr->IsOn(AliESDtrack::kITSin) ||
419 esdTr->GetKinkIndex(0)>0) return kFALSE;
421 if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
424 esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz);
425 // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
426 if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation())
427 && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
428 Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
430 // if (esdTr->Pt()<3) return kFALSE;//RS
434 //_________________________________________________________________________
435 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
437 // find prolongaion candidates finding for single seed
439 AliITSUSeed seedUC; // copy of the seed from the upper layer
440 AliITSUSeed seedT; // transient seed between the seedUC and new prolongation hypothesis
442 AliITSUTrackHyp* hypTr = InitHypothesis(esdTr,esdID); // initialize prolongations hypotheses tree
443 if (!hypTr || hypTr->GetSkip()) return;
446 fCurrHyp->InitFrom(hypTr);
448 AliITSURecoSens *hitSens[AliITSURecoSens::kNNeighbors+1];
450 TObjArray clArr; // container for transfer of clusters matching to seed
452 for (int ila=fNLrActive;ila--;) {
454 fCurrLayer = fITS->GetLayerActive(ila);
455 Bool_t noClSharing = fCurrTrackCond->GetClSharing(ila)==0;
456 int ilaUp = ila+1; // prolong seeds from layer above
458 // for the outermost layer the seed is created from the ESD track
459 int nSeedsUp = (ilaUp==fNLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
460 int maxNBranches = fCurrTrackCond->GetMaxBranches(ila);
461 int maxNCandidates = fCurrTrackCond->GetMaxCandidates(ila);
463 for (int isd=0;isd<nSeedsUp;isd++) {
465 if (ilaUp==fNLrActive) {
467 seedUC.InitFromSeed(fCurrHyp->GetTPCSeed()); // init with TPC seed from ref.R
470 seedU = fCurrHyp->GetSeed(ilaUp,isd); // seed on prev.active layer to prolong
471 if (seedU->IsKilled()) continue;
472 seedUC = *seedU; // its copy will be prolonged
473 seedUC.SetParent(seedU);
475 seedUC.ResetFMatrix(); // reset the matrix for propagation to next layer
476 // go till next active layer
477 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()));
478 if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) { // external seed already prolonged
480 AliDebug(2,Form("Transport failed | esdID=%d (MClb:%d)",esdID,fCurrESDtrMClb));
481 // Check if the seed satisfies to track definition
482 if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
483 continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
485 if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
486 if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
492 if (!seedUC.ContainsFake() && AliDebugLevelClass()>2) {
493 mcquest = fCurrESDtrMClb;
495 printf("FSParams: "); for (int ip=0;ip<kNTrImpData;ip++) printf("%+e ",fTrImpData[ip]); printf("\n");
498 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens, mcquest); // find detectors which may be hit by the track
500 int nsens = fCurrLayer->FindSensors(&fTrImpData[kTrPhi0], hitSens); // find detectors which may be hit by the track
501 AliDebug(2,Form("Will check %d sensors on lr:%d | esdID=%d (MClb:%d)",nsens,ila,esdID,fCurrESDtrMClb));
506 for (int isn=nsens;isn--;) {
507 AliITSURecoSens* sens = hitSens[isn];
508 int ncl = sens->GetNClusters();
512 // We need to propagate the seed to sensor on fCurrLayer staying in the frame of the sensor from prev.layer,
513 // since the transport matrix should be defined in this frame.
514 double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
515 if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) {
516 if (AliDebugLevelClass()>2) {
517 printf("Failed on GetTrackingXAtXAlpha: X=%.4f alp=%.3f\n",sens->GetXTF(),sens->GetPhiTF());
522 if (xs<seedT.GetX()) {
523 if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
525 else { // some low precision tracks may hit the sensor plane outside of the layer radius
526 if (AliDebugLevelClass()>2) {
527 if (!seedT.ContainsFake()) {
528 printf("WRONG PROP on L%d, sens%d of %d: %.4f > %.4f\n",ila,isn,nsens,xs,seedT.GetX());
533 if (!seedT.PropagateParamOnlyTo(xs,bz)) continue;
535 // if (!seedT.PropagateToX(xs,bz)) continue;
536 // if (!seedT.Rotate(sens->GetPhiTF())) continue;
537 if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
539 int clID0 = sens->GetFirstClusterId();
540 for (int icl=ncl;icl--;) { // don't change the order, clusters are sorted
541 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo",10);
542 int clID = clID0 + icl;
543 if (noClSharing && fCurrLayer->GetCluster(clID)->IsClusterUsed()) continue;
544 int res = CheckCluster(&seedT,ila,clID0+icl);
545 //AliLog::SetClassDebugLevel("AliITSUTrackerGlo", 0);
547 if (res==kStopSearchOnSensor) break; // stop looking on this sensor
548 if (res==kClusterNotMatching) continue; // cluster does not match
549 // cluster is matching and it was added to the hypotheses tree
552 // cluster search is done. Do we need to have a version of this seed skipping current layer
553 if (!NeedToKill(&seedUC,kMissingCluster)) {
554 AliITSUSeed* seedSkp = NewSeedFromPool(&seedUC);
555 double penalty = -fCurrTrackCond->GetMissPenalty(ila);
556 // to do: make penalty to account for probability to miss the cluster for good reason
557 seedSkp->SetChi2Cl(penalty);
558 AddSeedBranch(seedSkp);
559 // AddProlongationHypothesis(seedSkp,ila);
561 // transfer the new branches of the seed to the hypothesis container
562 if (fNBranchesAdded) ValidateAllowedBranches(maxNBranches);
565 if (fNCandidatesAdded) ValidateAllowedCandidates(ila,maxNCandidates);
566 // ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
567 if (AliDebugLevelClass()>2) { //RS
568 printf(">>> All hypotheses on lr %d: \n",ila);
569 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();}
573 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
574 if (ila!=0) continue;
575 double vecL[5] = {0};
576 double matL[15] = {0};
577 AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
578 while(sp->GetParent()) {
579 sp->Smooth(vecL,matL);
580 if (sp->GetLayerID()>=fNLrActive-1) break;
581 sp = (AliITSUSeed*)sp->GetParent();
586 SaveReducedHypothesesTree(hypTr);
587 if (AliDebugLevelClass()>1) {
588 printf("\nSaved hypotheses for esdTrack %d (MCLab:%d)\n",esdID,fCurrESDtrMClb);
595 //_________________________________________________________________________
596 AliITSUTrackHyp* AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
598 // init prolongaion candidates finding for single seed
599 AliITSUTrackHyp* hyp = GetTrackHyp(esdID);
602 fCountProlongationTrials++;
604 fCurrMass = esdTr->GetMass();
605 if (fCurrMass<kPionMass*0.9) fCurrMass = kPionMass; // don't trust to mu, e identification from TPCin
607 hyp = new AliITSUTrackHyp(fNLrActive);
608 hyp->SetESDTrack(esdTr);
609 hyp->SetUniqueID(esdID);
610 hyp->SetMass(fCurrMass);
611 hyp->SetTPCSeed( new AliExternalTrackParam(*esdTr) );
612 SetTrackHyp(hyp,esdID);
614 if (!TransportToLayer(hyp->GetTPCSeed(),fITS->GetNLayers(), fITS->GetLrIDActive(fNLrActive-1), fITS->GetRITSTPCRef())) hyp->SetSkip(); // propagate to outer R of ITS
620 //_________________________________________________________________________
621 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
623 // transport seed from layerFrom to the entrance (along the direction of the track) of layerTo or to rLim (if>0), wathever is closer
626 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
628 int dir = lTo > lFrom ? 1:-1;
629 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
630 Bool_t checkFirst = kTRUE;
631 Bool_t limReached = kFALSE;
634 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
637 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
638 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
640 // go the entrance of the layer, assuming no materials in between
641 double xToGo = lrTo->GetR(-dir);
644 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
647 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
650 // double xts = xToGo;
651 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
652 //printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
653 //seed->Print("etp");
657 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
658 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
659 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
660 //seed->Print("etp");
664 if (limReached) break;
670 //_________________________________________________________________________
671 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t rLim)
673 // transport track from layerFrom to the entrance of layerTo or to rLim (if>0), wathever is closer
675 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
677 int dir = lTo > lFrom ? 1:-1;
678 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
679 Bool_t checkFirst = kTRUE;
680 Bool_t limReached = kFALSE;
683 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
686 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
687 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
689 // go the entrance of the layer, assuming no materials in between
690 double xToGo = lrTo->GetR(-dir);
693 if (rLim<xToGo) {xToGo = rLim; limReached = kTRUE;}
696 if (rLim>xToGo) {xToGo = rLim; limReached = kTRUE;}
699 // double xts = xToGo;
700 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
701 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
702 // seed->Print("etp");
705 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
706 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
707 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
708 //seed->Print("etp");
712 if (limReached) break;
718 //_________________________________________________________________________
719 Bool_t AliITSUTrackerGlo::TransportToLayerX(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo, Double_t xStop)
721 // transport track from layerFrom to the entrance of layerTo but do not pass control parameter X
723 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
725 int dir = lTo > lFrom ? 1:-1;
726 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
727 Bool_t checkFirst = kTRUE;
730 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
733 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
734 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
736 // go the entrance of the layer, assuming no materials in between
737 double xToGo = lrTo->GetR(-dir); // R of the entrance to layer
739 // double xts = xToGo;
740 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) {
741 // printf("FailHere1: %f %f %d\n",xts,xToGo,dir);
742 // seed->Print("etp");
745 if ( (dir>0&&xToGo>xStop) || (dir<0&&xToGo<xStop) ) xToGo = xStop;
747 AliDebug(2,Form("go in dir=%d to R=%.4f(X:%.4f)",dir,lrTo->GetR(-dir), xToGo));
748 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) {
749 //printf("FailHere2: %f %f %d\n",xts,xToGo,dir);
750 //seed->Print("etp");
759 //_________________________________________________________________________
760 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
762 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
763 // If check is requested, do this only provided the track has not exited the layer already
764 double xToGo = lr->GetR(dir);
765 if (check) { // do we need to track till the surface of the current layer ?
766 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
767 AliDebug(3,Form(" dir:%d Cur: %e Tgt: %e",dir,Sqrt(curR2),xToGo));
768 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
769 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
771 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
772 // go via layer to its boundary, applying material correction.
773 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
779 //_________________________________________________________________________
780 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
782 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
783 // If check is requested, do this only provided the track has not exited the layer already
784 double xToGo = lr->GetR(dir);
785 if (check) { // do we need to track till the surface of the current layer ?
786 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
787 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // on the surface or outside of the layer
788 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // on the surface or outside of the layer
790 AliDebug(2,Form(" dir=%d : from R=%.4f -> R=%.4f",dir,Sqrt(seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY()), xToGo));
791 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
792 // go via layer to its boundary, applying material correction.
793 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
800 //_________________________________________________________________________
801 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
803 // go to the entrance of lr in direction dir, w/o applying material corrections.
804 // If check is requested, do this only provided the track did not reach the layer already
805 double xToGo = lr->GetR(-dir);
806 if (check) { // do we need to track till the surface of the current layer ?
807 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
808 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
809 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
811 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
812 // go via layer to its boundary, applying material correction.
813 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
818 //_________________________________________________________________________
819 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
821 // go to the entrance of lr in direction dir, w/o applying material corrections.
822 // If check is requested, do this only provided the track did not reach the layer already
823 double xToGo = lr->GetR(-dir);
824 if (check) { // do we need to track till the surface of the current layer ?
825 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
826 if (dir>0) { if (curR2-xToGo*xToGo>-fgkToler) return kTRUE; } // already passed
827 else if (dir<0) { if (xToGo*xToGo-curR2>-fgkToler) return kTRUE; } // already passed
829 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
830 // go via layer to its boundary, applying material correction.
831 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
836 //_________________________________________________________________________
837 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
839 // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
840 // as well as some aux info
842 AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
843 seed->GetXYZ(&fTrImpData[kTrXIn]); // lab position at the entrance from above
844 static AliExternalTrackParam sc; // seed copy for manipulations
847 fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
848 if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
849 double dr = lrA->GetDR(); // approximate X dist at the inner radius
850 if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
851 // special case: track does not reach inner radius, might be tangential
852 double r = sc.GetD(0,0,bz);
854 if (!sc.GetXatLabR(r,x,bz,-1)) return kFALSE;
855 dr = Abs(sc.GetX() - x);
856 if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) return kFALSE;
859 fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
860 double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
861 double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
862 sgy = Sqrt(sgy)*fCurrTrackCond->GetNSigmaRoadY(ilrA);
863 sgz = Sqrt(sgz)*fCurrTrackCond->GetNSigmaRoadZ(ilrA);
864 double dphi0 = 0.5*Abs(fTrImpData[kTrPhiOut]-fTrImpData[kTrPhiIn]);
865 double phi0 = 0.5*(fTrImpData[kTrPhiOut]+fTrImpData[kTrPhiIn]);
866 if ( dphi0>(0.5*Pi()) ) {
867 // special case of angles around pi
872 fTrImpData[kTrPhi0] = phi0;
873 fTrImpData[kTrZ0] = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
874 dphi0 += sgy/lrA->GetR();
875 fTrImpData[kTrDPhi] = dphi0<PiOver2() ? dphi0 : PiOver2();
876 fTrImpData[kTrDZ] = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn]) + sgz;
881 //________________________________________
882 void AliITSUTrackerGlo::ResetSeedsPool()
884 // mark all seeds in the pool as unused
885 AliDebug(-1,Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
887 fSeedsPool.Clear(); // seeds don't allocate memory
891 //________________________________________
892 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd)
894 // account that this seed is "deleted"
895 int id = sd->GetPoolID();
897 AliError(Form("Freeing of seed %p NOT from the pool is requested",sd));
900 // AliInfo(Form("%d %p",id, seed));
901 fSeedsPool.RemoveAt(id);
902 if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
903 fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
906 //_________________________________________________________________________
907 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID)
909 // Check if the cluster (in tracking frame!) is matching to track.
910 // The track must be already propagated to sensor tracking frame.
911 // Returns: kStopSearchOnSensor if the search on given sensor should be stopped,
912 // kClusterMatching if the cluster is matching
913 // kClusterMatching otherwise
915 const double kTolerX = 5e-4;
916 // AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
917 AliCluster *cl = fCurrLayer->GetCluster(clID);
919 Bool_t goodCl = kFALSE;
920 int currLabel = Abs(fCurrESDtrMClb);
922 if (cl->GetLabel(0)>=0) {
923 for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
925 Bool_t goodSeed = !track->ContainsFake();
927 if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
928 if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
929 if (goodCl&&goodSeed && AliDebugLevelClass()>2 ) {
930 AliDebug(2,Form("Lost good cl on L:%d failed propagation. |ESDtrack#%d (MClb:%d)",lr,fCurrESDtrack->GetID(),fCurrESDtrMClb));
934 return kStopSearchOnSensor; // propagation failed, seedT is intact
937 double dy = cl->GetY()-track->GetY();
938 double dz = cl->GetZ()-track->GetZ();
940 #ifdef _FILL_CONTROL_HISTOS_
941 int hcOffs = (1+fTrackPhase)*kHistosPhase + lr;
943 if (goodCl && (((AliITSUClusterPix*)cl)->GetNPix()>1 || !((AliITSUClusterPix*)cl)->IsSplit()) && goodSeed && fCHistoArr /* && track->GetChi2Penalty()<1e-5*/) {
945 ((TH2*)fCHistoArr->At(kHResY+hcOffs))->Fill(htrPt,dy);
946 ((TH2*)fCHistoArr->At(kHResZ+hcOffs))->Fill(htrPt,dz);
947 double errY = track->GetSigmaY2();
948 double errZ = track->GetSigmaZ2();
949 if (errY>0) ((TH2*)fCHistoArr->At(kHResYP+hcOffs))->Fill(htrPt,dy/Sqrt(errY));
950 if (errZ>0) ((TH2*)fCHistoArr->At(kHResZP+hcOffs))->Fill(htrPt,dz/Sqrt(errZ));
955 double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
956 fCurrTrackCond->GetNSigmaRoadY(lr)*fCurrTrackCond->GetNSigmaRoadY(lr); // RS TOOPTIMIZE
957 if (dy2>tol2) { // the clusters are sorted in Z(col) then in Y(row).
958 if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
959 AliDebug(2,Form("Lost good cl: dy2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",dy2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
962 PrintSeedClusters(track);
964 // clusters are sorted in increasing Y and the loop goes from last (highers Y) to 1st.
965 // if the track has too large y for this cluster (dy<0) it will be so for all other clusters of the sensor
966 if (dy<0) return kStopSearchOnSensor; // No chance that other cluster of this sensor will match (all Y's will be even larger)
967 else return kClusterNotMatching; // Other clusters may match
970 tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
971 fCurrTrackCond->GetNSigmaRoadZ(lr)*fCurrTrackCond->GetNSigmaRoadZ(lr); // RS TOOPTIMIZE
973 if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
974 AliDebug(2,Form("Lost good cl on L:%d : dz2=%e > tol2=%e |ESDtrack#%d (MClb:%d)",lr,dz2,tol2,fCurrESDtrack->GetID(),fCurrESDtrMClb));
977 PrintSeedClusters(track);
979 return kClusterNotMatching; // Other clusters may match
983 Double_t p[2]={cl->GetY(), cl->GetZ()};
984 Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
985 double chi2 = track->GetPredictedChi2(p,cov);
987 #ifdef _FILL_CONTROL_HISTOS_
989 ((TH2*)fCHistoArr->At(kHChi2Cl+hcOffs))->Fill(htrPt,chi2);
993 if (chi2>fCurrTrackCond->GetMaxTr2ClChi2(lr)) {
994 if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
995 AliDebug(2,Form("Lost good cl on L:%d : Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e |ESDtrack#%d (MClb:%d)\n",
996 lr,chi2,fCurrTrackCond->GetMaxTr2ClChi2(lr),dy,dz,fCurrESDtrack->GetID(),fCurrESDtrMClb));
999 PrintSeedClusters(track);
1001 return kClusterNotMatching;
1004 track = NewSeedFromPool(track); // input track will be reused, use its clone for updates
1005 if (!track->Update()) {
1006 if (goodCl&&goodSeed && AliDebugLevelClass()>2) {
1007 AliDebug(2,Form("Lost good cluster on L:%d : failed to update",lr));
1008 track->Print("etp");
1010 PrintSeedClusters(track);
1012 MarkSeedFree(track);
1013 return kClusterNotMatching;
1015 track->SetChi2Cl(chi2);
1016 track->SetLrClusterID(lr,clID);
1017 // cl->IncreaseClusterUsage(); // do this only for winners
1019 track->SetFake(!goodCl);
1021 AliDebug(2,Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)| ",
1022 goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
1024 AddSeedBranch(track);
1025 #ifdef _FILL_CONTROL_HISTOS_
1027 ((TH2*)fCHistoArr->At(kHChi2Nrm+hcOffs))->Fill(htrPt,track->GetChi2GloNrm());
1031 return kClusterMatching;
1034 //_________________________________________________________________________
1035 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
1037 // check if the seed should not be discarded
1038 const UShort_t kMask = 0xffff;
1039 if (flag==kMissingCluster) {
1040 int lastChecked = seed->GetLayerID();
1041 UShort_t patt = seed->GetHitsPattern();
1042 if (lastChecked) patt |= ~(kMask<<lastChecked); // not all layers were checked, complete unchecked ones by potential hits
1043 Bool_t seedOK = fCurrTrackCond->CheckPattern(patt);
1049 //______________________________________________________________________________
1050 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1052 // propagate seed to given x applying material correction if requested
1053 const Double_t kEpsilon = 1e-5;
1054 Double_t xpos = seed->GetX();
1055 Int_t dir = (xpos<xToGo) ? 1:-1;
1056 Double_t xyz0[3],xyz1[3],param[7];
1058 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1059 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1060 while ( (xToGo-xpos)*dir > kEpsilon){
1061 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1062 Double_t x = xpos+step;
1063 Double_t bz=GetBz(); // getting the local Bz
1064 if (!seed->PropagateToX(x,bz)) return kFALSE;
1066 if (matCorr || updTime) {
1067 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1070 seed->GetXYZ(xyz1); // // global pos at the end of step
1072 MeanMaterialBudget(xyz0,xyz1,param);
1073 Double_t xrho=param[0]*param[4], xx0=param[1];
1074 if (dir>0) xrho = -xrho; // outward should be negative
1075 if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) {AliInfo("Fail2"); seed->Print("etp"); return kFALSE;}
1078 else { // matCorr is not requested but time integral is
1079 double d0 = xyz1[0]-xyz0[0];
1080 double d1 = xyz1[1]-xyz0[1];
1081 double d2 = xyz1[2]-xyz0[2];
1082 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1085 if (updTime) seed->AddTimeStep(ds);
1086 xpos = seed->GetX();
1091 //______________________________________________________________________________
1092 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
1094 // propagate seed to given x applying material correction if requested
1095 const Double_t kEpsilon = 1e-5;
1096 Double_t xpos = seed->GetX();
1097 Int_t dir = (xpos<xToGo) ? 1:-1;
1098 Double_t xyz0[3],xyz1[3],param[7];
1100 if (AliDebugLevelClass()>1) {
1101 AliDebug(2,Form("Before Propagate to X=%f with M=%.3f MaxStep=%.4f MatCorr=%d",xToGo,mass,maxStep,matCorr));
1102 seed->AliExternalTrackParam::Print();
1104 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
1105 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
1106 while ( (xToGo-xpos)*dir > kEpsilon){
1107 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1108 Double_t x = xpos+step;
1109 Double_t bz=GetBz(); // getting the local Bz
1110 if (!seed->PropagateTo(x,bz)) return kFALSE;
1112 if (matCorr || updTime) {
1113 xyz0[0]=xyz1[0]; // global pos at the beginning of step
1116 seed->GetXYZ(xyz1); // // global pos at the end of step
1119 MeanMaterialBudget(xyz0,xyz1,param);
1120 Double_t xrho=param[0]*param[4], xx0=param[1];
1121 if (dir>0) xrho = -xrho; // outward should be negative
1122 if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
1125 else { // matCorr is not requested but time integral is
1126 double d0 = xyz1[0]-xyz0[0];
1127 double d1 = xyz1[1]-xyz0[1];
1128 double d2 = xyz1[2]-xyz0[2];
1129 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
1132 if (updTime) seed->AddTimeStep(ds);
1134 xpos = seed->GetX();
1136 if (AliDebugLevelClass()>1) {
1137 AliDebug(2,Form("After Propagate to X=%f",xToGo));
1138 seed->AliExternalTrackParam::Print();
1143 //______________________________________________________________________________
1144 Bool_t AliITSUTrackerGlo::FinalizeHypothesis(AliITSUTrackHyp* hyp)
1146 // finalize single TPC track hypothesis
1147 if (!hyp || hyp->GetESDTrack()->IsOn(AliESDtrack::kITSin)) return kFALSE;
1148 AliITSUSeed* winner = 0;
1150 fCurrMass = hyp->GetMass();
1151 if (!(winner=hyp->DefineWinner())) return kFALSE;
1152 FlagSeedClusters(winner,kTRUE);
1155 #ifdef _FILL_CONTROL_HISTOS_
1159 int lb = Abs(hyp->GetESDTrack()->GetTPCLabel());
1160 if ( hyp->GetITSLabel()==lb ) {
1161 hchimt = (TH2*)fCHistoArr->At(kHChiMatchCorr);
1162 hchiSA = (TH2*)fCHistoArr->At(kHChiITSSACorr);
1164 else if ( Abs(hyp->GetITSLabel()) == lb ) hchimt = (TH2*)fCHistoArr->At(kHChiMatchFake);
1165 else if ( hyp->GetITSLabel()>=0 ) hchimt = (TH2*)fCHistoArr->At(kHChiMatchCorrMiss);
1166 else hchimt = (TH2*)fCHistoArr->At(kHChiMatchFakeMiss);
1167 if (!hchiSA) hchiSA = (TH2*)fCHistoArr->At(kHChiITSSAFake);
1168 // printf("MTStatus: ITS:%+5d TPC:%+5d chimt:%e chi2SA:%e-> %s\n",hyp->GetITSLabel(),lb,winner->GetChi2ITSTPC(),winner->GetChi2ITSSA(),hchimt->GetName());
1169 if (hchimt) hchimt->Fill(hyp->GetTPCSeed()->Pt(),winner->GetChi2ITSTPC());
1170 if (hchiSA) hchiSA->Fill(hyp->GetTPCSeed()->Pt(),winner->GetChi2ITSSA());
1174 UpdateESDTrack(hyp,AliESDtrack::kITSin);
1178 //______________________________________________________________________________
1179 void AliITSUTrackerGlo::FinalizeHypotheses()
1181 // select winner for each hypothesis, remove cl. sharing conflicts
1182 AliCodeTimerAuto("",0);
1185 int nh = fHypStore.GetEntriesFast();
1186 for (int ih=0;ih<nh;ih++) {
1187 if (FinalizeHypothesis( GetTrackHyp(ih) )) fCountITSin++;
1190 AliITSUSeed* winner = 0;
1192 #ifdef _FILL_CONTROL_HISTOS_
1193 // if requested, collect cluster sharing statistics
1195 if (fCHistoArr && (hShare=(TH2*)fCHistoArr->At(kHClShare))) {
1196 for (int ih=0;ih<nh;ih++) {
1197 AliITSUTrackHyp* hyp = GetTrackHyp(ih);
1198 if (!hyp || !(winner=hyp->GetWinner())) continue;
1200 double pt = hyp->Pt();
1202 int clID = winner->GetLrCluster(lrID);
1203 if (clID<0) continue;
1204 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1205 if (!cl->IsClusterShared()) continue;
1206 hShare->Fill(pt,winner->IsFake() ? lrID+fNLrActive : lrID);
1207 } while ((winner=(AliITSUSeed*)winner->GetParent()));
1212 if (AliDebugLevelClass()>+2) {
1214 AliRefArray** refArr = new AliRefArray*[fNLrActive];
1215 for (int ilr=0;ilr<fNLrActive;ilr++) refArr[ilr] = new AliRefArray(1000);
1216 for (int ih=0;ih<nh;ih++) {
1217 AliITSUTrackHyp* hyp = GetTrackHyp(ih);
1218 if (!hyp || !(winner=hyp->GetWinner())) continue;
1221 int clID = winner->GetLrCluster(lrID);
1222 if (clID<0) continue;
1223 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1224 if (!cl->IsClusterShared()) continue;
1225 refArr[lrID]->AddReference(clID,ih);
1226 } while ((winner=(AliITSUSeed*)winner->GetParent()));
1229 for (int ilr=0;ilr<fNLrActive;ilr++) {
1230 int ncl = fITS->GetLayerActive(ilr)->GetNClusters();
1231 printf("\nClusterSharingDump: Lr %d (%d cl)\n",ilr,ncl);
1233 for (int icl=0;icl<ncl;icl++) {
1234 if (!refArr[ilr]->HasReference(icl)) continue;
1235 int nref = refArr[ilr]->GetReferences(icl,refs,100);
1236 printf("--- cl%3d(#%d): NShare=%4d\n",cnt++,icl,nref);
1237 for (int ir=0;ir<nref;ir++) {
1238 AliITSUTrackHyp* hyp = GetTrackHyp(refs[ir]);
1239 winner = hyp->GetWinner();
1240 AliESDtrack* esdTr = hyp->GetESDTrack();
1241 printf("#%4d Pt:%.3f Chi:%6.2f Ncl:%d MCits%+5d MCtpc:%+5d ESD:%4d |",
1242 refs[ir],winner->Pt(),winner->GetChi2GloNrm(),winner->GetNLayersHit(),
1243 hyp->GetITSLabel(),esdTr->GetTPCLabel(),esdTr->GetID());
1247 int clID = winner->GetLrCluster(lrs);
1248 if (clID<0) continue;
1249 while( lrs>++prevL ) printf("%4s ","----");
1250 printf("%4d (%5.1f)",clID,winner->GetChi2Cl());
1251 } while ((winner=(AliITSUSeed*)winner->GetParent()));
1262 //______________________________________________________________________________
1263 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp,Int_t flag)
1265 // update ESD track with current best hypothesis
1266 AliESDtrack* esdTr = hyp->GetESDTrack();
1268 AliITSUSeed* win = hyp->GetWinner();
1273 case AliESDtrack::kITSin:
1274 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1275 // TODO: set cluster info
1278 case AliESDtrack::kITSout:
1279 // here the stored chi2 will correspond to backward ITS-SA fit
1280 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1281 // TODO: avoid setting friend
1284 case AliESDtrack::kITSrefit:
1286 // at the moment fill as a chi2 the TPC-ITS matching chi2
1287 chiSav = hyp->GetChi2();
1288 hyp->SetChi2(win->GetChi2ITSTPC());
1289 esdTr->UpdateTrackParams(hyp,flag); // update kinematics
1290 hyp->SetChi2(chiSav);
1291 // TODO: avoid setting cluster info
1294 AliFatal(Form("Unknown flag %d",flag));
1297 esdTr->SetITSLabel(hyp->GetITSLabel());
1298 // transfer module indices
1302 //______________________________________________________________________________
1303 Double_t AliITSUTrackerGlo::RefitTrack(AliExternalTrackParam* trc, Double_t rDest, Int_t stopCond)
1305 // refit track till radius rDest. The cluster,mass info is taken from fCurrHyp (and its winner seed)
1306 // if stopCond<0 : propagate till last cluster then stop
1307 // if stopCond==0: propagate till last cluster then try to go till limiting rDest, don't mind if fail
1308 // if stopCond>0 : rDest must be reached
1310 double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
1311 int dir,lrStart,lrStop;
1313 dir = rCurr<rDest ? 1 : -1;
1314 lrStart = fITS->FindFirstLayerID(rCurr,dir);
1315 lrStop = fITS->FindLastLayerID(rDest,dir); // lr id before which we have to stop
1317 if (AliDebugLevelClass()>2) {
1318 printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
1319 printf("Before refit: "); trc->Print();
1321 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));
1323 Int_t clInfo[2*AliITSUAux::kMaxLayers];
1324 Int_t nCl = fCurrHyp->FetchClusterInfo(clInfo);
1325 fCurrMass = fCurrHyp->GetMass();
1326 AliExternalTrackParam tmpTr(*trc);
1328 int iclLr[2],nclLr,clCount=0;
1330 int lrStop1 = lrStop+1;
1331 for (int ilr=lrStart;ilr!=lrStop1;ilr+=dir) {
1332 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
1333 if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
1334 int ilrA2,ilrA = lr->GetActiveID();
1335 // passive layer or active w/o hits will be traversed on the way to next cluster
1336 if (!lr->IsActive() || clInfo[ilrA2=(ilrA<<1)]<0) continue;
1338 // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
1340 if (dir>0) { // clusters are stored in increasing radius order
1341 iclLr[nclLr++]=clInfo[ilrA2++];
1342 if (clInfo[ilrA2]>=0) iclLr[nclLr++]=clInfo[ilrA2];
1345 if ( clInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=clInfo[ilrA2+1];
1346 iclLr[nclLr++]=clInfo[ilrA2];
1349 Bool_t transportedToLayer = kFALSE;
1350 for (int icl=0;icl<nclLr;icl++) {
1351 AliITSUClusterPix* clus = (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
1352 AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
1353 if (!tmpTr.Rotate(sens->GetPhiTF())) {
1354 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1355 AliDebug(2,Form("Failed on rotate to %f | ESDtrack#%d (MClb:%d)",sens->GetPhiTF(),trESD->GetID(),trESD->GetTPCLabel()));
1359 double xClus = sens->GetXTF()+clus->GetX();
1360 if (!transportedToLayer) {
1361 if (ilr!=lrStart && !TransportToLayerX(&tmpTr,lrStart,ilr,xClus)) {
1362 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1363 AliDebug(2,Form("Failed to transport %d -> %d | ESDtrack#%d (MClb:%d)\n",lrStart,ilr,trESD->GetID(),trESD->GetTPCLabel()));
1364 return -1; // go to the entrance to the layer
1367 transportedToLayer = kTRUE;
1370 if (AliDebugLevelClass()>1) {
1371 AliDebug(2,Form("Propagate to cl:%d on lr %d Need to go %.4f -> %.4f",icl,ilrA,tmpTr.GetX(),xClus));
1374 if (!PropagateSeed(&tmpTr,xClus,fCurrMass)) {
1375 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1376 AliDebug(2,Form("Failed on propagate to %f (dir=%d) | ESDtrack#%d (MClb:%d)",xClus,dir,trESD->GetID(),trESD->GetTPCLabel()));
1380 Double_t p[2]={clus->GetY(), clus->GetZ()};
1381 Double_t cov[3]={clus->GetSigmaY2(), clus->GetSigmaYZ(), clus->GetSigmaZ2()};
1382 double chi2cl = tmpTr.GetPredictedChi2(p,cov);
1385 #ifdef _FILL_CONTROL_HISTOS_
1386 int hcOffs = (1+fTrackPhase)*kHistosPhase + ilrA;
1388 if (fCHistoArr && fTrackPhase>kClus2Tracks && trc->GetLabel()>=0/* && trc->Charge()>0*/) {
1390 double dy = p[0]-tmpTr.GetY();
1391 double dz = p[1]-tmpTr.GetZ();
1392 ((TH2*)fCHistoArr->At(kHResY+hcOffs))->Fill(htrPt,dy);
1393 ((TH2*)fCHistoArr->At(kHResZ+hcOffs))->Fill(htrPt,dz);
1394 double errY = tmpTr.GetSigmaY2();
1395 double errZ = tmpTr.GetSigmaZ2();
1396 if (errY>0) ((TH2*)fCHistoArr->At(kHResYP+hcOffs))->Fill(htrPt,dy/Sqrt(errY));
1397 if (errZ>0) ((TH2*)fCHistoArr->At(kHResZP+hcOffs))->Fill(htrPt,dz/Sqrt(errZ));
1398 ((TH2*)fCHistoArr->At(kHChi2Cl+hcOffs))->Fill(htrPt,chi2cl);
1402 if ( !tmpTr.Update(p,cov) ) {
1403 AliESDtrack* trESD = fCurrHyp->GetESDTrack();
1404 AliDebug(2,Form("Failed on Update | ESDtrack#%d (MClb:%d)",trESD->GetID(),trESD->GetTPCLabel()));
1407 if (AliDebugLevelClass()>1) {
1408 printf("AfterRefit: "); tmpTr.Print();
1410 if (++clCount==nCl && stopCond<0) {
1412 return chi2; // it was requested to not propagate after last update
1417 // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
1418 // Still, try to go as close as possible to rDest.
1420 if (lrStart!=lrStop) {
1421 if (!TransportToLayer(&tmpTr,lrStart,lrStop)) {
1422 AliDebug(-1,Form("Failed on TransportToLayer %d->%d | ESDtrack#%d (MClb:%d), X=%f",lrStart,lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1423 if (stopCond>0) return -chi2; // rDest was obligatory
1425 if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) {
1426 AliDebug(-1,Form("Failed on GoToExitFromLayer %d | ESDtrack#%d (MClb:%d), X=%f",lrStop,fCurrESDtrack->GetID(),fCurrESDtrMClb,tmpTr.GetX()));
1427 if (stopCond>0) return -chi2; // rDest was obligatory
1430 // go to the destination radius. Note that here we don't select direction to avoid precision problems
1431 if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),0) || !PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) {
1432 if (stopCond>0) return -chi2; // rDest was obligatory
1435 if (AliDebugLevelClass()>2) {
1436 printf("After refit (now at lr %d): ",lrStop); trc->Print();
1441 //__________________________________________________________________
1442 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp)
1446 const int kMaxLbPerCl = 3;
1447 int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
1448 Int_t lr,nLab=0,nCl=0;
1449 AliITSUSeed *seed = hyp->GetWinner();
1451 int clID = seed->GetLrCluster(lr);
1453 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
1455 for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
1456 int trLb = cl->GetLabel(imc);
1458 // search this mc track in already accounted ones
1460 for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
1461 if (iLab<nLab) lbStat[iLab]++;
1466 } // loop over given cluster's labels
1468 seed = (AliITSUSeed*)seed->GetParent();
1469 } // loop over clusters
1471 AliESDtrack* esdTr = hyp->GetESDTrack();
1472 int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
1474 int maxLab=0,nTPCok=0;
1475 for (int ilb=nLab;ilb--;) {
1476 int st = lbStat[ilb];
1477 if (lbStat[maxLab]<st) maxLab=ilb;
1478 if (tpcLab==lbID[ilb]) nTPCok += st;
1480 hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
1481 hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
1482 hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbID[maxLab] : -lbID[maxLab]); // winner label
1486 hyp->SetFakeRatio(-1.);
1487 hyp->SetLabel( -tpcLab );
1488 hyp->SetITSLabel( kDummyLabel );
1491 //__________________________________________________________________
1492 Bool_t AliITSUTrackerGlo::AddSeedBranch(AliITSUSeed* seed)
1494 // add new prolongation branch for the seed to temporary array. The seeds are sorted according to Compare f-n
1495 if (fNCandidatesAdded+fNBranchesAdded>=fLayerMaxCandidates ) {
1496 AliInfo(Form("Number of candidates at layer %d for current seed reached %d, increasing buffer",fCurrActLrID,fLayerMaxCandidates));
1497 fLayerMaxCandidates = 2*(fLayerMaxCandidates+1);
1498 AliITSUSeed** tmpArr = fLayerCandidates;
1499 fLayerCandidates = new AliITSUSeed*[fLayerMaxCandidates];
1500 memcpy(fLayerCandidates,tmpArr,(fNCandidatesAdded+fNBranchesAdded)*sizeof(AliITSUSeed*));
1501 delete tmpArr; // delete only array, not objects
1503 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded]; // note: fNCandidatesAdded is incremented after adding all branches of current seed
1504 int slot=fNBranchesAdded++;
1505 for (int slotF=slot;slotF--;) { // slotF is always slot-1
1506 AliITSUSeed* si = branches[slotF];
1507 if (si->Compare(seed)<0) break; // found the last seed with better quality
1508 // otherwise, shift the worse seed to the next slot
1509 branches[slot] = si;
1510 slot = slotF; // slot should be slotF+1
1512 // if needed, move worse seeds
1513 branches[slot] = seed;
1518 //__________________________________________________________________
1519 void AliITSUTrackerGlo::ValidateAllowedBranches(Int_t acceptMax)
1521 // keep allowed number of branches for current seed and disable extras
1522 AliCodeTimerAuto("",0);
1523 int nb = Min(fNBranchesAdded,acceptMax);
1524 // if (nb<fNBranchesAdded) printf("ValidateAllowedBranches: %d of %d (%d) on lr %d\n",nb,fNBranchesAdded,fNCandidatesAdded,ilr);
1525 // disable unused branches
1526 AliITSUSeed** branches = &fLayerCandidates[fNCandidatesAdded];
1527 #ifdef _FILL_CONTROL_HISTOS_
1529 for (int ib=0;ib<fNBranchesAdded;ib++) {
1530 AliITSUSeed* sd = branches[ib];
1531 if (!sd->ContainsFake() && (bestID<0 || sd->Compare(branches[bestID])<0) ) bestID = ib;
1534 TH2* hb = (TH2*)fCHistoArr->At(kHBestInBranch + (1+fTrackPhase)*kHistosPhase + fCurrActLrID);
1535 if (hb) hb->Fill(branches[bestID]->Pt(), bestID);
1539 for (int ib=nb;ib<fNBranchesAdded;ib++) {
1541 #ifdef _FILL_CONTROL_HISTOS_
1542 if (AliDebugLevelClass()>-2 && !branches[ib]->ContainsFake() /*&& branches[ib]->GetNLayersHit()*/
1543 && (bestID<0 || branches[ib]->Compare(branches[bestID])<0 ) ) {
1544 printf("Suppress good branch as %d of %d |",ib,fNBranchesAdded); branches[ib]->Print();
1545 // printf("Survivors : \n");
1546 // for (int j=0;j<nb;j++) branches[j]->Print();
1549 MarkSeedFree(branches[ib]);
1551 fNCandidatesAdded += nb; // update total candidates counter
1552 fNBranchesAdded = 0; // reset branches counter
1556 //__________________________________________________________________
1557 void AliITSUTrackerGlo::ValidateAllowedCandidates(Int_t ilr, Int_t acceptMax)
1559 // transfer allowed number of branches to hypothesis container
1561 AliCodeTimerAuto("",0);
1562 // sort candidates in increasing order of chi2
1563 static int lastSize = 0;
1564 static int *index = 0;
1565 static float *chi2 = 0;
1566 if (fLayerMaxCandidates>lastSize) {
1567 lastSize = fLayerMaxCandidates;
1570 index = new int[lastSize];
1571 chi2 = new float[lastSize];
1573 for (int i=fNCandidatesAdded;i--;) chi2[i] = fLayerCandidates[i]->GetChi2GloNrm();
1574 Sort(fNCandidatesAdded,chi2,index,kFALSE);
1576 #ifdef _FILL_CONTROL_HISTOS_
1578 for (int ib=0;ib<fNCandidatesAdded;ib++) {
1579 AliITSUSeed* sd = fLayerCandidates[index[ib]];
1580 if (!sd->ContainsFake() && (bestID<0 || sd->Compare(fLayerCandidates[index[bestID]])<0) ) bestID = ib;
1583 TH2* hb = (TH2*)fCHistoArr->At(kHBestInCand + (1+fTrackPhase)*kHistosPhase + fCurrActLrID);
1584 if (hb) hb->Fill(fLayerCandidates[index[bestID]]->Pt(), bestID);
1589 if (ilr>0) { // just take 1st acceptMax candidates
1590 nb = Min(fNCandidatesAdded,acceptMax);
1591 for (int ib=0;ib<nb;ib++) AddProlongationHypothesis(fLayerCandidates[index[ib]],ilr);
1593 else { // for innermost layer test ITS SA backward chi2 and matching to TPC
1594 AliITSUSeed* wn0 = fCurrHyp->GetWinner(); // in principle, should be NULL
1595 for (nb=0;nb<fNCandidatesAdded;nb++) {
1596 AliITSUSeed* sd = fLayerCandidates[index[nb]];
1597 fCurrHyp->SetWinner(sd);
1598 if (!CheckBackwardMatching(sd)) MarkSeedFree(sd); // discard bad backward tracks
1600 AddProlongationHypothesis(sd,ilr);
1601 if (++nacc==acceptMax) {nb++; break;} // the rest will be discarded
1604 fCurrHyp->SetWinner(wn0); // restore original winner (NULL?)
1607 // disable unused candidates
1608 for (int ib=nb;ib<fNCandidatesAdded;ib++) {
1610 #ifdef _FILL_CONTROL_HISTOS_
1611 if (AliDebugLevelClass()>-2 && !fLayerCandidates[index[ib]]->ContainsFake() /*&& fLayerCandidates[index[ib]]->GetNLayersHit()*/
1612 && (bestID<0 || fLayerCandidates[index[ib]]->Compare(fLayerCandidates[index[bestID]])<0 ) ) {
1613 printf("Suppress good candidate as %d of %d |",index[ib],fNCandidatesAdded); fLayerCandidates[index[ib]]->Print();
1616 MarkSeedFree(fLayerCandidates[index[ib]]);
1618 fNCandidatesAdded = 0; // reset candidates counter
1622 //__________________________________________________________________
1623 Bool_t AliITSUTrackerGlo::CheckBackwardMatching(AliITSUSeed* seed)
1625 // check seed backward propagation chi2 and matching to TPC
1626 double bz0 = GetBz();
1627 double rDest = fITS->GetRITSTPCRef(); // reference radius for comparison
1628 static AliExternalTrackParam trback;
1630 trback.ResetCovariance(10000);
1631 int ndf = seed->GetNLayersHit()*2-5;
1632 double chi2sa = RefitTrack(&trback,rDest,1);
1633 if (chi2sa<0) return kFALSE;
1634 if (ndf>0) chi2sa /= ndf;
1635 if (chi2sa>fCurrTrackCond->GetMaxITSSAChi2()) return kFALSE;
1637 // relate to TPC track at outer layer
1638 AliExternalTrackParam* tpcSeed = fCurrHyp->GetTPCSeed();
1639 if (!trback.Rotate(tpcSeed->GetAlpha()) || !trback.PropagateParamOnlyTo(tpcSeed->GetX(),bz0)) {
1640 if (AliDebugLevelClass()>+1 && !seed->ContainsFake()) {
1641 AliInfo(Form("Failed to propagate seed to TPC track @ X:%.3f Alpha:%.3f",tpcSeed->GetX(),tpcSeed->GetAlpha()));
1647 double chi2Match = trback.GetPredictedChi2(tpcSeed)/5;
1648 if (chi2Match>fCurrTrackCond->GetMaxITSTPCMatchChi2()) return kFALSE;
1650 seed->SetChi2ITSSA(chi2sa);
1651 seed->SetChi2ITSTPC(chi2Match);
1655 //__________________________________________________________________
1656 void AliITSUTrackerGlo::SaveReducedHypothesesTree(AliITSUTrackHyp* dest)
1658 // remove those hypothesis seeds which dont lead to candidates at final layer
1660 // 1st, flag the seeds to save
1662 for (int isd=0;isd<fCurrHyp->GetNSeeds(lr0);isd++) {
1663 AliITSUSeed* seed = fCurrHyp->RemoveSeed(lr0,isd);
1664 if (!seed) continue;
1665 seed->FlagTree(AliITSUSeed::kSave);
1666 dest->AddSeed(seed,lr0);
1668 for (int ilr=1;ilr<fNLrActive;ilr++) {
1669 int nsd = fCurrHyp->GetNSeeds(ilr);
1670 for (int isd=0;isd<nsd;isd++) {
1671 AliITSUSeed* seed = fCurrHyp->RemoveSeed(ilr,isd);
1672 if (!seed) continue; // already discarded or saved
1673 if (seed->IsSaved()) dest->AddSeed(seed,ilr);
1674 else MarkSeedFree(seed);
1678 // AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
1681 //__________________________________________________________________
1682 void AliITSUTrackerGlo::FlagSplitClusters()
1684 // set special bit on split clusters using MC info
1685 for (int ilr=fNLrActive;ilr--;) {
1687 AliITSURecoLayer* lr = fITS->GetLayerActive(ilr);
1688 for (int isn=lr->GetNSensors();isn--;) {
1689 AliITSURecoSens* sens = lr->GetSensor(isn);
1690 int nCl = sens->GetNClusters();
1692 int cl0 = sens->GetFirstClusterId();
1693 for (int icl=nCl;icl--;) {
1694 AliITSUClusterPix *cl = (AliITSUClusterPix*)lr->GetCluster(cl0+icl);
1695 for (int icl1=icl;icl1--;) {
1696 AliITSUClusterPix *cl1 = (AliITSUClusterPix*)lr->GetCluster(cl0+icl1);
1697 if (cl->HasCommonTrack(cl1)) {
1698 if (!cl->IsSplit()) nsplit++;
1699 if (!cl1->IsSplit()) nsplit++;
1706 AliInfo(Form("%4d out of %4d clusters are split",nsplit,lr->GetNClusters()));
1711 //__________________________________________________________________
1712 Bool_t AliITSUTrackerGlo::ContainsSplitCluster(const AliITSUSeed* seed, Int_t maxSize)
1714 // check if the seed contains split cluster with size < maxSize
1716 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1717 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1718 if (cl->IsSplit() && cl->GetNPix()<maxSize ) return kTRUE;
1720 return seed->GetParent() ? ContainsSplitCluster((AliITSUSeed*)seed->GetParent(),maxSize) : kFALSE;
1724 //__________________________________________________________________
1725 void AliITSUTrackerGlo::PrintSeedClusters(const AliITSUSeed* seed, Option_t* option)
1727 // print seeds clusters
1729 if ( (clID=seed->GetLrCluster(lrID))>=0 ) {
1730 AliITSUClusterPix* cl = (AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID);
1733 if (seed->GetParent()) PrintSeedClusters((AliITSUSeed*)seed->GetParent(), option);
1737 //__________________________________________________________________
1738 void AliITSUTrackerGlo::FlagSeedClusters(const AliITSUSeed* seed, Bool_t flg)
1740 // mark used clusters
1743 if ( (clID=seed->GetLrCluster(lrID))>=0 ) ((AliITSUClusterPix*)fITS->GetLayerActive(lrID)->GetCluster(clID))->ModClUsage(flg);
1744 seed = (AliITSUSeed*)seed->GetParent();
1751 #ifdef _FILL_CONTROL_HISTOS_
1752 //__________________________________________________________________
1753 void AliITSUTrackerGlo::BookControlHistos()
1755 // book special control histos
1756 if (fCHistoArr) return;
1757 const int kNResDef=7;
1758 const double kResDef[kNResDef]={0.05,0.05,0.3, 0.05,1,0.5,1.5};
1759 fCHistoArr = new TObjArray();
1760 fCHistoArr->SetOwner(kTRUE);
1761 const double ptMax=10;
1762 const double plMax=10;
1763 const double chiMax=100;
1764 const int nptbins=50;
1765 const int nresbins=400;
1766 const int nplbins=50;
1767 const int nchbins=200;
1768 const int maxBr = 15;
1769 const int maxCand = 200;
1771 for (int stp=0;stp<kNTrackingPhases;stp++) {
1772 for (int ilr=0;ilr<fNLrActive;ilr++) {
1773 int hoffs = (1+stp)*kHistosPhase + ilr;
1774 double mxdf = ilr>=kNResDef ? kResDef[kNResDef-1] : kResDef[ilr];
1775 ttl = Form("S%d_residY%d",stp,ilr);
1776 TH2F* hdy = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1777 fCHistoArr->AddAtAndExpand(hdy,hoffs + kHResY);
1778 hdy->SetDirectory(0);
1780 ttl = Form("S%d_residYPull%d",stp,ilr);
1781 TH2F* hdyp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1782 fCHistoArr->AddAtAndExpand(hdyp,hoffs + kHResYP);
1783 hdyp->SetDirectory(0);
1785 ttl = Form("S%d_residZ%d",stp,ilr);
1786 TH2F* hdz = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nresbins,-mxdf,mxdf);
1787 fCHistoArr->AddAtAndExpand(hdz,hoffs + kHResZ);
1788 hdz->SetDirectory(0);
1790 ttl = Form("S%d_residZPull%d",stp,ilr);
1791 TH2F* hdzp = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax,nplbins,-plMax,plMax);
1792 hdzp->SetDirectory(0);
1793 fCHistoArr->AddAtAndExpand(hdzp,hoffs + kHResZP);
1795 ttl = Form("S%d_chi2Cl%d",stp,ilr);
1796 TH2F* hchi = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1797 hchi->SetDirectory(0);
1798 fCHistoArr->AddAtAndExpand(hchi,hoffs + kHChi2Cl);
1800 ttl = Form("S%d_chi2Nrm%d",stp,ilr);
1801 TH2F* hchiN = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1802 hchiN->SetDirectory(0);
1803 fCHistoArr->AddAtAndExpand(hchiN,hoffs + kHChi2Nrm);
1805 if (stp==0) { // these histos make sense only for clusters2tracks stage
1806 ttl = Form("S%d_bestInBranch%d",stp,ilr);
1807 TH2* hnbr = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxBr,-0.5,maxBr-0.5);
1808 hnbr->SetDirectory(0);
1809 fCHistoArr->AddAtAndExpand(hnbr,hoffs + kHBestInBranch);
1811 ttl = Form("S%d_bestInCands%d",stp,ilr);
1812 TH2* hncn = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, maxCand,-0.5,maxCand-0.5);
1813 hncn->SetDirectory(0);
1814 fCHistoArr->AddAtAndExpand(hncn,hoffs + kHBestInCand);
1820 ttl = Form("ClSharing");
1821 TH2* hclShare = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, 2*fNLrActive,0,2*fNLrActive);
1822 hclShare->SetDirectory(0);
1823 fCHistoArr->AddAtAndExpand(hclShare,kHClShare);
1826 ttl = Form("Chi2MatchCorr");
1827 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1828 hchiMatch->SetDirectory(0);
1829 fCHistoArr->AddAtAndExpand(hchiMatch,kHChiMatchCorr);
1831 ttl = Form("Chi2MatchFake");
1832 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1833 hchiMatch->SetDirectory(0);
1834 fCHistoArr->AddAtAndExpand(hchiMatch,kHChiMatchFake);
1836 ttl = Form("Chi2MatchCorrMiss");
1837 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1838 hchiMatch->SetDirectory(0);
1839 fCHistoArr->AddAtAndExpand(hchiMatch,kHChiMatchCorrMiss);
1841 ttl = Form("Chi2MatchFakeMiss");
1842 hchiMatch = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins,0.,chiMax);
1843 hchiMatch->SetDirectory(0);
1844 fCHistoArr->AddAtAndExpand(hchiMatch,kHChiMatchFakeMiss);
1847 ttl = Form("Chi2ITSSACorr");
1848 hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
1849 hchiSA->SetDirectory(0);
1850 fCHistoArr->AddAtAndExpand(hchiSA,kHChiITSSACorr);
1852 ttl = Form("Chi2ITSSAFake");
1853 hchiSA = new TH2F(ttl.Data(),ttl.Data(),nptbins,0,ptMax, nchbins/2,0.,chiMax/2);
1854 hchiSA->SetDirectory(0);
1855 fCHistoArr->AddAtAndExpand(hchiSA,kHChiITSSAFake);