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>
29 #include "AliITSUTrackerGlo.h"
30 #include "AliESDEvent.h"
31 #include "AliESDtrack.h"
32 #include "AliITSURecoDet.h"
33 #include "AliITSURecoSens.h"
34 #include "AliITSUReconstructor.h"
35 #include "AliITSReconstructor.h"
36 #include "AliITSUSeed.h"
37 #include "AliITSUAux.h"
38 #include "AliITSUClusterPix.h"
39 using namespace AliITSUAux;
40 using namespace TMath;
44 //----------------- tmp stuff -----------------
46 ClassImp(AliITSUTrackerGlo)
48 const Double_t AliITSUTrackerGlo::fgkToler = 1e-6;// tolerance for layer on-surface check
51 //_________________________________________________________________________
52 AliITSUTrackerGlo::AliITSUTrackerGlo(AliITSUReconstructor* rec)
59 ,fSeedsPool("AliITSUSeed",0)
67 // Default constructor
71 //_________________________________________________________________________
72 AliITSUTrackerGlo::~AliITSUTrackerGlo()
81 //_________________________________________________________________________
82 void AliITSUTrackerGlo::Init(AliITSUReconstructor* rec)
84 // init with external reconstructor
86 fITS = new AliITSURecoDet(rec->GetGeom(),"ITSURecoInterface");
87 int nLr = fITS->GetNLayersActive();
88 fClInfo = new Int_t[nLr<<1];
89 for (int ilr=nLr;ilr--;) {
90 fITS->GetLayerActive(ilr)->SetClusters(rec->GetClusters(ilr));
93 fSeedsPool.ExpandCreateFast(1000); // RS TOCHECK
94 fFreeSeedsID.Set(1000);
96 fTrCond.SetNLayers(fITS->GetNLayersActive());
97 fTrCond.AddNewCondition(5);
98 fTrCond.AddGroupPattern( (0x1<<0)|(0x1<<1) );
99 fTrCond.AddGroupPattern( (0x1<<3)|(0x1<<4) );
100 fTrCond.AddGroupPattern( (0x1<<5)|(0x1<<6) );
102 fTrCond.AddNewCondition(5);
103 fTrCond.AddGroupPattern( (0x1<<0)|(0x1<<2) );
104 fTrCond.AddGroupPattern( (0x1<<3)|(0x1<<4) );
105 fTrCond.AddGroupPattern( (0x1<<5)|(0x1<<6) );
107 fTrCond.AddNewCondition(5);
108 fTrCond.AddGroupPattern( (0x1<<1)|(0x1<<2) );
109 fTrCond.AddGroupPattern( (0x1<<3)|(0x1<<4) );
110 fTrCond.AddGroupPattern( (0x1<<5)|(0x1<<6) );
112 printf("Tracking Conditions: ");
116 //_________________________________________________________________________
117 Int_t AliITSUTrackerGlo::Clusters2Tracks(AliESDEvent *esdEv)
122 fTrackPhase = kClus2Tracks;
123 int nTrESD = esdEv->GetNumberOfTracks();
124 AliInfo(Form("Will try to find prolongations for %d tracks",nTrESD));
125 AliITSUReconstructor::GetRecoParam()->Print();
127 if (fHypStore.GetSize()<nTrESD) fHypStore.Expand(nTrESD+100);
129 fITS->ProcessClusters();
130 // select ESD tracks to propagate
131 for (int itr=0;itr<nTrESD;itr++) {
132 AliESDtrack *esdTr = esdEv->GetTrack(itr);
133 AliInfo(Form("Processing track %d | MCLabel: %d",itr,esdTr->GetTPCLabel()));
134 if (!NeedToProlong(esdTr)) continue; // are we interested in this track?
135 FindTrack(esdTr, itr);
138 AliInfo(Form("SeedsPool: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
140 FinalizeHypotheses();
145 //_________________________________________________________________________
146 Int_t AliITSUTrackerGlo::PropagateBack(AliESDEvent *esdEv)
149 // Do outward fits in ITS
152 fTrackPhase = kPropBack;
153 int nTrESD = esdEv->GetNumberOfTracks();
154 AliInfo(Form("Will fit %d tracks",nTrESD));
156 double bz0 = GetBz();
157 Double_t xyzTrk[3],xyzVtx[3]={GetX(),GetY(),GetZ()};
158 AliITSUTrackHyp dummyTr,*currTr=0;
159 const double kWatchStep=10.; // for larger steps watch arc vs segment difference
160 Double_t times[AliPID::kSPECIES];
163 for (int itr=0;itr<nTrESD;itr++) {
164 AliESDtrack *esdTr = esdEv->GetTrack(itr);
165 // Start time integral and add distance from current position to vertex
166 if (esdTr->IsOn(AliESDtrack::kITSout)) continue; // already done
168 esdTr->GetXYZ(xyzTrk);
169 Double_t dst = 0.; // set initial track lenght, tof
171 double dxs = xyzTrk[0] - xyzVtx[0];
172 double dys = xyzTrk[1] - xyzVtx[1];
173 double dzs = xyzTrk[2] - xyzVtx[2];
174 // RS: for large segment steps d use approximation of cicrular arc s by
175 // s = 2R*asin(d/2R) = d/p asin(p) \approx d/p (p + 1/6 p^3) = d (1+1/6 p^2)
176 // where R is the track radius, p = d/2R = 2C*d (C is the abs curvature)
177 // Hence s^2/d^2 = (1+1/6 p^2)^2
178 dst = dxs*dxs + dys*dys;
179 if (dst > kWatchStep*kWatchStep) { // correct circular part for arc/segment factor
180 double crv = Abs(esdTr->GetC(bz0));
181 double fcarc = 1.+crv*crv*dst/6.;
188 esdTr->SetStatus(AliESDtrack::kTIME);
190 if (!esdTr->IsOn(AliESDtrack::kITSin)) { // case of tracks w/o ITS prolongation: just set the time integration
191 dummyTr.AliExternalTrackParam::operator=(*esdTr);
192 dummyTr.StartTimeIntegral();
193 dummyTr.AddTimeStep(dst);
194 dummyTr.GetIntegratedTimes(times);
195 esdTr->SetIntegratedTimes(times);
196 esdTr->SetIntegratedLength(dummyTr.GetIntegratedLength());
200 currTr = GetTrackHyp(itr);
201 currTr->StartTimeIntegral();
202 currTr->AddTimeStep(dst);
203 printf("Before resetCov: "); currTr->AliExternalTrackParam::Print();
204 currTr->ResetCovariance(10000);
205 RefitTrack(currTr,fITS->GetRMax()); // propagate to exit from the ITS/TPC screen
212 //_________________________________________________________________________
213 Int_t AliITSUTrackerGlo::RefitInward(AliESDEvent * /*event*/)
219 Info("RefitInward","To be implemented");
223 //_________________________________________________________________________
224 Int_t AliITSUTrackerGlo::LoadClusters(TTree * treeRP)
226 // read from tree (if pointer provided) or directly from the ITS reco interface
228 return fReconstructor->LoadClusters(treeRP);
231 //_________________________________________________________________________
232 void AliITSUTrackerGlo::UnloadClusters()
238 Info("UnloadClusters","To be implemented");
240 //_________________________________________________________________________
241 AliCluster * AliITSUTrackerGlo::GetCluster(Int_t /*index*/) const
247 Info("GetCluster","To be implemented");
251 //_________________________________________________________________________
252 Bool_t AliITSUTrackerGlo::NeedToProlong(AliESDtrack* esdTr)
254 // do we need to match this track to ITS?
256 static double bz = GetBz();
257 if (!esdTr->IsOn(AliESDtrack::kTPCin) ||
258 esdTr->IsOn(AliESDtrack::kTPCout) ||
259 esdTr->IsOn(AliESDtrack::kITSin) ||
260 esdTr->GetKinkIndex(0)>0) return kFALSE;
262 if (esdTr->Pt()<AliITSUReconstructor::GetRecoParam()->GetMinPtForProlongation()) return kFALSE;
265 esdTr->GetDZ(GetX(),GetY(),GetZ(),bz,dtz);
266 // if track is not V0 candidata but has large offset wrt IP, reject it. RS TOCHECK
267 if ( !(esdTr->GetV0Index(0)>0 && dtz[0]>AliITSUReconstructor::GetRecoParam()->GetMaxDforV0dghtrForProlongation())
268 && (Abs(dtz[0])>AliITSUReconstructor::GetRecoParam()->GetMaxDForProlongation() ||
269 Abs(dtz[1])>AliITSUReconstructor::GetRecoParam()->GetMaxDZForProlongation())) return kFALSE;
274 //_________________________________________________________________________
275 void AliITSUTrackerGlo::FindTrack(AliESDtrack* esdTr, Int_t esdID)
277 // find prolongaion candidates finding for single seed
279 AliITSUSeed seedUC; // copy of the seed from the upper layer
280 AliITSUSeed seedT; // transient seed between the seedUC and new prolongation hypothesis
282 if (!InitHypothesis(esdTr,esdID)) return; // initialize prolongations hypotheses tree
284 AliITSURecoSens *hitSens[AliITSURecoSens::kNNeighbors+1];
286 TObjArray clArr; // container for transfer of clusters matching to seed
288 int nLrActive = fITS->GetNLayersActive();
289 for (int ila=nLrActive;ila--;) {
290 int ilaUp = ila+1; // prolong seeds from layer above
292 // for the outermost layer the seed is created from the ESD track
293 int nSeedsUp = (ilaUp==nLrActive) ? 1 : fCurrHyp->GetNSeeds(ilaUp);
295 for (int isd=0;isd<nSeedsUp;isd++) {
297 if (ilaUp==nLrActive) {
299 seedUC.InitFromESDTrack(esdTr);
302 seedU = fCurrHyp->GetSeed(ilaUp,isd); // seed on prev.active layer to prolong
303 seedUC = *seedU; // its copy will be prolonged
304 seedUC.SetParent(seedU);
306 seedUC.ResetFMatrix(); // reset the matrix for propagation to next layer
307 // go till next active layer
308 AliInfo(Form("working on Lr:%d Seed:%d of %d",ila,isd,nSeedsUp));
309 if (!TransportToLayer(&seedUC, fITS->GetLrIDActive(ilaUp), fITS->GetLrIDActive(ila)) ) {
311 AliInfo("Transport failed");
312 // Check if the seed satisfies to track definition
313 if (NeedToKill(&seedUC,kTransportFailed) && seedU) KillSeed(seedU,kTRUE);
314 continue; // RS TODO: decide what to do with tracks stopped on higher layers w/o killing
316 AliITSURecoLayer* lrA = fITS->GetLayerActive(ila);
317 if (!GetRoadWidth(&seedUC, ila)) { // failed to find road width on the layer
318 if (NeedToKill(&seedUC,kRWCheckFailed) && seedU) KillSeed(seedU,kTRUE);
321 int nsens = lrA->FindSensors(&fTrImpData[kTrPhi0], hitSens); // find detectors which may be hit by the track
322 AliInfo(Form("Will check %d sensors on lr:%d ",nsens,ila));
325 for (int isn=nsens;isn--;) {
327 AliITSURecoSens* sens = hitSens[isn];
329 // We need to propagate the seed to sensor on lrA staying the frame of the sensor from prev.layer,
330 // since the transport matrix should be defined in this frame.
331 double xs; // X in the TF of current seed, corresponding to intersection with sensor plane
332 if (!seedT.GetTrackingXAtXAlpha(sens->GetXTF(),sens->GetPhiTF(),bz, xs)) continue;
333 if (!PropagateSeed(&seedT,xs,fCurrMass)) continue;
334 // if (!seedT.PropagateToX(xs,bz)) continue;
335 // if (!seedT.Rotate(sens->GetPhiTF())) continue;
336 if (!seedT.RotateToAlpha(sens->GetPhiTF())) continue;
338 int clID0 = sens->GetFirstClusterId();
339 for (int icl=sens->GetNClusters();icl--;) {
340 int res = CheckCluster(&seedT,ila,clID0+icl);
342 if (res==kStopSearchOnSensor) break; // stop looking on this sensor
343 if (res==kClusterNotMatching) continue; // cluster does not match
344 // cluster is matching and it was added to the hypotheses tree
347 // cluster search is done. Do we need ta have a version of this seed skipping current layer
349 if (!NeedToKill(&seedT,kMissingCluster)) {
350 AliITSUSeed* seedSkp = NewSeedFromPool(&seedT);
351 double penalty = -AliITSUReconstructor::GetRecoParam()->GetMissPenalty(ila);
352 // to do: make penalty to account for probability to miss the cluster for good reason
353 seedSkp->SetChi2Cl(penalty);
354 AddProlongationHypothesis(seedSkp,ila);
357 ((TObjArray*)fCurrHyp->GetLayerSeeds(ila))->Sort();
358 printf(">>> All hypotheses on lr %d: \n",ila);
359 for (int ih=0;ih<fCurrHyp->GetNSeeds(ila);ih++) {
360 printf(" #%3d ",ih); fCurrHyp->GetSeed(ila,ih)->Print();
363 if (ila!=0) continue;
364 double vecL[5] = {0};
365 double matL[15] = {0};
366 AliITSUSeed* sp = fCurrHyp->GetSeed(ila,ih);
367 while(sp->GetParent()) {
368 sp->Smooth(vecL,matL);
369 if (sp->GetLayerID()>=fITS->GetNLayersActive()-1) break;
370 sp = (AliITSUSeed*)sp->GetParent();
376 SaveCurrentTrackHypotheses();
380 //_________________________________________________________________________
381 Bool_t AliITSUTrackerGlo::InitHypothesis(AliESDtrack *esdTr, Int_t esdID)
383 // init prolongaion candidates finding for single seed
384 fCurrHyp = GetTrackHyp(esdID);
385 if (fCurrHyp) return kTRUE;
387 fCurrMass = esdTr->GetMass();
388 fCurrESDtrack = esdTr;
389 if (fCurrMass<kPionMass*0.9) fCurrMass = kPionMass; // don't trust to mu, e identification from TPCin
391 fCurrHyp = new AliITSUTrackHyp(fITS->GetNLayersActive());
392 fCurrHyp->SetESDTrack(esdTr);
393 fCurrHyp->SetUniqueID(esdID);
394 fCurrHyp->SetMass(fCurrMass);
395 SetTrackHyp(fCurrHyp,esdID);
401 //_________________________________________________________________________
402 Bool_t AliITSUTrackerGlo::TransportToLayer(AliITSUSeed* seed, Int_t lFrom, Int_t lTo)
404 // transport seed from layerFrom to the entrance of layerTo
407 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
409 int dir = lTo > lFrom ? 1:-1;
410 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
411 Bool_t checkFirst = kTRUE;
414 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
417 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
418 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
420 // go the entrance of the layer, assuming no materials in between
421 double xToGo = lrTo->GetR(-dir);
422 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
423 if (!PropagateSeed(seed,xToGo,fCurrMass,100, kFALSE )) return kFALSE;
430 //_________________________________________________________________________
431 Bool_t AliITSUTrackerGlo::TransportToLayer(AliExternalTrackParam* seed, Int_t lFrom, Int_t lTo)
433 // transport track from layerFrom to the entrance of layerTo
435 if (lTo==lFrom) AliFatal(Form("was called with lFrom=%d lTo=%d",lFrom,lTo));
437 int dir = lTo > lFrom ? 1:-1;
438 AliITSURecoLayer* lrFr = fITS->GetLayer(lFrom); // this can be 0 when extrapolation from TPC to ITS is requested
439 Bool_t checkFirst = kTRUE;
442 if (!GoToExitFromLayer(seed,lrFr,dir,checkFirst)) return kFALSE; // go till the end of current layer
445 AliITSURecoLayer* lrTo = fITS->GetLayer( (lFrom+=dir) );
446 if (!lrTo) AliFatal(Form("Layer %d does not exist",lFrom));
448 if (!GoToExitFromLayer(seed,lrTo,dir)) return kFALSE; // go the entrance of the layer, assuming no materials in between
455 //_________________________________________________________________________
456 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
458 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
459 // If check is requested, do this only provided the track has not exited the layer already
460 double xToGo = lr->GetR(dir);
461 if (check) { // do we need to track till the surface of the current layer ?
462 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
463 if (dir>0) { if (curR2-xToGo*xToGo>fgkToler) return kTRUE; } // on the surface or outside of the layer
464 else if (dir<0) { if (xToGo*xToGo-curR2>fgkToler) return kTRUE; } // on the surface or outside of the layer
466 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
467 // go via layer to its boundary, applying material correction.
468 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
473 //_________________________________________________________________________
474 Bool_t AliITSUTrackerGlo::GoToExitFromLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
476 // go to the exit from lr in direction dir, applying material corrections in steps specific for this layer
477 // If check is requested, do this only provided the track has not exited the layer already
478 double xToGo = lr->GetR(dir);
479 if (check) { // do we need to track till the surface of the current layer ?
480 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
481 if (dir>0) { if (curR2-xToGo*xToGo>fgkToler) return kTRUE; } // on the surface or outside of the layer
482 else if (dir<0) { if (xToGo*xToGo-curR2>fgkToler) return kTRUE; } // on the surface or outside of the layer
484 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
485 // go via layer to its boundary, applying material correction.
486 if (!PropagateSeed(seed,xToGo,fCurrMass, lr->GetMaxStep())) return kFALSE;
492 //_________________________________________________________________________
493 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliITSUSeed* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
495 // go to the entrance of lr in direction dir, w/o applying material corrections.
496 // If check is requested, do this only provided the track did not reach the layer already
497 double xToGo = lr->GetR(-dir);
498 if (check) { // do we need to track till the surface of the current layer ?
499 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
500 if (dir>0) { if (curR2-xToGo*xToGo>fgkToler) return kTRUE; } // already passed
501 else if (dir<0) { if (xToGo*xToGo-curR2>fgkToler) return kTRUE; } // already passed
503 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
504 // go via layer to its boundary, applying material correction.
505 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
510 //_________________________________________________________________________
511 Bool_t AliITSUTrackerGlo::GoToEntranceToLayer(AliExternalTrackParam* seed, AliITSURecoLayer* lr, Int_t dir, Bool_t check)
513 // go to the entrance of lr in direction dir, w/o applying material corrections.
514 // If check is requested, do this only provided the track did not reach the layer already
515 double xToGo = lr->GetR(-dir);
516 if (check) { // do we need to track till the surface of the current layer ?
517 double curR2 = seed->GetX()*seed->GetX() + seed->GetY()*seed->GetY(); // current radius
518 if (dir>0) { if (curR2-xToGo*xToGo>fgkToler) return kTRUE; } // already passed
519 else if (dir<0) { if (xToGo*xToGo-curR2>fgkToler) return kTRUE; } // already passed
521 if (!seed->GetXatLabR(xToGo,xToGo,GetBz(),dir)) return kFALSE;
522 // go via layer to its boundary, applying material correction.
523 if (!PropagateSeed(seed,xToGo,fCurrMass, 100, kFALSE)) return kFALSE;
528 //_________________________________________________________________________
529 Bool_t AliITSUTrackerGlo::GetRoadWidth(AliITSUSeed* seed, int ilrA)
531 // calculate road width in terms of phi and z for the track which MUST be on the external radius of the layer
532 // as well as some aux info
534 AliITSURecoLayer* lrA = fITS->GetLayerActive(ilrA);
535 seed->GetXYZ(&fTrImpData[kTrXIn]); // lab position at the entrance from above
536 static AliExternalTrackParam sc; // seed copy for manipalitions
539 fTrImpData[kTrPhiIn] = ATan2(fTrImpData[kTrYIn],fTrImpData[kTrXIn]);
540 if (!sc.Rotate(fTrImpData[kTrPhiIn])) return kFALSE; // go to the frame of the entry point into the layer
541 double dr = lrA->GetDR(); // approximate X dist at the inner radius
542 if (!sc.GetXYZAt(sc.GetX()-dr, bz, fTrImpData + kTrXOut)) {
543 // special case: track does not reach inner radius, might be tangential
544 double r = sc.GetD(0,0,bz);
546 if (!sc.GetXatLabR(r,x,bz,-1)) {
548 AliFatal(Form("This should not happen: r=%f",r));
550 dr = Abs(sc.GetX() - x);
551 if (!sc.GetXYZAt(x, bz, fTrImpData + kTrXOut)) {
553 AliFatal(Form("This should not happen: x=%f",x));
557 fTrImpData[kTrPhiOut] = ATan2(fTrImpData[kTrYOut],fTrImpData[kTrXOut]);
558 double sgy = sc.GetSigmaY2() + dr*dr*sc.GetSigmaSnp2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(ilrA);
559 double sgz = sc.GetSigmaZ2() + dr*dr*sc.GetSigmaTgl2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(ilrA);
560 sgy = Sqrt(sgy)*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY();
561 sgz = Sqrt(sgz)*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ();
562 fTrImpData[kTrPhi0] = 0.5*(fTrImpData[kTrPhiOut]+fTrImpData[kTrPhiIn]);
563 fTrImpData[kTrZ0] = 0.5*(fTrImpData[kTrZOut]+fTrImpData[kTrZIn]);
564 fTrImpData[kTrDPhi] = 0.5*Abs(fTrImpData[kTrPhiOut]-fTrImpData[kTrPhiIn]) + sgy/lrA->GetR();
565 fTrImpData[kTrDZ] = 0.5*Abs(fTrImpData[kTrZOut]-fTrImpData[kTrZIn]) + sgz;
570 //________________________________________
571 void AliITSUTrackerGlo::ResetSeedsPool()
573 // mark all seeds in the pool as unused
574 AliInfo(Form("CurrentSize: %d, BookedUpTo: %d, free: %d",fSeedsPool.GetSize(),fSeedsPool.GetEntriesFast(),fNFreeSeeds));
576 fSeedsPool.Clear(); // seeds don't allocate memory
580 //________________________________________
581 void AliITSUTrackerGlo::MarkSeedFree(AliITSUSeed *sd)
583 // account that this seed is "deleted"
584 int id = sd->GetPoolID();
586 AliError(Form("Freeing of seed %p NOT from the pool is requested",sd));
589 // AliInfo(Form("%d %p",id, seed));
590 fSeedsPool.RemoveAt(id);
591 if (fFreeSeedsID.GetSize()<=fNFreeSeeds) fFreeSeedsID.Set( 2*fNFreeSeeds + 100 );
592 fFreeSeedsID.GetArray()[fNFreeSeeds++] = id;
595 //_________________________________________________________________________
596 Int_t AliITSUTrackerGlo::CheckCluster(AliITSUSeed* track, Int_t lr, Int_t clID)
598 // Check if the cluster (in tracking frame!) is matching to track.
599 // The track must be already propagated to sensor tracking frame.
600 // Returns: kStopSearchOnSensor if the search on given sensor should be stopped,
601 // kClusterMatching if the cluster is matching
602 // kClusterMatching otherwise
604 const double kTolerX = 5e-4;
605 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
607 Bool_t goodCl = kFALSE;
608 int currLabel = Abs(fCurrESDtrack->GetTPCLabel());
610 if (cl->GetLabel(0)>=0) {
611 for (int i=0;i<3;i++) if (cl->GetLabel(i)>=0 && cl->GetLabel(i)==currLabel) {goodCl = kTRUE; break;}
615 if (TMath::Abs(cl->GetX())>kTolerX) { // if due to the misalingment X is large, propagate track only
616 if (!track->PropagateParamOnlyTo(track->GetX()+cl->GetX(),GetBz())) {
617 if (goodCl) {printf("Loose good cl: Failed propagation. |"); cl->Print();}
618 return kStopSearchOnSensor; // propagation failed, seedT is intact
621 double dy = cl->GetY()-track->GetY();
622 double dz = cl->GetZ()-track->GetZ();
625 double tol2 = (track->GetSigmaY2() + AliITSUReconstructor::GetRecoParam()->GetSigmaY2(lr))*
626 AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY()*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadY(); // RS TOOPTIMIZE
627 if (dy2>tol2) { // the clusters are sorted in Z(col) then in Y(row).
628 if (goodCl) {printf("Loose good cl: dy2=%e > tol2=%e |",dy2,tol2); cl->Print();}
629 if (dy>0) return kStopSearchOnSensor; // No chance that other cluster of this sensor will match (all Y's will be even larger)
630 else return kClusterNotMatching; // Other clusters may match
633 tol2 = (track->GetSigmaZ2() + AliITSUReconstructor::GetRecoParam()->GetSigmaZ2(lr))*
634 AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ()*AliITSUReconstructor::GetRecoParam()->GetNSigmaRoadZ(); // RS TOOPTIMIZE
636 if (goodCl) {printf("Loose good cl: dz2=%e > tol2=%e |",dz2,tol2); cl->Print();}
637 return kClusterNotMatching; // Other clusters may match
641 Double_t p[2]={cl->GetY(), cl->GetZ()};
642 Double_t cov[3]={cl->GetSigmaY2(), cl->GetSigmaYZ(), cl->GetSigmaZ2()};
643 double chi2 = track->GetPredictedChi2(p,cov);
644 if (chi2>AliITSUReconstructor::GetRecoParam()->GetMaxTr2ClChi2(lr)) {
646 printf("Loose good cl: Chi2=%e > Chi2Max=%e |dy: %+.3e dz: %+.3e\n",
647 chi2,AliITSUReconstructor::GetRecoParam()->GetMaxTr2ClChi2(lr),dy,dz);
651 return kClusterNotMatching;
654 track = NewSeedFromPool(track); // input track will be reused, use its clone for updates
655 if (!track->Update()) {
656 if (goodCl) {printf("Loose good cl: Failed update |"); cl->Print();}
658 return kClusterNotMatching;
660 track->SetChi2Cl(chi2);
661 track->SetLrClusterID(lr,clID);
662 cl->IncreaseClusterUsage();
664 track->SetFake(!goodCl);
666 AliInfo(Form("AddCl(%d) Cl%d lr:%d: dY:%+8.4f dZ:%+8.4f (MC: %5d %5d %5d) |Chi2=%f(%c)",
667 goodCl,clID,lr,dy,dz2,cl->GetLabel(0),cl->GetLabel(1),cl->GetLabel(2), chi2, track->IsFake() ? '-':'+'));
669 AddProlongationHypothesis(track,lr);
671 return kClusterMatching;
674 //_________________________________________________________________________
675 Bool_t AliITSUTrackerGlo::NeedToKill(AliITSUSeed *seed, Int_t flag)
677 // check if the seed should not be discarded
678 const UShort_t kMask = 0xffff;
679 if (flag==kMissingCluster) {
680 int lastChecked = seed->GetLayerID();
681 UShort_t patt = seed->GetHitsPattern();
682 if (lastChecked) patt |= ~(kMask<<lastChecked); // not all layers were checked, complete unchecked once by potential hits
683 Bool_t seedOK = fTrCond.CheckPattern(patt);
689 //______________________________________________________________________________
690 Bool_t AliITSUTrackerGlo::PropagateSeed(AliITSUSeed *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
692 // propagate seed to given x applying material correction if requested
693 const Double_t kEpsilon = 1e-5;
694 Double_t xpos = seed->GetX();
695 Int_t dir = (xpos<xToGo) ? 1:-1;
696 Double_t xyz0[3],xyz1[3],param[7];
698 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
699 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
700 while ( (xToGo-xpos)*dir > kEpsilon){
701 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
702 Double_t x = xpos+step;
703 Double_t bz=GetBz(); // getting the local Bz
704 if (!seed->PropagateToX(x,bz)) return kFALSE;
706 if (matCorr || updTime) {
707 xyz0[0]=xyz1[0]; // global pos at the beginning of step
710 seed->GetXYZ(xyz1); // // global pos at the end of step
712 MeanMaterialBudget(xyz0,xyz1,param);
713 Double_t xrho=param[0]*param[4], xx0=param[1];
714 if (dir>0) xrho = -xrho; // outward should be negative
715 if (!seed->ApplyMaterialCorrection(xx0,xrho,mass,kFALSE)) return kFALSE;
718 else { // matCorr is not requested but time integral is
719 double d0 = xyz1[0]-xyz0[0];
720 double d1 = xyz1[1]-xyz0[1];
721 double d2 = xyz1[2]-xyz0[2];
722 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
725 if (updTime) seed->AddTimeStep(ds);
731 //______________________________________________________________________________
732 Bool_t AliITSUTrackerGlo::PropagateSeed(AliExternalTrackParam *seed, Double_t xToGo, Double_t mass, Double_t maxStep, Bool_t matCorr)
734 // propagate seed to given x applying material correction if requested
735 const Double_t kEpsilon = 1e-5;
736 Double_t xpos = seed->GetX();
737 Int_t dir = (xpos<xToGo) ? 1:-1;
738 Double_t xyz0[3],xyz1[3],param[7];
740 Bool_t updTime = dir>0 && seed->IsStartedTimeIntegral();
741 if (matCorr || updTime) seed->GetXYZ(xyz1); //starting global position
742 while ( (xToGo-xpos)*dir > kEpsilon){
743 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
744 Double_t x = xpos+step;
745 Double_t bz=GetBz(); // getting the local Bz
746 if (!seed->PropagateTo(x,bz)) return kFALSE;
748 if (matCorr || updTime) {
749 xyz0[0]=xyz1[0]; // global pos at the beginning of step
752 seed->GetXYZ(xyz1); // // global pos at the end of step
755 MeanMaterialBudget(xyz0,xyz1,param);
756 Double_t xrho=param[0]*param[4], xx0=param[1];
757 if (dir>0) xrho = -xrho; // outward should be negative
758 if (!seed->CorrectForMeanMaterial(xx0,xrho,mass)) return kFALSE;
761 else { // matCorr is not requested but time integral is
762 double d0 = xyz1[0]-xyz0[0];
763 double d1 = xyz1[1]-xyz0[1];
764 double d2 = xyz1[2]-xyz0[2];
765 ds = TMath::Sqrt(d0*d0+d1*d1+d2*d2);
768 if (updTime) seed->AddTimeStep(ds);
775 //______________________________________________________________________________
776 void AliITSUTrackerGlo::SaveCurrentTrackHypotheses()
778 // RS: shall we clean up killed seeds?
784 //______________________________________________________________________________
785 void AliITSUTrackerGlo::FinalizeHypotheses()
787 // select winner for each hypothesis, remove cl. sharing conflicts
790 int nh = fHypStore.GetEntriesFast();
791 for (int ih=0;ih<nh;ih++) {
792 AliITSUTrackHyp* hyp = (AliITSUTrackHyp*) fHypStore.UncheckedAt(ih);
794 hyp->DefineWinner(); // TODO
800 //______________________________________________________________________________
801 void AliITSUTrackerGlo::UpdateESDTrack(AliITSUTrackHyp* hyp)
803 // update ESD track with current best hypothesis
804 AliESDtrack* esdTr = hyp->GetESDTrack();
806 AliITSUSeed* win = hyp->GetWinner();
808 esdTr->UpdateTrackParams(hyp,AliESDtrack::kITSin);
810 // transfer module indices
814 //______________________________________________________________________________
815 Bool_t AliITSUTrackerGlo::RefitTrack(AliITSUTrackHyp* trc, Double_t rDest)
817 // refit track till radius rDest
818 AliITSUTrackHyp tmpTr;
820 double rCurr = Sqrt(trc->GetX()*trc->GetX() + trc->GetY()*trc->GetY());
821 int dir,lrStart,lrStop;
823 dir = rCurr<rDest ? 1 : -1;
824 lrStart = fITS->FindFirstLayerID(rCurr,dir);
825 lrStop = fITS->FindLastLayerID(rDest,dir);
826 printf("Refit %d: Lr:%d (%f) -> Lr:%d (%f)\n",dir,lrStart,rCurr, lrStop,rDest);
827 printf("Before refit: "); trc->AliExternalTrackParam::Print();
828 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));
830 trc->FetchClusterInfo(fClInfo);
831 fCurrMass = trc->GetMass();
832 tmpTr.AliKalmanTrack::operator=(*trc);
836 for (int ilr=lrStart;ilr!=lrStop;ilr+=dir) {
837 AliITSURecoLayer* lr = fITS->GetLayer(ilr);
838 if ( dir*(rCurr-lr->GetR(dir))>0) continue; // this layer is already passed
839 int ilrA2,ilrA = lr->GetActiveID();
840 // passive layer or active w/o hits will be traversed on the way to next cluster
841 if (!lr->IsActive() || fClInfo[ilrA2=(ilrA<<1)]<0) continue;
843 if (ilr!=lrStart && !TransportToLayer(&tmpTr,lrStart,ilr)) {
844 printf("Failed to transport %d -> %d\n",lrStart,ilr);
845 return kFALSE; // go to the entrance to the layer
849 // select the order in which possible 2 clusters (in case of the overlap) will be traversed and fitted
851 if (dir>0) { // clusters are stored in increasing radius order
852 iclLr[nclLr++]=fClInfo[ilrA2++];
853 if (fClInfo[ilrA2]>=0) iclLr[nclLr++]=fClInfo[ilrA2];
856 if ( fClInfo[ilrA2+1]>=0 ) iclLr[nclLr++]=fClInfo[ilrA2+1];
857 iclLr[nclLr++]=fClInfo[ilrA2];
860 for (int icl=0;icl<nclLr;icl++) {
861 AliITSUClusterPix* clus = (AliITSUClusterPix*)lr->GetCluster(iclLr[icl]);
862 AliITSURecoSens* sens = lr->GetSensorFromID(clus->GetVolumeId());
863 if (!tmpTr.Rotate(sens->GetPhiTF())) return kFALSE;
864 printf("Refit cl:%d on lr %d Need to go %.4f -> %.4f\n",icl,ilrA,tmpTr.GetX(),sens->GetXTF()+clus->GetX());
865 if (!PropagateSeed(&tmpTr,sens->GetXTF()+clus->GetX(),fCurrMass)) return kFALSE;
866 if (!tmpTr.Update(clus)) return kFALSE;
867 printf("AfterRefit: "); tmpTr.AliExternalTrackParam::Print();
871 // All clusters were succesfully fitted. Even if the track does not reach rDest, this is enough to validate it.
872 // Still, try to go as close as possible to rDest.
874 if (lrStart!=lrStop) {
875 printf("Going to last layer %d -> %d\n",lrStart,lrStop);
876 if (!TransportToLayer(&tmpTr,lrStart,lrStop)) return kTRUE;
877 if (!GoToExitFromLayer(&tmpTr,fITS->GetLayer(lrStop),dir)) return kFALSE; // go till the exit from layer
879 printf("On exit from last layer\n");
880 tmpTr.AliExternalTrackParam::Print();
881 // go to the destination radius
882 if (!tmpTr.GetXatLabR(rDest,rDest,GetBz(),dir)) return kTRUE;
883 if (!PropagateSeed(&tmpTr,rDest,fCurrMass, 100, kFALSE)) return kTRUE;
885 trc->AliKalmanTrack::operator=(tmpTr);
886 printf("After refit (now at lr %d): ",lrStart); trc->AliExternalTrackParam::Print();
890 //__________________________________________________________________
891 void AliITSUTrackerGlo::CookMCLabel(AliITSUTrackHyp* hyp)
895 const int kMaxLbPerCl = 3;
896 int lbID[kMaxLayers*6],lbStat[kMaxLayers*6];
897 Int_t lr,nLab=0,nCl=0;
898 AliITSUSeed *seed = hyp->GetWinner();
900 int clID = seed->GetLrCluster(lr);
902 AliCluster *cl = fITS->GetLayerActive(lr)->GetCluster(clID);
904 for (int imc=0;imc<kMaxLbPerCl;imc++) { // labels within single cluster
905 int trLb = cl->GetLabel(imc);
907 // search this mc track in already accounted ones
909 for (iLab=0;iLab<nLab;iLab++) if (lbID[iLab]==trLb) break;
910 if (iLab<nLab) lbStat[iLab]++;
915 } // loop over given cluster's labels
917 seed = (AliITSUSeed*)seed->GetParent();
918 } // loop over clusters
921 int maxLab=0,nTPCok=0;
922 AliESDtrack* esdTr = hyp->GetESDTrack();
923 int tpcLab = esdTr ? Abs(esdTr->GetTPCLabel()) : -kDummyLabel;
924 for (int ilb=nLab;ilb--;) {
925 int st = lbStat[ilb];
926 if (lbStat[maxLab]<st) maxLab=ilb;
927 if (tpcLab==lbID[ilb]) nTPCok += st;
929 hyp->SetFakeRatio(1.-float(nTPCok)/nCl);
930 hyp->SetLabel( nTPCok==nCl ? tpcLab : -tpcLab);
931 hyp->SetITSLabel( lbStat[maxLab]==nCl ? lbStat[maxLab] : -lbStat[maxLab]); // winner label
935 hyp->SetFakeRatio(-1.);
936 hyp->SetLabel( kDummyLabel );
937 hyp->SetITSLabel( kDummyLabel );