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 //-------------------------------------------------------------------------
25 #include <Riostream.h>
27 #include "AliITSgeom.h"
28 #include "AliITSRecPoint.h"
29 #include "AliTPCtrack.h"
30 #include "AliITSclusterV2.h"
31 #include "AliITStrackerV2.h"
39 AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[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 (fLayers+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 r =-x*rot[1] + y*rot[0]; if (i==1) r=-r;
74 Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927;
75 phi+=0.5*TMath::Pi(); if (phi<0) phi += 2*TMath::Pi();
76 AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1);
78 new(&det) AliITSdetector(r,phi);
87 fConstraint[0]=1; fConstraint[1]=0;
90 void AliITStrackerV2::LoadClusters() {
91 //--------------------------------------------------------------------
92 //This function loads ITS clusters
93 //--------------------------------------------------------------------
95 sprintf(cname,"TreeC_ITS_%d",GetEventNumber());
96 TTree *cTree=(TTree*)gDirectory->Get(cname);
98 cerr<<"AliITStrackerV2::LoadClusters can't get cTree !\n";
101 TBranch *branch=cTree->GetBranch("Clusters");
103 cerr<<"AliITStrackerV2::LoadClusters can't get Clusters branch !\n";
107 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
108 branch->SetAddress(&clusters);
111 for (Int_t i=0; i<kMaxLayer; i++) {
112 Int_t jmax = j + fLayers[i].GetNladders()*fLayers[i].GetNdetectors();
113 for (; j<jmax; j++) {
114 if (!cTree->GetEvent(j)) continue;
115 Int_t ncl=clusters->GetEntriesFast();
117 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
120 if (c->GetLabel(2)!=LAB)
121 if (c->GetLabel(1)!=LAB)
122 if (c->GetLabel(0)!=LAB) continue;
123 cout<<lay-1<<' '<<lad-1<<' '<<det-1<<' '<<c->GetY()<<' '<<c->GetZ()<<endl;
126 fLayers[i].InsertCluster(new AliITSclusterV2(*c));
131 delete cTree; //Thanks to Mariana Bondila
134 void AliITStrackerV2::UnloadClusters() {
135 //--------------------------------------------------------------------
136 //This function unloads ITS clusters
137 //--------------------------------------------------------------------
138 for (Int_t i=0; i<kMaxLayer; i++) fLayers[i].ResetClusters();
145 Int_t AliITStrackerV2::Clusters2Tracks(const TFile *inp, TFile *out) {
146 //--------------------------------------------------------------------
147 //This functions reconstructs ITS tracks
148 //--------------------------------------------------------------------
149 TFile *in=(TFile*)inp;
150 TDirectory *savedir=gDirectory;
155 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
156 cerr<<"file with TPC tracks is not open !\n";
160 if (!out->IsOpen()) {
161 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
162 cerr<<"file for ITS tracks is not open !\n";
169 Int_t nentr=0; TObjArray itsTracks(15000);
171 {/* Read TPC tracks */
172 sprintf(tname,"TreeT_TPC_%d",GetEventNumber());
173 TTree *tpcTree=(TTree*)in->Get(tname);
175 cerr<<"AliITStrackerV2::Clusters2Tracks(): "
176 "can't get a tree with TPC tracks !\n";
179 AliTPCtrack *itrack=new AliTPCtrack;
180 tpcTree->SetBranchAddress("tracks",&itrack);
181 nentr=(Int_t)tpcTree->GetEntries();
182 for (Int_t i=0; i<nentr; i++) {
183 tpcTree->GetEvent(i);
186 t=new AliITStrackV2(*itrack);
187 } catch (const Char_t *msg) {
192 if (TMath::Abs(t->GetD())>4) continue;
194 t->PropagateTo(80.,0.0053,30);
195 if (TMath::Abs(t->GetY())>12.8) t->CorrectForMaterial(0.03);
196 if (TMath::Abs(t->GetZ())<0.2) t->CorrectForMaterial(0.40);
197 t->PropagateTo(61.,0.0053,30);
198 //Double_t xk=52.,x,y,z; t->GetGlobalXYZat(xk,x,y,z);
199 //if (TMath::Abs(y)<7.77) t->PropagateTo(xk,0.19,24.);
200 t->PropagateTo(50.,0.001);
202 itsTracks.AddLast(t);
204 delete tpcTree; //Thanks to Mariana Bondila
211 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
212 TTree itsTree(tname,"Tree with ITS tracks");
213 AliITStrackV2 *otrack=&fBestTrack;
214 itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
216 for (fPass=0; fPass<2; fPass++) {
217 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
218 for (Int_t i=0; i<nentr; i++) {
219 if (i%10==0) cerr<<nentr-i<<" \r";
220 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
221 if (t==0) continue; //this track has been already tracked
222 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
227 if (TMath::Abs(tpcLabel)!=LAB) continue;
228 cout<<tpcLabel<<" *****************\n";
231 ResetTrackToFollow(*t);
234 for (FollowProlongation(); fI<kMaxLayer; fI++) {
235 while (TakeNextProlongation()) FollowProlongation();
239 cout<<fBestTrack.GetNumberOfClusters()<<" number of clusters\n\n";
242 if (fBestTrack.GetNumberOfClusters() < kMaxLayer-kLayersToSkip) continue;
244 if (fConstraint[fPass]) {
245 Int_t index[kMaxLayer];
247 for (k=0; k<kMaxLayer; k++) index[k]=-1;
248 Int_t nc=fBestTrack.GetNumberOfClusters();
249 for (k=0; k<nc; k++) {
250 Int_t idx=fBestTrack.GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
253 fBestTrack.~AliITStrackV2(); new(&fBestTrack) AliITStrackV2(*t);
254 if (!RefitAt(3.7, &fBestTrack, index)) continue;
257 fBestTrack.SetLabel(tpcLabel);
258 fBestTrack.CookdEdx();
259 CookLabel(&fBestTrack,0.); //For comparison only
261 UseClusters(&fBestTrack);
262 delete itsTracks.RemoveAt(i);
274 cerr<<"Number of TPC tracks: "<<nentr<<endl;
275 cerr<<"Number of prolonged tracks: "<<itsTree.GetEntries()<<endl;
280 Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
281 //--------------------------------------------------------------------
282 //This functions propagates reconstructed ITS tracks back
283 //--------------------------------------------------------------------
284 TFile *in=(TFile*)inp;
285 TDirectory *savedir=gDirectory;
290 cerr<<"AliITStrackerV2::PropagateBack(): ";
291 cerr<<"file with ITS tracks is not open !\n";
295 if (!out->IsOpen()) {
296 cerr<<"AliITStrackerV2::PropagateBack(): ";
297 cerr<<"file for back propagated ITS tracks is not open !\n";
304 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
305 TTree *itsTree=(TTree*)in->Get(tname);
307 cerr<<"AliITStrackerV2::PropagateBack() ";
308 cerr<<"can't get a tree with ITS tracks !\n";
311 AliITStrackV2 *itrack=new AliITStrackV2;
312 itsTree->SetBranchAddress("tracks",&itrack);
316 sprintf(tname,"TreeT_ITSb_%d",GetEventNumber());
317 TTree backTree(tname,"Tree with back propagated ITS tracks");
318 AliTPCtrack *otrack=0;
319 backTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
323 Int_t nentr=(Int_t)itsTree->GetEntries();
324 for (Int_t i=0; i<nentr; i++) {
325 itsTree->GetEvent(i);
326 ResetTrackToFollow(*itrack);
327 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
328 Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
331 if (TMath::Abs(itsLabel)!=LAB) continue;
332 cout<<itsLabel<<" *****************\n";
336 Int_t nc=itrack->GetNumberOfClusters();
338 for (Int_t k=0; k<nc; k++) {
339 Int_t index=itrack->GetClusterIndex(k);
340 AliITSclusterV2 *c=(AliITSclusterV2*)GetCluster(index);
341 cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
345 const AliITSclusterV2 *c=0;
347 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
348 c=(AliITSclusterV2*)GetCluster(idx);
350 for (fI=0; fI<kMaxLayer; fI++) {
351 AliITSlayer &layer=fLayers[fI];
352 Double_t r=layer.GetR();
353 if (fI==2 || fI==4) {
354 Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
355 Double_t d=0.0034, x0=38.6;
356 if (fI==2) {rs=9.; d=0.0097; x0=42.;}
357 if (!fTrackToFollow.PropagateTo(rs,-d,x0)) throw "";
361 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z))
362 throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
363 Double_t phi=TMath::ATan2(y,x);
364 Int_t idet=layer.FindDetectorIndex(phi,z);
366 throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
367 const AliITSdetector &det=layer.GetDetector(idet);
368 r=det.GetR(); phi=det.GetPhi();
369 if (!fTrackToFollow.Propagate(phi,r)) throw "";
370 fTrackToFollow.SetDetectorIndex(idet);
372 const AliITSclusterV2 *cl=0;
374 Double_t maxchi2=kMaxChi2;
377 idet=c->GetDetectorIndex();
378 if (idet != fTrackToFollow.GetDetectorIndex()) {
379 const AliITSdetector &det=layer.GetDetector(idet);
380 r=det.GetR(); phi=det.GetPhi();
381 if (!fTrackToFollow.Propagate(phi,r)) throw "";
382 fTrackToFollow.SetDetectorIndex(idet);
384 Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
386 cl=c; maxchi2=chi2; index=idx;
389 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
390 c=(AliITSclusterV2*)GetCluster(idx);
394 if (fTrackToFollow.GetNumberOfClusters()>2) {
395 Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
396 Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
397 Double_t zmin=fTrackToFollow.GetZ() - dz;
398 Double_t zmax=fTrackToFollow.GetZ() + dz;
399 Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
400 Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
401 layer.SelectClusters(zmin,zmax,ymin,ymax);
403 const AliITSclusterV2 *cc=0; Int_t ci;
404 while ((cc=layer.GetNextCluster(ci))!=0) {
405 idet=cc->GetDetectorIndex();
406 if (idet != fTrackToFollow.GetDetectorIndex()) continue;
407 Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
409 cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
415 if (!fTrackToFollow.Update(cl,maxchi2,index))
416 cerr<<"AliITStrackerV2::PropagateBack: filtering failed !\n";
420 x=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
421 fTrackToFollow.CorrectForMaterial(-x,x0);
425 fTrackToFollow.PropagateTo(50.,-0.001);
426 //Double_t xk=52.,x,y,z; fTrackToFollow.GetGlobalXYZat(xk,x,y,z);
427 //if (TMath::Abs(y)<7.77) fTrackToFollow.PropagateTo(xk,-0.19,24.);
428 fTrackToFollow.PropagateTo(61,-0.0053,30);
429 fTrackToFollow.PropagateTo(80.,-0.0053,30);
431 fTrackToFollow.SetLabel(itsLabel);
432 otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha());
433 backTree.Fill(); delete otrack;
434 UseClusters(&fTrackToFollow);
435 cerr << ++ntrk << " \r";
437 catch (const Char_t *msg) {
444 cerr<<"Number of ITS tracks: "<<nentr<<endl;
445 cerr<<"Number of back propagated ITS tracks: "<<ntrk<<endl;
449 delete itsTree; //Thanks to Mariana Bondila
458 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
459 //--------------------------------------------------------------------
460 // Return pointer to a given cluster
461 //--------------------------------------------------------------------
462 Int_t l=(index & 0xf0000000) >> 28;
463 Int_t c=(index & 0x0fffffff) >> 00;
464 return fLayers[l].GetCluster(c);
468 void AliITStrackerV2::FollowProlongation() {
469 //--------------------------------------------------------------------
470 //This function finds a track prolongation
471 //--------------------------------------------------------------------
472 Int_t tryAgain=kLayersToSkip;
479 AliITSlayer &layer=fLayers[i];
480 AliITStrackV2 &track=fTracks[i];
482 Double_t r=layer.GetR();
485 Double_t rs=0.5*(fLayers[i+1].GetR() + r);
486 Double_t d=0.0034, x0=38.6;
487 if (i==1) {rs=9.; d=0.0097; x0=42;}
488 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
489 //cerr<<"AliITStrackerV2::FollowProlongation: "
490 //"propagation failed !\n";
497 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
498 //cerr<<"AliITStrackerV2::FollowProlongation: "
499 //"failed to estimate track !\n";
502 Double_t phi=TMath::ATan2(y,x);
503 Int_t idet=layer.FindDetectorIndex(phi,z);
505 //cerr<<"AliITStrackerV2::FollowProlongation: "
506 //"failed to find a detector !\n";
510 //propagate to the intersection
511 const AliITSdetector &det=layer.GetDetector(idet);
513 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
514 //cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
517 fTrackToFollow.SetDetectorIndex(idet);
519 //Select possible prolongations and store the current track estimation
520 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
521 Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
522 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
524 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
528 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) break;
530 Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
531 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
533 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
537 Double_t zmin=track.GetZ() - dz;
538 Double_t zmax=track.GetZ() + dz;
539 Double_t ymin=track.GetY() + r*phi - dy;
540 Double_t ymax=track.GetY() + r*phi + dy;
541 layer.SelectClusters(zmin,zmax,ymin,ymax);
544 //take another prolongation
545 if (!TakeNextProlongation()) if (!tryAgain--) break;
546 tryAgain=kLayersToSkip;
550 //deal with the best track
551 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
552 Int_t nclb=fBestTrack.GetNumberOfClusters();
555 Double_t chi2=fTrackToFollow.GetChi2();
556 if (chi2/ncl < kChi2PerCluster) {
557 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
565 Int_t AliITStrackerV2::TakeNextProlongation() {
566 //--------------------------------------------------------------------
567 // This function takes another track prolongation
569 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
570 //--------------------------------------------------------------------
571 AliITSlayer &layer=fLayers[fI];
572 ResetTrackToFollow(fTracks[fI]);
574 Double_t dz=4*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
575 Double_t dy=4*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
577 const AliITSclusterV2 *c=0; Int_t ci=-1;
578 Double_t chi2=12345.;
579 while ((c=layer.GetNextCluster(ci))!=0) {
580 //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue;
581 Int_t idet=c->GetDetectorIndex();
583 if (fTrackToFollow.GetDetectorIndex()!=idet) {
584 const AliITSdetector &det=layer.GetDetector(idet);
585 ResetTrackToFollow(fTracks[fI]);
586 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
587 //cerr<<"AliITStrackerV2::TakeNextProlongation: "
588 //"propagation failed !\n";
591 fTrackToFollow.SetDetectorIndex(idet);
592 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
595 cout<<fI<<" change detector !\n";
600 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
601 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
603 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
607 cout<<fI<<" chi2="<<chi2<<' '
608 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
611 Double_t phi=TMath::ATan(fTrackToFollow.GetY()/fTrackToFollow.GetX())
612 +fTrackToFollow.GetAlpha(); phi=phi*180./3.141593;
617 if (chi2>=kMaxChi2) return 0;
620 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
621 //cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
625 if (fTrackToFollow.GetNumberOfClusters()>1)
626 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
629 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
633 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
634 fTrackToFollow.CorrectForMaterial(d,x0);
637 if (fConstraint[fPass]) {
638 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
639 fTrackToFollow.Improve(d,GetY(),GetZ());
643 cout<<"accepted lab="<<c->GetLabel(0)<<' '
644 <<fTrackToFollow.GetNumberOfClusters()<<' '
645 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
646 <<fTrackToFollow.Get1Pt()<<endl<<endl;
653 AliITStrackerV2::AliITSlayer::AliITSlayer() {
654 //--------------------------------------------------------------------
655 //default AliITSlayer constructor
656 //--------------------------------------------------------------------
661 AliITStrackerV2::AliITSlayer::
662 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
663 //--------------------------------------------------------------------
664 //main AliITSlayer constructor
665 //--------------------------------------------------------------------
666 fR=r; fPhiOffset=p; fZOffset=z;
667 fNladders=nl; fNdetectors=nd;
668 fDetectors=new AliITSdetector[fNladders*fNdetectors];
674 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
675 //--------------------------------------------------------------------
676 // AliITSlayer destructor
677 //--------------------------------------------------------------------
679 for (Int_t i=0; i<fN; i++) delete fClusters[i];
682 void AliITStrackerV2::AliITSlayer::ResetClusters() {
683 //--------------------------------------------------------------------
684 // This function removes loaded clusters
685 //--------------------------------------------------------------------
686 for (Int_t i=0; i<fN; i++) delete fClusters[i];
691 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
692 //--------------------------------------------------------------------
693 //This function adds a cluster to this layer
694 //--------------------------------------------------------------------
695 if (fN==kMaxClusterPerLayer) {
696 cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): "
697 "Too many clusters !\n";
701 if (fN==0) {fClusters[fN++]=c; return 0;}
702 Int_t i=FindClusterIndex(c->GetZ());
703 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
704 fClusters[i]=c; fN++;
709 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
710 //--------------------------------------------------------------------
711 // This function returns the index of the nearest cluster
712 //--------------------------------------------------------------------
714 if (z <= fClusters[0]->GetZ()) return 0;
715 if (z > fClusters[fN-1]->GetZ()) return fN;
716 Int_t b=0, e=fN-1, m=(b+e)/2;
717 for (; b<e; m=(b+e)/2) {
718 if (z > fClusters[m]->GetZ()) b=m+1;
724 void AliITStrackerV2::AliITSlayer::
725 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
726 //--------------------------------------------------------------------
727 // This function sets the "window"
728 //--------------------------------------------------------------------
729 fI=FindClusterIndex(zmin); fZmax=zmax;
730 Double_t circle=2*TMath::Pi()*fR;
731 if (ymax>circle) { ymax-=circle; ymin-=circle; }
732 fYmin=ymin; fYmax=ymax;
735 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
736 //--------------------------------------------------------------------
737 // This function returns clusters within the "window"
738 //--------------------------------------------------------------------
739 const AliITSclusterV2 *cluster=0;
740 for (Int_t i=fI; i<fN; i++) {
741 const AliITSclusterV2 *c=fClusters[i];
742 if (c->GetZ() > fZmax) break;
743 if (c->IsUsed()) continue;
744 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
745 Double_t y=fR*det.GetPhi() + c->GetY();
747 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
748 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
750 if (y<fYmin) continue;
751 if (y>fYmax) continue;
760 Int_t AliITStrackerV2::AliITSlayer::
761 FindDetectorIndex(Double_t phi, Double_t z) const {
762 //--------------------------------------------------------------------
763 //This function finds the detector crossed by the track
764 //--------------------------------------------------------------------
765 Double_t dphi=phi-fPhiOffset;
766 if (dphi < 0) dphi += 2*TMath::Pi();
767 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
768 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
769 if (np>=fNladders) np-=fNladders;
770 if (np<0) np+=fNladders;
772 Double_t dz=fZOffset-z;
773 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
774 if (nz>=fNdetectors) return -1;
778 cout<<np<<' '<<nz<<endl;
781 return np*fNdetectors + nz;
785 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
787 //--------------------------------------------------------------------
788 //This function returns the layer thickness at this point (units X0)
789 //--------------------------------------------------------------------
793 if (43<fR&&fR<45) { //SSD2
796 if (TMath::Abs(y-0.00)>3.40) d+=dd;
797 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
798 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
799 for (Int_t i=0; i<12; i++) {
800 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
801 if (TMath::Abs(y-0.00)>3.40) d+=dd;
805 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
806 if (TMath::Abs(y-0.00)>3.40) d+=dd;
810 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
811 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
814 if (37<fR&&fR<41) { //SSD1
817 if (TMath::Abs(y-0.00)>3.40) d+=dd;
818 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
819 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
820 for (Int_t i=0; i<11; i++) {
821 if (TMath::Abs(z-3.9*i)<0.15) {
822 if (TMath::Abs(y-0.00)>3.40) d+=dd;
826 if (TMath::Abs(z+3.9*i)<0.15) {
827 if (TMath::Abs(y-0.00)>3.40) d+=dd;
831 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
832 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
835 if (13<fR&&fR<26) { //SDD
838 if (TMath::Abs(y-0.00)>3.30) d+=dd;
840 if (TMath::Abs(y-1.80)<0.55) {
842 for (Int_t j=0; j<20; j++) {
843 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
844 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
847 if (TMath::Abs(y+1.80)<0.55) {
849 for (Int_t j=0; j<20; j++) {
850 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
851 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
855 for (Int_t i=0; i<4; i++) {
856 if (TMath::Abs(z-7.3*i)<0.60) {
858 if (TMath::Abs(y-0.00)>3.30) d+=dd;
861 if (TMath::Abs(z+7.3*i)<0.60) {
863 if (TMath::Abs(y-0.00)>3.30) d+=dd;
868 if (6<fR&&fR<8) { //SPD2
869 Double_t dd=0.0063; x0=21.5;
871 if (TMath::Abs(y-3.08)>0.5) d+=dd;
872 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
873 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
875 if (3<fR&&fR<5) { //SPD1
876 Double_t dd=0.0063; x0=21.5;
878 if (TMath::Abs(y+0.21)>0.6) d+=dd;
879 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
880 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
890 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
892 //--------------------------------------------------------------------
893 //Returns the thickness between the current layer and the vertex (units X0)
894 //--------------------------------------------------------------------
895 Double_t d=0.0028*3*3; //beam pipe
898 Double_t xn=fLayers[fI].GetR();
899 for (Int_t i=0; i<fI; i++) {
900 Double_t xi=fLayers[i].GetR();
901 d+=fLayers[i].GetThickness(y,z,x0)*xi*xi;
910 Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
917 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
918 //--------------------------------------------------------------------
919 // This function returns number of clusters within the "window"
920 //--------------------------------------------------------------------
922 for (Int_t i=fI; i<fN; i++) {
923 const AliITSclusterV2 *c=fClusters[i];
924 if (c->GetZ() > fZmax) break;
925 if (c->IsUsed()) continue;
926 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
927 Double_t y=fR*det.GetPhi() + c->GetY();
929 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
930 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
932 if (y<fYmin) continue;
933 if (y>fYmax) continue;
939 Bool_t AliITStrackerV2::RefitAt(Double_t x, AliITStrackV2 *t, Int_t *index) {
940 //--------------------------------------------------------------------
941 // This function refits a track at a given position
942 //--------------------------------------------------------------------
943 Int_t from, to, step;
945 from=0; to=kMaxLayer;
948 from=kMaxLayer-1; to=-1;
952 for (Int_t i=from; i != to; i += step) {
953 AliITSlayer &layer=fLayers[i];
954 Double_t r=layer.GetR();
957 Double_t hI=i-0.5*step;
958 if (hI==1.5 || hI==3.5) {
959 Double_t rs=0.5*(fLayers[i-step].GetR() + r);
960 Double_t d=0.0034, x0=38.6;
961 if (hI==1.5) {rs=9.; d=0.0097; x0=42;}
962 if (!t->PropagateTo(rs,d,x0)) {
969 if (!t->GetGlobalXYZat(r,x,y,z)) {
972 Double_t phi=TMath::ATan2(y,x);
973 Int_t idet=layer.FindDetectorIndex(phi,z);
977 const AliITSdetector &det=layer.GetDetector(idet);
979 if (!t->Propagate(phi,det.GetR())) {
982 t->SetDetectorIndex(idet);
984 const AliITSclusterV2 *cl=0;
985 Double_t maxchi2=kMaxChi2;
989 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
990 if (idet != c->GetDetectorIndex()) {
991 idet=c->GetDetectorIndex();
992 const AliITSdetector &det=layer.GetDetector(idet);
993 if (!t->Propagate(det.GetPhi(),det.GetR())) {
996 t->SetDetectorIndex(idet);
998 Double_t chi2=t->GetPredictedChi2(c);
999 if (chi2<maxchi2) { cl=c; maxchi2=chi2; }
1005 if (t->GetNumberOfClusters()>2) {
1006 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1007 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1008 Double_t zmin=t->GetZ() - dz;
1009 Double_t zmax=t->GetZ() + dz;
1010 Double_t ymin=t->GetY() + phi*r - dy;
1011 Double_t ymax=t->GetY() + phi*r + dy;
1012 layer.SelectClusters(zmin,zmax,ymin,ymax);
1014 const AliITSclusterV2 *c=0; Int_t ci=-1;
1015 while ((c=layer.GetNextCluster(ci))!=0) {
1016 if (idet != c->GetDetectorIndex()) continue;
1017 Double_t chi2=t->GetPredictedChi2(c);
1018 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1024 if (!t->Update(cl,maxchi2,idx)) {
1031 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1032 t->CorrectForMaterial(-step*d,x0);
1037 if (!t->PropagateTo(x,0.,0.)) return kFALSE;