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 **************************************************************************/
16 //-------------------------------------------------------------------------
17 // Implementation of the ITS tracker class
18 // It reads AliITSclusterV2 clusters and creates AliITStrackV2 tracks
19 // and fills with them the ESD
20 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
21 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
22 //-------------------------------------------------------------------------
30 #include "AliITSgeom.h"
31 #include "AliITSRecPoint.h"
32 #include "AliTPCtrack.h"
34 #include "AliITSclusterV2.h"
35 #include "AliITStrackerV2.h"
37 ClassImp(AliITStrackerV2)
39 AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[kMaxLayer]; // ITS layers
41 AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
42 //--------------------------------------------------------------------
43 //This is the AliITStrackerV2 constructor
44 //--------------------------------------------------------------------
45 AliITSgeom *g=(AliITSgeom*)geom;
49 for (i=1; i<kMaxLayer+1; i++) {
50 Int_t nlad=g->GetNladders(i);
51 Int_t ndet=g->GetNdetectors(i);
53 g->GetTrans(i,1,1,x,y,z);
54 Double_t r=TMath::Sqrt(x*x + y*y);
55 Double_t poff=TMath::ATan2(y,x);
58 g->GetTrans(i,1,2,x,y,z);
59 r += TMath::Sqrt(x*x + y*y);
60 g->GetTrans(i,2,1,x,y,z);
61 r += TMath::Sqrt(x*x + y*y);
62 g->GetTrans(i,2,2,x,y,z);
63 r += TMath::Sqrt(x*x + y*y);
66 new (fgLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
68 for (Int_t j=1; j<nlad+1; j++) {
69 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
70 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
71 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
73 Double_t phi=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
75 if (i==1) phi+=TMath::Pi();
76 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
79 AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1);
80 new(&det) AliITSdetector(r,phi);
89 fConstraint[0]=1; fConstraint[1]=0;
91 Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV};
94 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=kLayersNotToSkip[i];
95 fLastLayerToTrackTo=kLastLayerToTrackTo;
99 void AliITStrackerV2::SetLayersNotToSkip(Int_t *l) {
100 //--------------------------------------------------------------------
101 //This function set masks of the layers which must be not skipped
102 //--------------------------------------------------------------------
103 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=l[i];
106 Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
107 //--------------------------------------------------------------------
108 //This function loads ITS clusters
109 //--------------------------------------------------------------------
110 TBranch *branch=cTree->GetBranch("Clusters");
112 Error("LoadClusters"," can't get the branch !\n");
116 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
117 branch->SetAddress(&clusters);
120 for (Int_t i=0; i<kMaxLayer; i++) {
121 Int_t ndet=fgLayers[i].GetNdetectors();
122 Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
123 for (; j<jmax; j++) {
124 if (!cTree->GetEvent(j)) continue;
125 Int_t ncl=clusters->GetEntriesFast();
127 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
128 fgLayers[i].InsertCluster(new AliITSclusterV2(*c));
132 fgLayers[i].ResetRoad(); //road defined by the cluster density
138 void AliITStrackerV2::UnloadClusters() {
139 //--------------------------------------------------------------------
140 //This function unloads ITS clusters
141 //--------------------------------------------------------------------
142 for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
145 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
146 //--------------------------------------------------------------------
147 // Correction for the material between the TPC and the ITS
148 // (should it belong to the TPC code ?)
149 //--------------------------------------------------------------------
150 Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ?
151 Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
152 Double_t yr=12.8, dr=0.03; // rods ?
153 Double_t zm=0.2, dm=0.40; // membrane
154 //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
155 Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
157 if (t->GetX() > riw) {
158 if (!t->PropagateTo(riw,diw,x0iw)) return 1;
159 if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
160 if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
161 if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
162 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
163 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r);
164 if (!t->PropagateTo(rs,ds)) return 1;
165 } else if (t->GetX() < rs) {
166 if (!t->PropagateTo(rs,-ds)) return 1;
167 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
168 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r);
169 if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
170 if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
172 ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
179 Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
180 //--------------------------------------------------------------------
181 // This functions reconstructs ITS tracks
182 // The clusters must be already loaded !
183 //--------------------------------------------------------------------
184 TObjArray itsTracks(15000);
186 {/* Read ESD tracks */
187 Int_t nentr=event->GetNumberOfTracks();
188 Info("Clusters2Tracks", "Number of ESD tracks: %d\n", nentr);
190 AliESDtrack *esd=event->GetTrack(nentr);
192 if ((esd->GetStatus()&AliESDtrack::kTPCin)==0) continue;
193 if (esd->GetStatus()&AliESDtrack::kTPCout) continue;
194 if (esd->GetStatus()&AliESDtrack::kITSin) continue;
198 t=new AliITStrackV2(*esd);
199 } catch (const Char_t *msg) {
200 Warning("Clusters2Tracks",msg);
204 if (TMath::Abs(t->GetD())>4) {
209 if (CorrectForDeadZoneMaterial(t)!=0) {
210 Warning("Clusters2Tracks",
211 "failed to correct for the material in the dead zone !\n");
215 itsTracks.AddLast(t);
217 } /* End Read ESD tracks */
220 Int_t nentr=itsTracks.GetEntriesFast();
221 fTrackHypothesys.Expand(nentr);
223 for (fPass=0; fPass<2; fPass++) {
224 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
225 for (Int_t i=0; i<nentr; i++) {
226 fCurrentEsdTrack = i;
227 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
228 if (t==0) continue; //this track has been already tracked
229 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
231 ResetTrackToFollow(*t);
234 for (FollowProlongation(); fI<kMaxLayer; fI++) {
235 while (TakeNextProlongation()) {
236 FollowProlongation();
240 CompressTrackHypothesys(fCurrentEsdTrack,5); //MI change
243 if (fBestTrack.GetNumberOfClusters() == 0) continue;
245 if (fConstraint[fPass]) {
246 ResetTrackToFollow(*t);
247 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
251 fBestTrack.SetLabel(tpcLabel);
252 fBestTrack.CookdEdx();
253 CookLabel(&fBestTrack,0.); //For comparison only
254 fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
255 //UseClusters(&fBestTrack);
256 delete itsTracks.RemoveAt(i);
259 AliITStrackV2 * besttrack = GetBestHypothesys(fCurrentEsdTrack,t);
260 if (!besttrack) continue;
261 besttrack->SetLabel(tpcLabel);
262 besttrack->CookdEdx();
263 CookLabel(besttrack,0.); //For comparison only
264 besttrack->UpdateESDtrack(AliESDtrack::kITSin);
265 delete itsTracks.RemoveAt(i);
272 Int_t entries = fTrackHypothesys.GetEntriesFast();
273 for (Int_t ientry=0;ientry<entries;ientry++){
274 TObjArray * array =(TObjArray*)fTrackHypothesys.UncheckedAt(ientry);
275 if (array) array->Delete();
276 delete fTrackHypothesys.RemoveAt(ientry);
279 fTrackHypothesys.Delete();
280 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
285 Int_t AliITStrackerV2::Clusters2Tracks(TTree *tpcTree, TTree *itsTree) {
286 //--------------------------------------------------------------------
287 // This functions reconstructs ITS tracks
288 // The clusters must be already loaded !
289 //--------------------------------------------------------------------
290 Int_t nentr=0; TObjArray itsTracks(15000);
292 Warning("Clusters2Tracks(TTree *, TTree *)",
293 "Will be removed soon ! Use Clusters2Tracks(AliESD *) instead.");
295 {/* Read TPC tracks */
296 AliTPCtrack *itrack=new AliTPCtrack;
297 TBranch *branch=tpcTree->GetBranch("tracks");
299 Error("Clusters2Tracks","Can't get the branch !");
302 tpcTree->SetBranchAddress("tracks",&itrack);
303 nentr=(Int_t)tpcTree->GetEntries();
305 Info("Clusters2Tracks","Number of TPC tracks: %d\n",nentr);
307 for (Int_t i=0; i<nentr; i++) {
308 tpcTree->GetEvent(i);
311 t=new AliITStrackV2(*itrack);
312 } catch (const Char_t *msg) {
313 Warning("Clusters2Tracks",msg);
317 if (TMath::Abs(t->GetD())>4) {
322 if (CorrectForDeadZoneMaterial(t)!=0) {
323 Warning("Clusters2Tracks",
324 "failed to correct for the material in the dead zone !\n");
329 itsTracks.AddLast(t);
334 nentr=itsTracks.GetEntriesFast();
337 AliITStrackV2 *otrack=&fBestTrack;
338 TBranch *branch=itsTree->GetBranch("tracks");
339 if (!branch) itsTree->Branch("tracks","AliITStrackV2",&otrack,32000,3);
340 else branch->SetAddress(&otrack);
342 for (fPass=0; fPass<2; fPass++) {
343 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
344 for (Int_t i=0; i<nentr; i++) {
345 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
346 if (t==0) continue; //this track has been already tracked
347 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
349 ResetTrackToFollow(*t);
352 for (FollowProlongation(); fI<kMaxLayer; fI++) {
353 while (TakeNextProlongation()) FollowProlongation();
356 if (fBestTrack.GetNumberOfClusters() == 0) continue;
358 if (fConstraint[fPass]) {
359 ResetTrackToFollow(*t);
360 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
364 fBestTrack.SetLabel(tpcLabel);
365 fBestTrack.CookdEdx();
366 CookLabel(&fBestTrack,0.); //For comparison only
368 //UseClusters(&fBestTrack);
369 delete itsTracks.RemoveAt(i);
373 nentr=(Int_t)itsTree->GetEntries();
374 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",nentr);
381 Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
382 //--------------------------------------------------------------------
383 // This functions propagates reconstructed ITS tracks back
384 // The clusters must be loaded !
385 //--------------------------------------------------------------------
386 Int_t nentr=event->GetNumberOfTracks();
387 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
390 for (Int_t i=0; i<nentr; i++) {
391 AliESDtrack *esd=event->GetTrack(i);
393 if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
394 if (esd->GetStatus()&AliESDtrack::kITSout) continue;
398 t=new AliITStrackV2(*esd);
399 } catch (const Char_t *msg) {
400 Warning("PropagateBack",msg);
405 ResetTrackToFollow(*t);
407 // propagete to vertex [SR, GSI 17.02.2003]
408 // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
409 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
410 if (fTrackToFollow.PropagateToVertex()) {
411 fTrackToFollow.StartTimeIntegral();
413 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
416 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
417 if (RefitAt(49.,&fTrackToFollow,t)) {
418 if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
419 Warning("PropagateBack",
420 "failed to correct for the material in the dead zone !\n");
424 fTrackToFollow.SetLabel(t->GetLabel());
425 //fTrackToFollow.CookdEdx();
426 CookLabel(&fTrackToFollow,0.); //For comparison only
427 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
428 //UseClusters(&fTrackToFollow);
434 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
439 Int_t AliITStrackerV2::RefitInward(AliESD *event) {
440 //--------------------------------------------------------------------
441 // This functions refits ITS tracks using the
442 // "inward propagated" TPC tracks
443 // The clusters must be loaded !
444 //--------------------------------------------------------------------
445 Int_t nentr=event->GetNumberOfTracks();
446 Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
449 for (Int_t i=0; i<nentr; i++) {
450 AliESDtrack *esd=event->GetTrack(i);
452 if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
453 if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
454 if (esd->GetStatus()&AliESDtrack::kTPCout)
455 if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
459 t=new AliITStrackV2(*esd);
460 } catch (const Char_t *msg) {
461 Warning("RefitInward",msg);
466 if (CorrectForDeadZoneMaterial(t)!=0) {
467 Warning("RefitInward",
468 "failed to correct for the material in the dead zone !\n");
473 ResetTrackToFollow(*t);
474 fTrackToFollow.ResetClusters();
476 if ((esd->GetStatus()&AliESDtrack::kTPCin)==0)
477 fTrackToFollow.ResetCovariance();
480 if (RefitAt(3.7, &fTrackToFollow, t)) {
481 fTrackToFollow.SetLabel(t->GetLabel());
482 fTrackToFollow.CookdEdx();
483 CookLabel(&fTrackToFollow,0.0); //For comparison only
485 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe
486 Double_t a=fTrackToFollow.GetAlpha();
487 Double_t cs=TMath::Cos(a),sn=TMath::Sin(a);
488 Double_t xv= GetX()*cs + GetY()*sn;
489 Double_t yv=-GetX()*sn + GetY()*cs;
491 Double_t c=fTrackToFollow.GetC(), snp=fTrackToFollow.GetSnp();
492 Double_t x=fTrackToFollow.GetX(), y=fTrackToFollow.GetY();
493 Double_t tgfv=-(c*(x-xv)-snp)/(c*(y-yv) + TMath::Sqrt(1.-snp*snp));
494 Double_t fv=TMath::ATan(tgfv);
496 cs=TMath::Cos(fv); sn=TMath::Sin(fv);
498 yv=-xv*sn + yv*cs; xv=x;
500 if (fTrackToFollow.Propagate(fv+a,xv)) {
501 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
502 //UseClusters(&fTrackToFollow);
504 AliITSclusterV2 c; c.SetY(yv); c.SetZ(GetZ());
505 c.SetSigmaY2(GetSigmaY()*GetSigmaY());
506 c.SetSigmaZ2(GetSigmaZ()*GetSigmaZ());
507 Double_t chi2=fTrackToFollow.GetPredictedChi2(&c);
509 if (fTrackToFollow.Update(&c,-chi2,0))
510 fTrackToFollow.SetConstrainedESDtrack(chi2);
519 Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
524 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
525 //--------------------------------------------------------------------
526 // Return pointer to a given cluster
527 //--------------------------------------------------------------------
528 Int_t l=(index & 0xf0000000) >> 28;
529 Int_t c=(index & 0x0fffffff) >> 00;
530 return fgLayers[l].GetCluster(c);
534 void AliITStrackerV2::FollowProlongation() {
535 //--------------------------------------------------------------------
536 //This function finds a track prolongation
537 //--------------------------------------------------------------------
539 while (fI>fLastLayerToTrackTo) {
542 AliITSlayer &layer=fgLayers[i];
543 AliITStrackV2 &track=fTracks[i];
545 Double_t r=layer.GetR();
548 Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
549 Double_t d=0.0034, x0=38.6;
550 if (i==1) {rs=9.; d=0.0097; x0=42;}
551 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
552 //Warning("FollowProlongation","propagation failed !\n");
559 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
560 //Warning("FollowProlongation","failed to estimate track !\n");
563 Double_t phi=TMath::ATan2(y,x);
565 Int_t idet=layer.FindDetectorIndex(phi,z);
567 //Warning("FollowProlongation","failed to find a detector !\n");
571 //propagate to the intersection
572 const AliITSdetector &det=layer.GetDetector(idet);
574 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
575 //Warning("FollowProlongation","propagation failed !\n");
578 fTrackToFollow.SetDetectorIndex(idet);
580 //Select possible prolongations and store the current track estimation
581 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
582 Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
583 Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
584 Double_t road=layer.GetRoad();
585 if (dz*dy>road*road) {
586 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
587 dz=road*scz; dy=road*scy;
590 //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
591 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
593 //Warning("FollowProlongation","too broad road in Z !\n");
597 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
599 //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
600 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
602 //Warning("FollowProlongation","too broad road in Y !\n");
606 Double_t zmin=track.GetZ() - dz;
607 Double_t zmax=track.GetZ() + dz;
608 Double_t ymin=track.GetY() + r*phi - dy;
609 Double_t ymax=track.GetY() + r*phi + dy;
610 layer.SelectClusters(zmin,zmax,ymin,ymax);
613 //take another prolongation
614 if (!TakeNextProlongation()){
616 if (fLayersNotToSkip[fI]||skipped>1) return;
620 //deal with the best track
621 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
622 Int_t nclb=fBestTrack.GetNumberOfClusters();
624 if ( (ncl>1) && (nclb>3) && ((fTrackToFollow.GetChi2()/ncl)<(2*fBestTrack.GetChi2()/(nclb))))
625 AddTrackHypothesys(new AliITStrackV2(fTrackToFollow), fCurrentEsdTrack);
626 else if (ncl>3) AddTrackHypothesys(new AliITStrackV2(fTrackToFollow), fCurrentEsdTrack);
628 Double_t chi2=fTrackToFollow.GetChi2();
629 if (chi2/ncl < kChi2PerCluster) {
633 if ( (ncl == nclb) && chi2 < fBestTrack.GetChi2()) {
641 Int_t AliITStrackerV2::TakeNextProlongation() {
642 //--------------------------------------------------------------------
643 // This function takes another track prolongation
645 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
646 //--------------------------------------------------------------------
647 AliITSlayer &layer=fgLayers[fI];
648 ResetTrackToFollow(fTracks[fI]);
650 Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
651 Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
652 Double_t road=layer.GetRoad();
653 if (dz*dy>road*road) {
654 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
655 dz=road*scz; dy=road*scy;
658 const AliITSclusterV2 *c=0; Int_t ci=-1;
659 Double_t chi2=12345.;
660 while ((c=layer.GetNextCluster(ci))!=0) {
661 Int_t idet=c->GetDetectorIndex();
663 if (fTrackToFollow.GetDetectorIndex()!=idet) {
664 const AliITSdetector &det=layer.GetDetector(idet);
665 ResetTrackToFollow(fTracks[fI]);
666 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
667 //Warning("TakeNextProlongation","propagation failed !\n");
670 fTrackToFollow.SetDetectorIndex(idet);
671 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
674 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
675 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
677 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
680 if (chi2>=kMaxChi2) return 0;
683 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
684 //Warning("TakeNextProlongation","filtering failed !\n");
688 if (fTrackToFollow.GetNumberOfClusters()>1)
689 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
692 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
696 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
697 fTrackToFollow.CorrectForMaterial(d,x0);
700 if (fConstraint[fPass]) {
701 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
702 Double_t xyz[]={GetX(),GetY(),GetZ()};
703 Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
704 fTrackToFollow.Improve(d,xyz,ers);
711 AliITStrackerV2::AliITSlayer::AliITSlayer() {
712 //--------------------------------------------------------------------
713 //default AliITSlayer constructor
714 //--------------------------------------------------------------------
719 AliITStrackerV2::AliITSlayer::
720 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
721 //--------------------------------------------------------------------
722 //main AliITSlayer constructor
723 //--------------------------------------------------------------------
724 fR=r; fPhiOffset=p; fZOffset=z;
725 fNladders=nl; fNdetectors=nd;
726 fDetectors=new AliITSdetector[fNladders*fNdetectors];
731 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
734 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
735 //--------------------------------------------------------------------
736 // AliITSlayer destructor
737 //--------------------------------------------------------------------
739 for (Int_t i=0; i<fN; i++) delete fClusters[i];
742 void AliITStrackerV2::AliITSlayer::ResetClusters() {
743 //--------------------------------------------------------------------
744 // This function removes loaded clusters
745 //--------------------------------------------------------------------
746 for (Int_t i=0; i<fN; i++) delete fClusters[i];
751 void AliITStrackerV2::AliITSlayer::ResetRoad() {
752 //--------------------------------------------------------------------
753 // This function calculates the road defined by the cluster density
754 //--------------------------------------------------------------------
756 for (Int_t i=0; i<fN; i++) {
757 if (TMath::Abs(fClusters[i]->GetZ())<fR) n++;
759 if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
762 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
763 //--------------------------------------------------------------------
764 //This function adds a cluster to this layer
765 //--------------------------------------------------------------------
766 if (fN==kMaxClusterPerLayer) {
767 ::Error("InsertCluster","Too many clusters !\n");
771 if (fN==0) {fClusters[fN++]=c; return 0;}
772 Int_t i=FindClusterIndex(c->GetZ());
773 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
774 fClusters[i]=c; fN++;
779 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
780 //--------------------------------------------------------------------
781 // This function returns the index of the nearest cluster
782 //--------------------------------------------------------------------
784 if (z <= fClusters[0]->GetZ()) return 0;
785 if (z > fClusters[fN-1]->GetZ()) return fN;
786 Int_t b=0, e=fN-1, m=(b+e)/2;
787 for (; b<e; m=(b+e)/2) {
788 if (z > fClusters[m]->GetZ()) b=m+1;
794 void AliITStrackerV2::AliITSlayer::
795 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
796 //--------------------------------------------------------------------
797 // This function sets the "window"
798 //--------------------------------------------------------------------
799 fI=FindClusterIndex(zmin); fZmax=zmax;
800 Double_t circle=2*TMath::Pi()*fR;
801 if (ymax>circle) { ymax-=circle; ymin-=circle; }
802 fYmin=ymin; fYmax=ymax;
805 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
806 //--------------------------------------------------------------------
807 // This function returns clusters within the "window"
808 //--------------------------------------------------------------------
809 const AliITSclusterV2 *cluster=0;
810 for (Int_t i=fI; i<fN; i++) {
811 const AliITSclusterV2 *c=fClusters[i];
812 if (c->GetZ() > fZmax) break;
813 if (c->IsUsed()) continue;
814 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
815 Double_t y=fR*det.GetPhi() + c->GetY();
817 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
818 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
820 if (y<fYmin) continue;
821 if (y>fYmax) continue;
830 Int_t AliITStrackerV2::AliITSlayer::
831 FindDetectorIndex(Double_t phi, Double_t z) const {
832 //--------------------------------------------------------------------
833 //This function finds the detector crossed by the track
834 //--------------------------------------------------------------------
835 Double_t dphi=-(phi-fPhiOffset);
836 if (dphi < 0) dphi += 2*TMath::Pi();
837 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
838 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
839 if (np>=fNladders) np-=fNladders;
840 if (np<0) np+=fNladders;
842 Double_t dz=fZOffset-z;
843 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
844 if (nz>=fNdetectors) return -1;
847 return np*fNdetectors + nz;
851 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
853 //--------------------------------------------------------------------
854 //This function returns the layer thickness at this point (units X0)
855 //--------------------------------------------------------------------
859 if (43<fR&&fR<45) { //SSD2
862 if (TMath::Abs(y-0.00)>3.40) d+=dd;
863 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
864 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
865 for (Int_t i=0; i<12; i++) {
866 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
867 if (TMath::Abs(y-0.00)>3.40) d+=dd;
871 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
872 if (TMath::Abs(y-0.00)>3.40) d+=dd;
876 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
877 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
880 if (37<fR&&fR<41) { //SSD1
883 if (TMath::Abs(y-0.00)>3.40) d+=dd;
884 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
885 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
886 for (Int_t i=0; i<11; i++) {
887 if (TMath::Abs(z-3.9*i)<0.15) {
888 if (TMath::Abs(y-0.00)>3.40) d+=dd;
892 if (TMath::Abs(z+3.9*i)<0.15) {
893 if (TMath::Abs(y-0.00)>3.40) d+=dd;
897 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
898 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
901 if (13<fR&&fR<26) { //SDD
904 if (TMath::Abs(y-0.00)>3.30) d+=dd;
906 if (TMath::Abs(y-1.80)<0.55) {
908 for (Int_t j=0; j<20; j++) {
909 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
910 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
913 if (TMath::Abs(y+1.80)<0.55) {
915 for (Int_t j=0; j<20; j++) {
916 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
917 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
921 for (Int_t i=0; i<4; i++) {
922 if (TMath::Abs(z-7.3*i)<0.60) {
924 if (TMath::Abs(y-0.00)>3.30) d+=dd;
927 if (TMath::Abs(z+7.3*i)<0.60) {
929 if (TMath::Abs(y-0.00)>3.30) d+=dd;
934 if (6<fR&&fR<8) { //SPD2
935 Double_t dd=0.0063; x0=21.5;
937 if (TMath::Abs(y-3.08)>0.5) d+=dd;
938 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
939 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
941 if (3<fR&&fR<5) { //SPD1
942 Double_t dd=0.0063; x0=21.5;
944 if (TMath::Abs(y+0.21)>0.6) d+=dd;
945 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
946 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
952 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
954 //--------------------------------------------------------------------
955 //Returns the thickness between the current layer and the vertex (units X0)
956 //--------------------------------------------------------------------
957 Double_t d=0.0028*3*3; //beam pipe
960 Double_t xn=fgLayers[fI].GetR();
961 for (Int_t i=0; i<fI; i++) {
962 Double_t xi=fgLayers[i].GetR();
963 d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
972 Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
979 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
980 //--------------------------------------------------------------------
981 // This function returns number of clusters within the "window"
982 //--------------------------------------------------------------------
984 for (Int_t i=fI; i<fN; i++) {
985 const AliITSclusterV2 *c=fClusters[i];
986 if (c->GetZ() > fZmax) break;
987 if (c->IsUsed()) continue;
988 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
989 Double_t y=fR*det.GetPhi() + c->GetY();
991 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
992 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
994 if (y<fYmin) continue;
995 if (y>fYmax) continue;
1002 AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
1003 //--------------------------------------------------------------------
1004 // This function refits the track "t" at the position "x" using
1005 // the clusters from "c"
1006 //--------------------------------------------------------------------
1007 Int_t index[kMaxLayer];
1009 for (k=0; k<kMaxLayer; k++) index[k]=-1;
1010 Int_t nc=c->GetNumberOfClusters();
1011 for (k=0; k<nc; k++) {
1012 Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
1016 Int_t from, to, step;
1017 if (xx > t->GetX()) {
1018 from=0; to=kMaxLayer;
1021 from=kMaxLayer-1; to=-1;
1025 for (Int_t i=from; i != to; i += step) {
1026 AliITSlayer &layer=fgLayers[i];
1027 Double_t r=layer.GetR();
1030 Double_t hI=i-0.5*step;
1031 if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {
1032 Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
1033 Double_t d=0.0034, x0=38.6;
1034 if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1035 if (!t->PropagateTo(rs,-step*d,x0)) {
1041 // remember old position [SR, GSI 18.02.2003]
1042 Double_t oldX=0., oldY=0., oldZ=0.;
1043 if (t->IsStartedTimeIntegral() && step==1) {
1044 t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1049 if (!t->GetGlobalXYZat(r,x,y,z)) {
1052 Double_t phi=TMath::ATan2(y,x);
1053 Int_t idet=layer.FindDetectorIndex(phi,z);
1057 const AliITSdetector &det=layer.GetDetector(idet);
1059 if (!t->Propagate(phi,det.GetR())) {
1062 t->SetDetectorIndex(idet);
1064 const AliITSclusterV2 *cl=0;
1065 Double_t maxchi2=kMaxChi2;
1069 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
1070 if (idet != c->GetDetectorIndex()) {
1071 idet=c->GetDetectorIndex();
1072 const AliITSdetector &det=layer.GetDetector(idet);
1073 if (!t->Propagate(det.GetPhi(),det.GetR())) {
1076 t->SetDetectorIndex(idet);
1078 Double_t chi2=t->GetPredictedChi2(c);
1088 if (t->GetNumberOfClusters()>2) {
1089 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1090 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1091 Double_t zmin=t->GetZ() - dz;
1092 Double_t zmax=t->GetZ() + dz;
1093 Double_t ymin=t->GetY() + phi*r - dy;
1094 Double_t ymax=t->GetY() + phi*r + dy;
1095 layer.SelectClusters(zmin,zmax,ymin,ymax);
1097 const AliITSclusterV2 *c=0; Int_t ci=-1;
1098 while ((c=layer.GetNextCluster(ci))!=0) {
1099 if (idet != c->GetDetectorIndex()) continue;
1100 Double_t chi2=t->GetPredictedChi2(c);
1101 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1106 if (!t->Update(cl,maxchi2,idx)) {
1109 t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1114 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1115 t->CorrectForMaterial(-step*d,x0);
1118 // track time update [SR, GSI 17.02.2003]
1119 if (t->IsStartedTimeIntegral() && step==1) {
1120 Double_t newX, newY, newZ;
1121 t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1122 Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) +
1123 (oldZ-newZ)*(oldZ-newZ);
1124 t->AddTimeStep(TMath::Sqrt(dL2));
1130 if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1134 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1135 //--------------------------------------------------------------------
1136 // This function marks clusters assigned to the track
1137 //--------------------------------------------------------------------
1138 AliTracker::UseClusters(t,from);
1140 AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1141 //if (c->GetQ()>2) c->Use();
1142 if (c->GetSigmaZ2()>0.1) c->Use();
1143 c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1144 //if (c->GetQ()>2) c->Use();
1145 if (c->GetSigmaZ2()>0.1) c->Use();
1150 void AliITStrackerV2::AddTrackHypothesys(AliITStrackV2 * track, Int_t esdindex)
1152 //------------------------------------------------------------------
1153 // add track to the list of hypothesys
1154 //------------------------------------------------------------------
1156 if (esdindex>=fTrackHypothesys.GetEntriesFast()) fTrackHypothesys.Expand(esdindex*2+10);
1158 TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1160 array = new TObjArray(10);
1161 fTrackHypothesys.AddAt(array,esdindex);
1163 array->AddLast(track);
1166 void AliITStrackerV2::CompressTrackHypothesys(Int_t esdindex, Int_t maxsize)
1168 //-------------------------------------------------------------------
1169 // compress array of track hypothesys
1170 // keep only maxsize best hypothesys
1171 //-------------------------------------------------------------------
1172 if (esdindex>fTrackHypothesys.GetEntriesFast()) return;
1173 if (! (fTrackHypothesys.At(esdindex)) ) return;
1174 TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1175 Int_t entries = array->GetEntries();
1176 if (entries<maxsize) return;
1177 Float_t * chi2 = new Float_t[entries];
1178 Int_t * index = new Int_t[entries];
1181 for (Int_t i=0;i<array->GetEntriesFast();i++){
1182 AliITStrackV2 * track = (AliITStrackV2*)array->At(i);
1184 chi2[current] = track->GetNumberOfClusters()+0.001+track->GetChi2()/track->GetNumberOfClusters();
1186 chi2[current] = 100000000;
1189 TMath::Sort(current,chi2,index,kFALSE);
1190 TObjArray * newarray = new TObjArray();
1191 maxsize = TMath::Min(current, maxsize);
1192 for (Int_t i=0;i<maxsize;i++){
1193 newarray->AddLast(array->RemoveAt(index[i]));
1196 delete fTrackHypothesys.RemoveAt(esdindex);
1197 fTrackHypothesys.AddAt(newarray,esdindex);
1206 AliITStrackV2 * AliITStrackerV2::GetBestHypothesys(Int_t esdindex, AliITStrackV2 * original)
1208 //-------------------------------------------------------------
1209 // try to find best hypothesy
1210 // currently - minimal chi2 of track+backpropagated track+matching to the tpc track
1211 //-------------------------------------------------------------
1212 TObjArray * array = (TObjArray*) fTrackHypothesys.At(esdindex);
1213 if (!array) return 0;
1214 Int_t entries = array->GetEntriesFast();
1215 if (!entries) return 0;
1216 Float_t minchi2 = 100000;
1218 AliITStrackV2 * besttrack=0;
1220 // calculate "weight of the cluster"
1222 Int_t clusterindex[6][100];
1223 Double_t clusterweight[6][100];
1224 for (Int_t ilayer=0;ilayer<6;ilayer++)
1225 for (Int_t icluster=0;icluster<100;icluster++){
1226 clusterindex[ilayer][icluster] = -1;
1227 clusterweight[ilayer][icluster] = 0;
1229 //printf("%d\t%d\n",esdindex, entries);
1232 for (Int_t itrack=0;itrack<entries;itrack++){
1233 AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1234 if (!track) continue;
1235 sumchi2 +=track->GetChi2();
1237 for (Int_t itrack=0;itrack<entries;itrack++){
1238 AliITStrackV2 * track = (AliITStrackV2*)array->At(itrack);
1239 if (!track) continue;
1240 for (Int_t icluster=0;icluster<track->GetNumberOfClusters();icluster++){
1241 Int_t tindex = track->GetClusterIndex(icluster);
1242 Int_t ilayer = (tindex & 0xf0000000) >> 28;
1243 if (tindex<0) continue;
1246 for (cindex=0;cindex<100;cindex++){
1247 if (clusterindex[ilayer][cindex]<0) break;
1248 if (clusterindex[ilayer][cindex]==tindex) break;
1250 if (cindex>100) break;
1251 if (clusterindex[ilayer][cindex]!=tindex) clusterindex[ilayer][cindex] = tindex;
1252 clusterweight[ilayer][cindex]+= track->GetChi2()/sumchi2;
1256 for (Int_t i=0;i<entries;i++){
1257 AliITStrackV2 * track = (AliITStrackV2*)array->At(i);
1258 if (!track) continue;
1259 AliITStrackV2 * backtrack = new AliITStrackV2(*track);
1260 backtrack->ResetCovariance();
1261 backtrack->ResetClusters();
1262 Double_t x = original->GetX();
1263 RefitAt(x,backtrack,track);
1264 if (track->GetNumberOfClusters()<maxn) {
1268 Double_t deltac = backtrack->GetC()-original->GetC();
1269 Double_t deltatgl = backtrack->GetTgl()-original->GetTgl();
1270 Double_t poolc = (deltac*deltac)/(original->fC44);
1271 Double_t pooltgl = (deltatgl*deltatgl)/(original->fC33);
1272 Double_t chi2 = 0.5*(track->GetChi2()+backtrack->GetChi2())+poolc+pooltgl;
1274 if (track->GetNumberOfClusters()>maxn){
1276 maxn = track->GetNumberOfClusters();
1281 if (chi2 < minchi2){
1289 AliITStrackV2 * track2 = new AliITStrackV2(*original);
1292 for (Int_t icluster=0;icluster<besttrack->GetNumberOfClusters();icluster++){
1293 Int_t index = besttrack->GetClusterIndex(icluster);
1294 Int_t ilayer = (index & 0xf0000000) >> 28;
1295 AliITSclusterV2 *c = (AliITSclusterV2*)GetCluster(index);
1299 track2->fIndex[current] = index;
1302 for (Int_t icluster=0;icluster<100;icluster++){
1303 // sign non "doubt" clusters as used
1304 if (clusterindex[ilayer][icluster]!=index) continue;
1305 if (clusterweight[ilayer][icluster]>0.9){
1307 track2->fIndex[current] = index;
1311 if (c->IsUsed()) continue;
1316 // besttrack = new AliITStrackV2(*original);
1317 //RefitAt(3.7, besttrack, track2);