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
19 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
20 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
21 //-------------------------------------------------------------------------
27 #include "AliITSgeom.h"
28 #include "AliITSRecPoint.h"
29 #include "AliTPCtrack.h"
30 #include "AliITSclusterV2.h"
31 #include "AliITStrackerV2.h"
33 ClassImp(AliITStrackerV2)
35 AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[kMaxLayer]; // ITS layers
37 AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
38 //--------------------------------------------------------------------
39 //This is the AliITStrackerV2 constructor
40 //--------------------------------------------------------------------
41 AliITSgeom *g=(AliITSgeom*)geom;
45 for (i=1; i<kMaxLayer+1; i++) {
46 Int_t nlad=g->GetNladders(i);
47 Int_t ndet=g->GetNdetectors(i);
49 g->GetTrans(i,1,1,x,y,z);
50 Double_t r=TMath::Sqrt(x*x + y*y);
51 Double_t poff=TMath::ATan2(y,x);
54 g->GetTrans(i,1,2,x,y,z);
55 r += TMath::Sqrt(x*x + y*y);
56 g->GetTrans(i,2,1,x,y,z);
57 r += TMath::Sqrt(x*x + y*y);
58 g->GetTrans(i,2,2,x,y,z);
59 r += TMath::Sqrt(x*x + y*y);
62 new (fLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
64 for (Int_t j=1; j<nlad+1; j++) {
65 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
66 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
67 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
69 Double_t r =-x*rot[1] + y*rot[0]; if (i==1) r=-r;
70 Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927;
71 phi+=0.5*TMath::Pi(); if (phi<0) phi += 2*TMath::Pi();
72 AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1);
74 new(&det) AliITSdetector(r,phi);
83 fConstraint[0]=1; fConstraint[1]=0;
85 Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV};
89 Int_t AliITStrackerV2::LoadClusters() {
90 //--------------------------------------------------------------------
91 //This function loads ITS clusters
92 //--------------------------------------------------------------------
94 sprintf(cname,"TreeC_ITS_%d",GetEventNumber());
95 TTree *cTree=(TTree*)gDirectory->Get(cname);
97 Error("LoadClusters"," can't get cTree !\n");
100 TBranch *branch=cTree->GetBranch("Clusters");
102 Error("LoadClusters"," can't get Clusters branch !\n");
106 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
107 branch->SetAddress(&clusters);
110 for (Int_t i=0; i<kMaxLayer; i++) {
111 Int_t ndet=fLayers[i].GetNdetectors();
112 Int_t jmax = j + fLayers[i].GetNladders()*ndet;
113 for (; j<jmax; j++) {
114 if (!cTree->GetEvent(j)) continue;
115 Int_t ncl=clusters->GetEntriesFast();
117 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
118 fLayers[i].InsertCluster(new AliITSclusterV2(*c));
122 fLayers[i].ResetRoad(); //road defined by the cluster density
124 delete cTree; //Thanks to Mariana Bondila
129 void AliITStrackerV2::UnloadClusters() {
130 //--------------------------------------------------------------------
131 //This function unloads ITS clusters
132 //--------------------------------------------------------------------
133 for (Int_t i=0; i<kMaxLayer; i++) fLayers[i].ResetClusters();
136 Int_t AliITStrackerV2::Clusters2Tracks(const TFile *inp, TFile *out) {
137 //--------------------------------------------------------------------
138 //This functions reconstructs ITS tracks
139 //--------------------------------------------------------------------
140 TFile *in=(TFile*)inp;
141 TDirectory *savedir=gDirectory;
143 if (LoadClusters()!=0) return 1;
146 Error("Clusters2Tracks","file with TPC tracks is not open !\n");
150 if (!out->IsOpen()) {
151 Error("Clusters2Tracks","file for ITS tracks is not open !\n");
158 Int_t nentr=0; TObjArray itsTracks(15000);
160 {/* Read TPC tracks */
161 sprintf(tname,"TreeT_TPC_%d",GetEventNumber());
162 TTree *tpcTree=(TTree*)in->Get(tname);
164 Error("Clusters2Tracks","can't get a tree with TPC tracks !\n");
167 AliTPCtrack *itrack=new AliTPCtrack;
168 tpcTree->SetBranchAddress("tracks",&itrack);
169 nentr=(Int_t)tpcTree->GetEntries();
170 for (Int_t i=0; i<nentr; i++) {
171 tpcTree->GetEvent(i);
174 t=new AliITStrackV2(*itrack);
175 } catch (const Char_t *msg) {
176 Warning("Clusters2Tracks",msg);
180 if (TMath::Abs(t->GetD())>4) continue;
182 t->PropagateTo(80.,0.0053,30);
183 if (TMath::Abs(t->GetY())>12.8) t->CorrectForMaterial(0.03);
184 if (TMath::Abs(t->GetZ())<0.2) t->CorrectForMaterial(0.40);
185 t->PropagateTo(61.,0.0053,30);
186 //Double_t xk=52.,x,y,z; t->GetGlobalXYZat(xk,x,y,z);
187 //if (TMath::Abs(y)<7.77) t->PropagateTo(xk,0.19,24.);
188 t->PropagateTo(50.,0.001);
190 itsTracks.AddLast(t);
192 delete tpcTree; //Thanks to Mariana Bondila
199 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
200 TTree itsTree(tname,"Tree with ITS tracks");
201 AliITStrackV2 *otrack=&fBestTrack;
202 itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
204 for (fPass=0; fPass<2; fPass++) {
205 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
206 for (Int_t i=0; i<nentr; i++) {
207 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
208 if (t==0) continue; //this track has been already tracked
209 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
211 ResetTrackToFollow(*t);
214 for (FollowProlongation(); fI<kMaxLayer; fI++) {
215 while (TakeNextProlongation()) FollowProlongation();
218 if (fBestTrack.GetNumberOfClusters() < kMaxLayer-kLayersToSkip)continue;
220 if (fConstraint[fPass]) {
221 Int_t index[kMaxLayer];
223 for (k=0; k<kMaxLayer; k++) index[k]=-1;
224 Int_t nc=fBestTrack.GetNumberOfClusters();
225 for (k=0; k<nc; k++) {
226 Int_t idx=fBestTrack.GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
229 fBestTrack.~AliITStrackV2(); new(&fBestTrack) AliITStrackV2(*t);
230 if (!RefitAt(3.7, &fBestTrack, index)) continue;
233 fBestTrack.SetLabel(tpcLabel);
234 fBestTrack.CookdEdx();
235 CookLabel(&fBestTrack,0.); //For comparison only
237 UseClusters(&fBestTrack);
238 delete itsTracks.RemoveAt(i);
243 Info("Clusters2Tracks","Number of TPC tracks: %d\n",nentr);
244 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",
245 (Int_t)itsTree.GetEntries());
257 Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
258 //--------------------------------------------------------------------
259 //This functions propagates reconstructed ITS tracks back
260 //--------------------------------------------------------------------
261 TFile *in=(TFile*)inp;
262 TDirectory *savedir=gDirectory;
264 if (LoadClusters()!=0) return 1;
267 Error("PropagateBack","file with ITS tracks is not open !\n");
271 if (!out->IsOpen()) {
272 Error("PropagateBack","file for back propagated ITS tracks is not open !\n");
279 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
280 TTree *itsTree=(TTree*)in->Get(tname);
282 Error("PropagateBack","can't get a tree with ITS tracks !\n");
285 AliITStrackV2 *itrack=new AliITStrackV2;
286 itsTree->SetBranchAddress("tracks",&itrack);
290 sprintf(tname,"TreeT_ITSb_%d",GetEventNumber());
291 TTree backTree(tname,"Tree with back propagated ITS tracks");
292 AliTPCtrack *otrack=0;
293 backTree.Branch("tracks","AliTPCtrack",&otrack,32000,2);
297 Int_t nentr=(Int_t)itsTree->GetEntries();
298 for (Int_t i=0; i<nentr; i++) {
299 itsTree->GetEvent(i);
300 ResetTrackToFollow(*itrack);
302 // propagete to vertex [SR, GSI 17.02.2003]
303 fTrackToFollow.PropagateTo(3.,0.0028,65.19);
304 fTrackToFollow.PropagateToVertex();
306 // Start Time measurement [SR, GSI 17.02.2003]
307 fTrackToFollow.StartTimeIntegral();
308 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
311 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
313 Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
316 Int_t nc=itrack->GetNumberOfClusters();
318 const AliITSclusterV2 *c=0;
320 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
321 c=(AliITSclusterV2*)GetCluster(idx);
323 for (fI=0; fI<kMaxLayer; fI++) {
324 AliITSlayer &layer=fLayers[fI];
325 Double_t r=layer.GetR();
326 if (fI==2 || fI==4) {
327 Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
328 Double_t d=0.0034, x0=38.6;
329 if (fI==2) {rs=9.; d=0.0097; x0=42.;}
330 if (!fTrackToFollow.PropagateTo(rs,-d,x0)) throw "";
333 // remember old position [SR, GSI 18.02.2003]
334 Double_t oldX, oldY, oldZ;
335 fTrackToFollow.GetGlobalXYZat(fTrackToFollow.GetX(),oldX,oldY,oldZ);
339 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z))
340 throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
341 Double_t phi=TMath::ATan2(y,x);
342 Int_t idet=layer.FindDetectorIndex(phi,z);
344 throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
345 const AliITSdetector &det=layer.GetDetector(idet);
346 r=det.GetR(); phi=det.GetPhi();
347 if (!fTrackToFollow.Propagate(phi,r)) throw "";
348 fTrackToFollow.SetDetectorIndex(idet);
350 const AliITSclusterV2 *cl=0;
352 Double_t maxchi2=kMaxChi2;
355 idet=c->GetDetectorIndex();
356 if (idet != fTrackToFollow.GetDetectorIndex()) {
357 const AliITSdetector &det=layer.GetDetector(idet);
358 r=det.GetR(); phi=det.GetPhi();
359 if (!fTrackToFollow.Propagate(phi,r)) throw "";
360 fTrackToFollow.SetDetectorIndex(idet);
362 Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
364 cl=c; maxchi2=chi2; index=idx;
367 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
368 c=(AliITSclusterV2*)GetCluster(idx);
372 if (fTrackToFollow.GetNumberOfClusters()>2) {
373 Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
374 Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
375 Double_t zmin=fTrackToFollow.GetZ() - dz;
376 Double_t zmax=fTrackToFollow.GetZ() + dz;
377 Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
378 Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
379 layer.SelectClusters(zmin,zmax,ymin,ymax);
381 const AliITSclusterV2 *cc=0; Int_t ci;
382 while ((cc=layer.GetNextCluster(ci))!=0) {
383 idet=cc->GetDetectorIndex();
384 if (idet != fTrackToFollow.GetDetectorIndex()) continue;
385 Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
387 cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
393 if (!fTrackToFollow.Update(cl,maxchi2,index))
394 Info("PropagateBack","filtering failed !\n");
398 x=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
399 fTrackToFollow.CorrectForMaterial(-x,x0);
402 // track time update [SR, GSI 17.02.2003]
403 Double_t newX, newY, newZ;
404 fTrackToFollow.GetGlobalXYZat(fTrackToFollow.GetX(),newX,newY,newZ);
405 Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) +
406 (oldZ-newZ)*(oldZ-newZ);
407 fTrackToFollow.AddTimeStep(TMath::Sqrt(dL2));
412 fTrackToFollow.PropagateTo(50.,-0.001);
413 //Double_t xk=52.,x,y,z; fTrackToFollow.GetGlobalXYZat(xk,x,y,z);
414 //if (TMath::Abs(y)<7.77) fTrackToFollow.PropagateTo(xk,-0.19,24.);
415 fTrackToFollow.PropagateTo(61,-0.0053,30);
416 fTrackToFollow.PropagateTo(80.,-0.0053,30);
418 fTrackToFollow.SetLabel(itsLabel);
419 otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha());
420 backTree.Fill(); delete otrack;
421 UseClusters(&fTrackToFollow);
424 catch (const Char_t *msg) {
425 Warning("PropagateBack",msg);
431 Info("PropagateBack","Number of ITS tracks: %d\n",nentr);
432 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
436 delete itsTree; //Thanks to Mariana Bondila
445 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
446 //--------------------------------------------------------------------
447 // Return pointer to a given cluster
448 //--------------------------------------------------------------------
449 Int_t l=(index & 0xf0000000) >> 28;
450 Int_t c=(index & 0x0fffffff) >> 00;
451 return fLayers[l].GetCluster(c);
455 void AliITStrackerV2::FollowProlongation() {
456 //--------------------------------------------------------------------
457 //This function finds a track prolongation
458 //--------------------------------------------------------------------
459 Int_t tryAgain=kLayersToSkip;
464 AliITSlayer &layer=fLayers[i];
465 AliITStrackV2 &track=fTracks[i];
467 Double_t r=layer.GetR();
470 Double_t rs=0.5*(fLayers[i+1].GetR() + r);
471 Double_t d=0.0034, x0=38.6;
472 if (i==1) {rs=9.; d=0.0097; x0=42;}
473 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
474 //Warning("FollowProlongation","propagation failed !\n");
481 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
482 //Warning("FollowProlongation","failed to estimate track !\n");
485 Double_t phi=TMath::ATan2(y,x);
486 Int_t idet=layer.FindDetectorIndex(phi,z);
488 //Warning("FollowProlongation","failed to find a detector !\n");
492 //propagate to the intersection
493 const AliITSdetector &det=layer.GetDetector(idet);
495 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
496 //Warning("FollowProlongation","propagation failed !\n");
499 fTrackToFollow.SetDetectorIndex(idet);
501 //Select possible prolongations and store the current track estimation
502 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
503 Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
504 Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
505 Double_t road=layer.GetRoad();
506 if (dz*dy>road*road) {
507 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
508 dz=road*scz; dy=road*scy;
511 //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
512 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
514 //Warning("FollowProlongation","too broad road in Z !\n");
518 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) break;
520 //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
521 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
523 //Warning("FollowProlongation","too broad road in Y !\n");
527 Double_t zmin=track.GetZ() - dz;
528 Double_t zmax=track.GetZ() + dz;
529 Double_t ymin=track.GetY() + r*phi - dy;
530 Double_t ymax=track.GetY() + r*phi + dy;
531 layer.SelectClusters(zmin,zmax,ymin,ymax);
534 //take another prolongation
535 if (!TakeNextProlongation()) if (!tryAgain--) break;
536 tryAgain=kLayersToSkip;
540 //deal with the best track
541 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
542 Int_t nclb=fBestTrack.GetNumberOfClusters();
545 Double_t chi2=fTrackToFollow.GetChi2();
546 if (chi2/ncl < kChi2PerCluster) {
547 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
555 Int_t AliITStrackerV2::TakeNextProlongation() {
556 //--------------------------------------------------------------------
557 // This function takes another track prolongation
559 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
560 //--------------------------------------------------------------------
561 AliITSlayer &layer=fLayers[fI];
562 ResetTrackToFollow(fTracks[fI]);
564 Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
565 Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
566 Double_t road=layer.GetRoad();
567 if (dz*dy>road*road) {
568 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
569 dz=road*scz; dy=road*scy;
572 const AliITSclusterV2 *c=0; Int_t ci=-1;
573 Double_t chi2=12345.;
574 while ((c=layer.GetNextCluster(ci))!=0) {
575 Int_t idet=c->GetDetectorIndex();
577 if (fTrackToFollow.GetDetectorIndex()!=idet) {
578 const AliITSdetector &det=layer.GetDetector(idet);
579 ResetTrackToFollow(fTracks[fI]);
580 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
581 //Warning("TakeNextProlongation","propagation failed !\n");
584 fTrackToFollow.SetDetectorIndex(idet);
585 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
588 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
589 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
591 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
594 if (chi2>=kMaxChi2) return 0;
597 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
598 //Warning("TakeNextProlongation","filtering failed !\n");
602 if (fTrackToFollow.GetNumberOfClusters()>1)
603 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
606 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
610 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
611 fTrackToFollow.CorrectForMaterial(d,x0);
614 if (fConstraint[fPass]) {
615 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
616 Double_t xyz[]={GetX(),GetY(),GetZ()};
617 Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
618 fTrackToFollow.Improve(d,xyz,ers);
625 AliITStrackerV2::AliITSlayer::AliITSlayer() {
626 //--------------------------------------------------------------------
627 //default AliITSlayer constructor
628 //--------------------------------------------------------------------
633 AliITStrackerV2::AliITSlayer::
634 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
635 //--------------------------------------------------------------------
636 //main AliITSlayer constructor
637 //--------------------------------------------------------------------
638 fR=r; fPhiOffset=p; fZOffset=z;
639 fNladders=nl; fNdetectors=nd;
640 fDetectors=new AliITSdetector[fNladders*fNdetectors];
645 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
648 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
649 //--------------------------------------------------------------------
650 // AliITSlayer destructor
651 //--------------------------------------------------------------------
653 for (Int_t i=0; i<fN; i++) delete fClusters[i];
656 void AliITStrackerV2::AliITSlayer::ResetClusters() {
657 //--------------------------------------------------------------------
658 // This function removes loaded clusters
659 //--------------------------------------------------------------------
660 for (Int_t i=0; i<fN; i++) delete fClusters[i];
665 void AliITStrackerV2::AliITSlayer::ResetRoad() {
666 //--------------------------------------------------------------------
667 // This function calculates the road defined by the cluster density
668 //--------------------------------------------------------------------
670 for (Int_t i=0; i<fN; i++) {
671 if (TMath::Abs(fClusters[i]->GetZ())<fR) n++;
673 if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
676 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
677 //--------------------------------------------------------------------
678 //This function adds a cluster to this layer
679 //--------------------------------------------------------------------
680 if (fN==kMaxClusterPerLayer) {
681 ::Error("InsertCluster","Too many clusters !\n");
685 if (fN==0) {fClusters[fN++]=c; return 0;}
686 Int_t i=FindClusterIndex(c->GetZ());
687 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
688 fClusters[i]=c; fN++;
693 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
694 //--------------------------------------------------------------------
695 // This function returns the index of the nearest cluster
696 //--------------------------------------------------------------------
698 if (z <= fClusters[0]->GetZ()) return 0;
699 if (z > fClusters[fN-1]->GetZ()) return fN;
700 Int_t b=0, e=fN-1, m=(b+e)/2;
701 for (; b<e; m=(b+e)/2) {
702 if (z > fClusters[m]->GetZ()) b=m+1;
708 void AliITStrackerV2::AliITSlayer::
709 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
710 //--------------------------------------------------------------------
711 // This function sets the "window"
712 //--------------------------------------------------------------------
713 fI=FindClusterIndex(zmin); fZmax=zmax;
714 Double_t circle=2*TMath::Pi()*fR;
715 if (ymax>circle) { ymax-=circle; ymin-=circle; }
716 fYmin=ymin; fYmax=ymax;
719 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
720 //--------------------------------------------------------------------
721 // This function returns clusters within the "window"
722 //--------------------------------------------------------------------
723 const AliITSclusterV2 *cluster=0;
724 for (Int_t i=fI; i<fN; i++) {
725 const AliITSclusterV2 *c=fClusters[i];
726 if (c->GetZ() > fZmax) break;
727 if (c->IsUsed()) continue;
728 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
729 Double_t y=fR*det.GetPhi() + c->GetY();
731 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
732 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
734 if (y<fYmin) continue;
735 if (y>fYmax) continue;
744 Int_t AliITStrackerV2::AliITSlayer::
745 FindDetectorIndex(Double_t phi, Double_t z) const {
746 //--------------------------------------------------------------------
747 //This function finds the detector crossed by the track
748 //--------------------------------------------------------------------
749 Double_t dphi=phi-fPhiOffset;
750 if (dphi < 0) dphi += 2*TMath::Pi();
751 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
752 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
753 if (np>=fNladders) np-=fNladders;
754 if (np<0) np+=fNladders;
756 Double_t dz=fZOffset-z;
757 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
758 if (nz>=fNdetectors) return -1;
761 return np*fNdetectors + nz;
765 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
767 //--------------------------------------------------------------------
768 //This function returns the layer thickness at this point (units X0)
769 //--------------------------------------------------------------------
773 if (43<fR&&fR<45) { //SSD2
776 if (TMath::Abs(y-0.00)>3.40) d+=dd;
777 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
778 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
779 for (Int_t i=0; i<12; i++) {
780 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
781 if (TMath::Abs(y-0.00)>3.40) d+=dd;
785 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
786 if (TMath::Abs(y-0.00)>3.40) d+=dd;
790 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
791 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
794 if (37<fR&&fR<41) { //SSD1
797 if (TMath::Abs(y-0.00)>3.40) d+=dd;
798 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
799 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
800 for (Int_t i=0; i<11; i++) {
801 if (TMath::Abs(z-3.9*i)<0.15) {
802 if (TMath::Abs(y-0.00)>3.40) d+=dd;
806 if (TMath::Abs(z+3.9*i)<0.15) {
807 if (TMath::Abs(y-0.00)>3.40) d+=dd;
811 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
812 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
815 if (13<fR&&fR<26) { //SDD
818 if (TMath::Abs(y-0.00)>3.30) d+=dd;
820 if (TMath::Abs(y-1.80)<0.55) {
822 for (Int_t j=0; j<20; j++) {
823 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
824 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
827 if (TMath::Abs(y+1.80)<0.55) {
829 for (Int_t j=0; j<20; j++) {
830 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
831 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
835 for (Int_t i=0; i<4; i++) {
836 if (TMath::Abs(z-7.3*i)<0.60) {
838 if (TMath::Abs(y-0.00)>3.30) d+=dd;
841 if (TMath::Abs(z+7.3*i)<0.60) {
843 if (TMath::Abs(y-0.00)>3.30) d+=dd;
848 if (6<fR&&fR<8) { //SPD2
849 Double_t dd=0.0063; x0=21.5;
851 if (TMath::Abs(y-3.08)>0.5) d+=dd;
852 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
853 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
855 if (3<fR&&fR<5) { //SPD1
856 Double_t dd=0.0063; x0=21.5;
858 if (TMath::Abs(y+0.21)>0.6) d+=dd;
859 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
860 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
866 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
868 //--------------------------------------------------------------------
869 //Returns the thickness between the current layer and the vertex (units X0)
870 //--------------------------------------------------------------------
871 Double_t d=0.0028*3*3; //beam pipe
874 Double_t xn=fLayers[fI].GetR();
875 for (Int_t i=0; i<fI; i++) {
876 Double_t xi=fLayers[i].GetR();
877 d+=fLayers[i].GetThickness(y,z,x0)*xi*xi;
886 Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
893 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
894 //--------------------------------------------------------------------
895 // This function returns number of clusters within the "window"
896 //--------------------------------------------------------------------
898 for (Int_t i=fI; i<fN; i++) {
899 const AliITSclusterV2 *c=fClusters[i];
900 if (c->GetZ() > fZmax) break;
901 if (c->IsUsed()) continue;
902 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
903 Double_t y=fR*det.GetPhi() + c->GetY();
905 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
906 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
908 if (y<fYmin) continue;
909 if (y>fYmax) continue;
915 Bool_t AliITStrackerV2::RefitAt(Double_t x, AliITStrackV2 *t, Int_t *index) {
916 //--------------------------------------------------------------------
917 // This function refits a track at a given position
918 //--------------------------------------------------------------------
919 Int_t from, to, step;
921 from=0; to=kMaxLayer;
924 from=kMaxLayer-1; to=-1;
928 for (Int_t i=from; i != to; i += step) {
929 AliITSlayer &layer=fLayers[i];
930 Double_t r=layer.GetR();
933 Double_t hI=i-0.5*step;
934 if (hI==1.5 || hI==3.5) {
935 Double_t rs=0.5*(fLayers[i-step].GetR() + r);
936 Double_t d=0.0034, x0=38.6;
937 if (hI==1.5) {rs=9.; d=0.0097; x0=42;}
938 if (!t->PropagateTo(rs,d,x0)) {
945 if (!t->GetGlobalXYZat(r,x,y,z)) {
948 Double_t phi=TMath::ATan2(y,x);
949 Int_t idet=layer.FindDetectorIndex(phi,z);
953 const AliITSdetector &det=layer.GetDetector(idet);
955 if (!t->Propagate(phi,det.GetR())) {
958 t->SetDetectorIndex(idet);
960 const AliITSclusterV2 *cl=0;
961 Double_t maxchi2=kMaxChi2;
965 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
966 if (idet != c->GetDetectorIndex()) {
967 idet=c->GetDetectorIndex();
968 const AliITSdetector &det=layer.GetDetector(idet);
969 if (!t->Propagate(det.GetPhi(),det.GetR())) {
972 t->SetDetectorIndex(idet);
974 Double_t chi2=t->GetPredictedChi2(c);
975 if (chi2<maxchi2) { cl=c; maxchi2=chi2; }
981 if (t->GetNumberOfClusters()>2) {
982 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
983 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
984 Double_t zmin=t->GetZ() - dz;
985 Double_t zmax=t->GetZ() + dz;
986 Double_t ymin=t->GetY() + phi*r - dy;
987 Double_t ymax=t->GetY() + phi*r + dy;
988 layer.SelectClusters(zmin,zmax,ymin,ymax);
990 const AliITSclusterV2 *c=0; Int_t ci=-1;
991 while ((c=layer.GetNextCluster(ci))!=0) {
992 if (idet != c->GetDetectorIndex()) continue;
993 Double_t chi2=t->GetPredictedChi2(c);
994 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1000 if (!t->Update(cl,maxchi2,idx)) {
1007 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1008 t->CorrectForMaterial(-step*d,x0);
1013 if (!t->PropagateTo(x,0.,0.)) return kFALSE;
1017 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1018 //--------------------------------------------------------------------
1019 // This function marks clusters assigned to the track
1020 //--------------------------------------------------------------------
1021 AliTracker::UseClusters(t,from);
1023 AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1024 //if (c->GetQ()>2) c->Use();
1025 if (c->GetSigmaZ2()>0.1) c->Use();
1026 c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1027 //if (c->GetQ()>2) c->Use();
1028 if (c->GetSigmaZ2()>0.1) c->Use();