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();
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 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
227 if (t==0) continue; //this track has been already tracked
228 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
230 ResetTrackToFollow(*t);
233 for (FollowProlongation(); fI<kMaxLayer; fI++) {
234 while (TakeNextProlongation()) FollowProlongation();
237 if (fBestTrack.GetNumberOfClusters() == 0) continue;
239 if (fConstraint[fPass]) {
240 ResetTrackToFollow(*t);
241 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
245 fBestTrack.SetLabel(tpcLabel);
246 fBestTrack.CookdEdx();
247 CookLabel(&fBestTrack,0.); //For comparison only
248 fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
249 UseClusters(&fBestTrack);
250 delete itsTracks.RemoveAt(i);
257 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
262 Int_t AliITStrackerV2::Clusters2Tracks(TTree *tpcTree, TTree *itsTree) {
263 //--------------------------------------------------------------------
264 // This functions reconstructs ITS tracks
265 // The clusters must be already loaded !
266 //--------------------------------------------------------------------
267 Int_t nentr=0; TObjArray itsTracks(15000);
269 Warning("Clusters2Tracks(TTree *, TTree *)",
270 "Will be removed soon ! Use Clusters2Tracks(AliESD *) instead.");
272 {/* Read TPC tracks */
273 AliTPCtrack *itrack=new AliTPCtrack;
274 TBranch *branch=tpcTree->GetBranch("tracks");
276 Error("Clusters2Tracks","Can't get the branch !");
279 tpcTree->SetBranchAddress("tracks",&itrack);
280 nentr=(Int_t)tpcTree->GetEntries();
282 Info("Clusters2Tracks","Number of TPC tracks: %d\n",nentr);
284 for (Int_t i=0; i<nentr; i++) {
285 tpcTree->GetEvent(i);
288 t=new AliITStrackV2(*itrack);
289 } catch (const Char_t *msg) {
290 Warning("Clusters2Tracks",msg);
294 if (TMath::Abs(t->GetD())>4) {
299 if (CorrectForDeadZoneMaterial(t)!=0) {
300 Warning("Clusters2Tracks",
301 "failed to correct for the material in the dead zone !\n");
306 itsTracks.AddLast(t);
311 nentr=itsTracks.GetEntriesFast();
314 AliITStrackV2 *otrack=&fBestTrack;
315 TBranch *branch=itsTree->GetBranch("tracks");
316 if (!branch) itsTree->Branch("tracks","AliITStrackV2",&otrack,32000,3);
317 else branch->SetAddress(&otrack);
319 for (fPass=0; fPass<2; fPass++) {
320 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
321 for (Int_t i=0; i<nentr; i++) {
322 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
323 if (t==0) continue; //this track has been already tracked
324 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
326 ResetTrackToFollow(*t);
329 for (FollowProlongation(); fI<kMaxLayer; fI++) {
330 while (TakeNextProlongation()) FollowProlongation();
333 if (fBestTrack.GetNumberOfClusters() == 0) continue;
335 if (fConstraint[fPass]) {
336 ResetTrackToFollow(*t);
337 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
341 fBestTrack.SetLabel(tpcLabel);
342 fBestTrack.CookdEdx();
343 CookLabel(&fBestTrack,0.); //For comparison only
345 UseClusters(&fBestTrack);
346 delete itsTracks.RemoveAt(i);
350 nentr=(Int_t)itsTree->GetEntries();
351 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",nentr);
358 Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
359 //--------------------------------------------------------------------
360 // This functions propagates reconstructed ITS tracks back
361 // The clusters must be loaded !
362 //--------------------------------------------------------------------
363 Int_t nentr=event->GetNumberOfTracks();
364 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
367 for (Int_t i=0; i<nentr; i++) {
368 AliESDtrack *esd=event->GetTrack(i);
370 if ((esd->GetStatus()&AliESDtrack::kITSin)==0) continue;
371 if (esd->GetStatus()&AliESDtrack::kITSout) continue;
375 t=new AliITStrackV2(*esd);
376 } catch (const Char_t *msg) {
377 Warning("PropagateBack",msg);
382 ResetTrackToFollow(*t);
384 // propagete to vertex [SR, GSI 17.02.2003]
385 // Start Time measurement [SR, GSI 17.02.2003], corrected by I.Belikov
386 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {
387 if (fTrackToFollow.PropagateToVertex()) {
388 fTrackToFollow.StartTimeIntegral();
390 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
393 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
394 if (RefitAt(49.,&fTrackToFollow,t)) {
395 if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
396 Warning("PropagateBack",
397 "failed to correct for the material in the dead zone !\n");
401 fTrackToFollow.SetLabel(t->GetLabel());
402 //fTrackToFollow.CookdEdx();
403 CookLabel(&fTrackToFollow,0.); //For comparison only
404 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
405 UseClusters(&fTrackToFollow);
411 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
416 Int_t AliITStrackerV2::RefitInward(AliESD *event) {
417 //--------------------------------------------------------------------
418 // This functions refits ITS tracks using the
419 // "inward propagated" TPC tracks
420 // The clusters must be loaded !
421 //--------------------------------------------------------------------
422 Int_t nentr=event->GetNumberOfTracks();
423 Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
426 for (Int_t i=0; i<nentr; i++) {
427 AliESDtrack *esd=event->GetTrack(i);
429 if ((esd->GetStatus()&AliESDtrack::kITSout) == 0) continue;
430 if (esd->GetStatus()&AliESDtrack::kITSrefit) continue;
431 if (esd->GetStatus()&AliESDtrack::kTPCout)
432 if ((esd->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
436 t=new AliITStrackV2(*esd);
437 } catch (const Char_t *msg) {
438 Warning("RefitInward",msg);
443 if (CorrectForDeadZoneMaterial(t)!=0) {
444 Warning("RefitInward",
445 "failed to correct for the material in the dead zone !\n");
450 ResetTrackToFollow(*t);
451 fTrackToFollow.ResetClusters();
454 if (RefitAt(3.7, &fTrackToFollow, t)) {
455 fTrackToFollow.SetLabel(t->GetLabel());
456 fTrackToFollow.CookdEdx();
457 CookLabel(&fTrackToFollow,0.); //For comparison only
459 if (fTrackToFollow.PropagateTo(3.,0.0028,65.19)) {//The beam pipe
460 Double_t a=fTrackToFollow.GetAlpha();
461 Double_t cs=TMath::Cos(a),sn=TMath::Sin(a);
462 Double_t xv= GetX()*cs + GetY()*sn;
463 Double_t yv=-GetX()*sn + GetY()*cs;
465 Double_t c=fTrackToFollow.GetC(), snp=fTrackToFollow.GetSnp();
466 Double_t x=fTrackToFollow.GetX(), y=fTrackToFollow.GetY();
467 Double_t tgfv=-(c*(x-xv)-snp)/(c*(y-yv) + TMath::Sqrt(1.-snp*snp));
468 Double_t fv=TMath::ATan(tgfv);
470 cs=TMath::Cos(fv); sn=TMath::Sin(fv);
472 yv=-xv*sn + yv*cs; xv=x;
474 if (fTrackToFollow.Propagate(fv+a,xv)) {
475 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
476 UseClusters(&fTrackToFollow);
478 AliITSclusterV2 c; c.SetY(yv); c.SetZ(GetZ());
479 c.SetSigmaY2(GetSigmaY()*GetSigmaY());
480 c.SetSigmaZ2(GetSigmaZ()*GetSigmaZ());
481 Double_t chi2=fTrackToFollow.GetPredictedChi2(&c);
483 if (fTrackToFollow.Update(&c,-chi2,0))
484 fTrackToFollow.SetConstrainedESDtrack(chi2);
493 Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
498 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
499 //--------------------------------------------------------------------
500 // Return pointer to a given cluster
501 //--------------------------------------------------------------------
502 Int_t l=(index & 0xf0000000) >> 28;
503 Int_t c=(index & 0x0fffffff) >> 00;
504 return fgLayers[l].GetCluster(c);
508 void AliITStrackerV2::FollowProlongation() {
509 //--------------------------------------------------------------------
510 //This function finds a track prolongation
511 //--------------------------------------------------------------------
512 while (fI>fLastLayerToTrackTo) {
515 AliITSlayer &layer=fgLayers[i];
516 AliITStrackV2 &track=fTracks[i];
518 Double_t r=layer.GetR();
521 Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
522 Double_t d=0.0034, x0=38.6;
523 if (i==1) {rs=9.; d=0.0097; x0=42;}
524 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
525 //Warning("FollowProlongation","propagation failed !\n");
532 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
533 //Warning("FollowProlongation","failed to estimate track !\n");
536 Double_t phi=TMath::ATan2(y,x);
538 Int_t idet=layer.FindDetectorIndex(phi,z);
540 //Warning("FollowProlongation","failed to find a detector !\n");
544 //propagate to the intersection
545 const AliITSdetector &det=layer.GetDetector(idet);
547 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
548 //Warning("FollowProlongation","propagation failed !\n");
551 fTrackToFollow.SetDetectorIndex(idet);
553 //Select possible prolongations and store the current track estimation
554 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
555 Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
556 Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
557 Double_t road=layer.GetRoad();
558 if (dz*dy>road*road) {
559 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
560 dz=road*scz; dy=road*scy;
563 //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
564 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
566 //Warning("FollowProlongation","too broad road in Z !\n");
570 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
572 //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
573 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
575 //Warning("FollowProlongation","too broad road in Y !\n");
579 Double_t zmin=track.GetZ() - dz;
580 Double_t zmax=track.GetZ() + dz;
581 Double_t ymin=track.GetY() + r*phi - dy;
582 Double_t ymax=track.GetY() + r*phi + dy;
583 layer.SelectClusters(zmin,zmax,ymin,ymax);
586 //take another prolongation
587 if (!TakeNextProlongation())
588 if (fLayersNotToSkip[fI]) return;
592 //deal with the best track
593 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
594 Int_t nclb=fBestTrack.GetNumberOfClusters();
597 Double_t chi2=fTrackToFollow.GetChi2();
598 if (chi2/ncl < kChi2PerCluster) {
599 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
607 Int_t AliITStrackerV2::TakeNextProlongation() {
608 //--------------------------------------------------------------------
609 // This function takes another track prolongation
611 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
612 //--------------------------------------------------------------------
613 AliITSlayer &layer=fgLayers[fI];
614 ResetTrackToFollow(fTracks[fI]);
616 Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
617 Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
618 Double_t road=layer.GetRoad();
619 if (dz*dy>road*road) {
620 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
621 dz=road*scz; dy=road*scy;
624 const AliITSclusterV2 *c=0; Int_t ci=-1;
625 Double_t chi2=12345.;
626 while ((c=layer.GetNextCluster(ci))!=0) {
627 Int_t idet=c->GetDetectorIndex();
629 if (fTrackToFollow.GetDetectorIndex()!=idet) {
630 const AliITSdetector &det=layer.GetDetector(idet);
631 ResetTrackToFollow(fTracks[fI]);
632 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
633 //Warning("TakeNextProlongation","propagation failed !\n");
636 fTrackToFollow.SetDetectorIndex(idet);
637 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
640 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
641 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
643 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
646 if (chi2>=kMaxChi2) return 0;
649 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
650 //Warning("TakeNextProlongation","filtering failed !\n");
654 if (fTrackToFollow.GetNumberOfClusters()>1)
655 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
658 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
662 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
663 fTrackToFollow.CorrectForMaterial(d,x0);
666 if (fConstraint[fPass]) {
667 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
668 Double_t xyz[]={GetX(),GetY(),GetZ()};
669 Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
670 fTrackToFollow.Improve(d,xyz,ers);
677 AliITStrackerV2::AliITSlayer::AliITSlayer() {
678 //--------------------------------------------------------------------
679 //default AliITSlayer constructor
680 //--------------------------------------------------------------------
685 AliITStrackerV2::AliITSlayer::
686 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
687 //--------------------------------------------------------------------
688 //main AliITSlayer constructor
689 //--------------------------------------------------------------------
690 fR=r; fPhiOffset=p; fZOffset=z;
691 fNladders=nl; fNdetectors=nd;
692 fDetectors=new AliITSdetector[fNladders*fNdetectors];
697 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
700 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
701 //--------------------------------------------------------------------
702 // AliITSlayer destructor
703 //--------------------------------------------------------------------
705 for (Int_t i=0; i<fN; i++) delete fClusters[i];
708 void AliITStrackerV2::AliITSlayer::ResetClusters() {
709 //--------------------------------------------------------------------
710 // This function removes loaded clusters
711 //--------------------------------------------------------------------
712 for (Int_t i=0; i<fN; i++) delete fClusters[i];
717 void AliITStrackerV2::AliITSlayer::ResetRoad() {
718 //--------------------------------------------------------------------
719 // This function calculates the road defined by the cluster density
720 //--------------------------------------------------------------------
722 for (Int_t i=0; i<fN; i++) {
723 if (TMath::Abs(fClusters[i]->GetZ())<fR) n++;
725 if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
728 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
729 //--------------------------------------------------------------------
730 //This function adds a cluster to this layer
731 //--------------------------------------------------------------------
732 if (fN==kMaxClusterPerLayer) {
733 ::Error("InsertCluster","Too many clusters !\n");
737 if (fN==0) {fClusters[fN++]=c; return 0;}
738 Int_t i=FindClusterIndex(c->GetZ());
739 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
740 fClusters[i]=c; fN++;
745 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
746 //--------------------------------------------------------------------
747 // This function returns the index of the nearest cluster
748 //--------------------------------------------------------------------
750 if (z <= fClusters[0]->GetZ()) return 0;
751 if (z > fClusters[fN-1]->GetZ()) return fN;
752 Int_t b=0, e=fN-1, m=(b+e)/2;
753 for (; b<e; m=(b+e)/2) {
754 if (z > fClusters[m]->GetZ()) b=m+1;
760 void AliITStrackerV2::AliITSlayer::
761 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
762 //--------------------------------------------------------------------
763 // This function sets the "window"
764 //--------------------------------------------------------------------
765 fI=FindClusterIndex(zmin); fZmax=zmax;
766 Double_t circle=2*TMath::Pi()*fR;
767 if (ymax>circle) { ymax-=circle; ymin-=circle; }
768 fYmin=ymin; fYmax=ymax;
771 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
772 //--------------------------------------------------------------------
773 // This function returns clusters within the "window"
774 //--------------------------------------------------------------------
775 const AliITSclusterV2 *cluster=0;
776 for (Int_t i=fI; i<fN; i++) {
777 const AliITSclusterV2 *c=fClusters[i];
778 if (c->GetZ() > fZmax) break;
779 if (c->IsUsed()) continue;
780 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
781 Double_t y=fR*det.GetPhi() + c->GetY();
783 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
784 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
786 if (y<fYmin) continue;
787 if (y>fYmax) continue;
796 Int_t AliITStrackerV2::AliITSlayer::
797 FindDetectorIndex(Double_t phi, Double_t z) const {
798 //--------------------------------------------------------------------
799 //This function finds the detector crossed by the track
800 //--------------------------------------------------------------------
801 Double_t dphi=-(phi-fPhiOffset);
802 if (dphi < 0) dphi += 2*TMath::Pi();
803 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
804 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
805 if (np>=fNladders) np-=fNladders;
806 if (np<0) np+=fNladders;
808 Double_t dz=fZOffset-z;
809 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
810 if (nz>=fNdetectors) return -1;
813 return np*fNdetectors + nz;
817 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
819 //--------------------------------------------------------------------
820 //This function returns the layer thickness at this point (units X0)
821 //--------------------------------------------------------------------
825 if (43<fR&&fR<45) { //SSD2
828 if (TMath::Abs(y-0.00)>3.40) d+=dd;
829 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
830 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
831 for (Int_t i=0; i<12; i++) {
832 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
833 if (TMath::Abs(y-0.00)>3.40) d+=dd;
837 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
838 if (TMath::Abs(y-0.00)>3.40) d+=dd;
842 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
843 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
846 if (37<fR&&fR<41) { //SSD1
849 if (TMath::Abs(y-0.00)>3.40) d+=dd;
850 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
851 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
852 for (Int_t i=0; i<11; i++) {
853 if (TMath::Abs(z-3.9*i)<0.15) {
854 if (TMath::Abs(y-0.00)>3.40) d+=dd;
858 if (TMath::Abs(z+3.9*i)<0.15) {
859 if (TMath::Abs(y-0.00)>3.40) d+=dd;
863 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
864 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
867 if (13<fR&&fR<26) { //SDD
870 if (TMath::Abs(y-0.00)>3.30) d+=dd;
872 if (TMath::Abs(y-1.80)<0.55) {
874 for (Int_t j=0; j<20; j++) {
875 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
876 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
879 if (TMath::Abs(y+1.80)<0.55) {
881 for (Int_t j=0; j<20; j++) {
882 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
883 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
887 for (Int_t i=0; i<4; i++) {
888 if (TMath::Abs(z-7.3*i)<0.60) {
890 if (TMath::Abs(y-0.00)>3.30) d+=dd;
893 if (TMath::Abs(z+7.3*i)<0.60) {
895 if (TMath::Abs(y-0.00)>3.30) d+=dd;
900 if (6<fR&&fR<8) { //SPD2
901 Double_t dd=0.0063; x0=21.5;
903 if (TMath::Abs(y-3.08)>0.5) d+=dd;
904 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
905 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
907 if (3<fR&&fR<5) { //SPD1
908 Double_t dd=0.0063; x0=21.5;
910 if (TMath::Abs(y+0.21)>0.6) d+=dd;
911 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
912 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
918 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
920 //--------------------------------------------------------------------
921 //Returns the thickness between the current layer and the vertex (units X0)
922 //--------------------------------------------------------------------
923 Double_t d=0.0028*3*3; //beam pipe
926 Double_t xn=fgLayers[fI].GetR();
927 for (Int_t i=0; i<fI; i++) {
928 Double_t xi=fgLayers[i].GetR();
929 d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
938 Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
945 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
946 //--------------------------------------------------------------------
947 // This function returns number of clusters within the "window"
948 //--------------------------------------------------------------------
950 for (Int_t i=fI; i<fN; i++) {
951 const AliITSclusterV2 *c=fClusters[i];
952 if (c->GetZ() > fZmax) break;
953 if (c->IsUsed()) continue;
954 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
955 Double_t y=fR*det.GetPhi() + c->GetY();
957 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
958 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
960 if (y<fYmin) continue;
961 if (y>fYmax) continue;
968 AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
969 //--------------------------------------------------------------------
970 // This function refits the track "t" at the position "x" using
971 // the clusters from "c"
972 //--------------------------------------------------------------------
973 Int_t index[kMaxLayer];
975 for (k=0; k<kMaxLayer; k++) index[k]=-1;
976 Int_t nc=c->GetNumberOfClusters();
977 for (k=0; k<nc; k++) {
978 Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
982 Int_t from, to, step;
983 if (xx > t->GetX()) {
984 from=0; to=kMaxLayer;
987 from=kMaxLayer-1; to=-1;
991 for (Int_t i=from; i != to; i += step) {
992 AliITSlayer &layer=fgLayers[i];
993 Double_t r=layer.GetR();
996 Double_t hI=i-0.5*step;
997 if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {
998 Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
999 Double_t d=0.0034, x0=38.6;
1000 if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1001 if (!t->PropagateTo(rs,-step*d,x0)) {
1007 // remember old position [SR, GSI 18.02.2003]
1008 Double_t oldX=0., oldY=0., oldZ=0.;
1009 if (t->IsStartedTimeIntegral() && step==1) {
1010 t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1015 if (!t->GetGlobalXYZat(r,x,y,z)) {
1018 Double_t phi=TMath::ATan2(y,x);
1019 Int_t idet=layer.FindDetectorIndex(phi,z);
1023 const AliITSdetector &det=layer.GetDetector(idet);
1025 if (!t->Propagate(phi,det.GetR())) {
1028 t->SetDetectorIndex(idet);
1030 const AliITSclusterV2 *cl=0;
1031 Double_t maxchi2=kMaxChi2;
1035 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
1036 if (idet != c->GetDetectorIndex()) {
1037 idet=c->GetDetectorIndex();
1038 const AliITSdetector &det=layer.GetDetector(idet);
1039 if (!t->Propagate(det.GetPhi(),det.GetR())) {
1042 t->SetDetectorIndex(idet);
1044 Double_t chi2=t->GetPredictedChi2(c);
1054 if (t->GetNumberOfClusters()>2) {
1055 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1056 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1057 Double_t zmin=t->GetZ() - dz;
1058 Double_t zmax=t->GetZ() + dz;
1059 Double_t ymin=t->GetY() + phi*r - dy;
1060 Double_t ymax=t->GetY() + phi*r + dy;
1061 layer.SelectClusters(zmin,zmax,ymin,ymax);
1063 const AliITSclusterV2 *c=0; Int_t ci=-1;
1064 while ((c=layer.GetNextCluster(ci))!=0) {
1065 if (idet != c->GetDetectorIndex()) continue;
1066 Double_t chi2=t->GetPredictedChi2(c);
1067 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1072 if (!t->Update(cl,maxchi2,idx)) {
1075 t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1080 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1081 t->CorrectForMaterial(-step*d,x0);
1084 // track time update [SR, GSI 17.02.2003]
1085 if (t->IsStartedTimeIntegral() && step==1) {
1086 Double_t newX, newY, newZ;
1087 t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1088 Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) +
1089 (oldZ-newZ)*(oldZ-newZ);
1090 t->AddTimeStep(TMath::Sqrt(dL2));
1096 if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1100 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1101 //--------------------------------------------------------------------
1102 // This function marks clusters assigned to the track
1103 //--------------------------------------------------------------------
1104 AliTracker::UseClusters(t,from);
1106 AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1107 //if (c->GetQ()>2) c->Use();
1108 if (c->GetSigmaZ2()>0.1) c->Use();
1109 c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1110 //if (c->GetQ()>2) c->Use();
1111 if (c->GetSigmaZ2()>0.1) c->Use();