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
42 AliITStrackerV2(const AliITSgeom *geom, Int_t eventn) throw (const Char_t *) :
44 //--------------------------------------------------------------------
45 //This is the AliITStracker constructor
46 //It also reads clusters from a root file
47 //--------------------------------------------------------------------
48 fEventN=eventn; //MI change add event number - used to generate identifier
50 AliITSgeom *g=(AliITSgeom*)geom;
54 for (i=1; i<kMaxLayer+1; i++) {
55 Int_t nlad=g->GetNladders(i);
56 Int_t ndet=g->GetNdetectors(i);
58 g->GetTrans(i,1,1,x,y,z);
59 Double_t r=TMath::Sqrt(x*x + y*y);
60 Double_t poff=TMath::ATan2(y,x);
63 g->GetTrans(i,1,2,x,y,z);
64 r += TMath::Sqrt(x*x + y*y);
65 g->GetTrans(i,2,1,x,y,z);
66 r += TMath::Sqrt(x*x + y*y);
67 g->GetTrans(i,2,2,x,y,z);
68 r += TMath::Sqrt(x*x + y*y);
71 new (fLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
73 for (Int_t j=1; j<nlad+1; j++) {
74 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
75 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
76 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
78 Double_t r =-x*rot[1] + y*rot[0]; if (i==1) r=-r;
79 Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927;
81 AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1);
83 if (phi<0) phi += 2*TMath::Pi();
85 new(&det) AliITSdetector(r,phi);
95 sprintf(cname,"TreeC_ITS_%d",eventn);
96 TTree *cTree=(TTree*)gDirectory->Get(cname);
99 ("AliITStrackerV2::AliITStrackerV2 can't get cTree !\n");
101 TBranch *branch=cTree->GetBranch("Clusters");
103 ("AliITStrackerV2::AliITStrackerV2 can't get Clusters branch !\n");
105 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
106 branch->SetAddress(&clusters);
108 Int_t nentr=(Int_t)cTree->GetEntries();
109 for (i=0; i<nentr; i++) {
110 if (!cTree->GetEvent(i)) continue;
111 Int_t lay,lad,det; g->GetModuleId(i,lay,lad,det);
112 Int_t ncl=clusters->GetEntriesFast();
114 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
117 if (c->GetLabel(2)!=LAB)
118 if (c->GetLabel(1)!=LAB)
119 if (c->GetLabel(0)!=LAB) continue;
120 cout<<lay-1<<' '<<lad-1<<' '<<det-1<<' '<<c->GetY()<<' '<<c->GetZ()<<endl;
123 fLayers[lay-1].InsertCluster(new AliITSclusterV2(*c));
127 delete cTree; //Thanks to Mariana Bondila
129 catch (const Char_t *msg) {
137 fConstraint[0]=1; fConstraint[1]=0;
144 Int_t AliITStrackerV2::Clusters2Tracks(const TFile *inp, TFile *out) {
145 //--------------------------------------------------------------------
146 //This functions reconstructs ITS tracks
147 //--------------------------------------------------------------------
148 TFile *in=(TFile*)inp;
149 TDirectory *savedir=gDirectory;
152 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
153 cerr<<"file with TPC tracks is not open !\n";
157 if (!out->IsOpen()) {
158 cerr<<"AliITStrackerV2::Clusters2Tracks(): ";
159 cerr<<"file for ITS tracks is not open !\n";
166 Int_t nentr=0; TObjArray itsTracks(15000);
168 {/* Read TPC tracks */
169 sprintf(tname,"TreeT_TPC_%d",fEventN);
170 TTree *tpcTree=(TTree*)in->Get(tname);
172 cerr<<"AliITStrackerV2::Clusters2Tracks(): "
173 "can't get a tree with TPC tracks !\n";
176 AliTPCtrack *itrack=new AliTPCtrack;
177 tpcTree->SetBranchAddress("tracks",&itrack);
178 nentr=(Int_t)tpcTree->GetEntries();
179 for (Int_t i=0; i<nentr; i++) {
180 tpcTree->GetEvent(i);
183 t=new AliITStrackV2(*itrack);
184 } catch (const Char_t *msg) {
189 if (TMath::Abs(t->GetD())>4) continue;
191 t->PropagateTo(80.,0.0053);
192 if (TMath::Abs(t->GetY())>13.) t->CorrectForMaterial(0.03);
193 if (TMath::Abs(t->GetZ())<0.2) t->CorrectForMaterial(0.40);
194 t->PropagateTo(61.,0.0052);
195 Double_t xk=52.,x,y,z; t->GetGlobalXYZat(xk,x,y,z);
196 if (TMath::Abs(y)<7.77) t->PropagateTo(xk,0.19,24.);
197 t->PropagateTo(50.,0.001);
199 itsTracks.AddLast(t);
201 delete tpcTree; //Thanks to Mariana Bondila
208 sprintf(tname,"TreeT_ITS_%d",fEventN);
209 TTree itsTree(tname,"Tree with ITS tracks");
210 AliITStrackV2 *otrack=&fBestTrack;
211 itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
213 for (fPass=0; fPass<2; fPass++) {
214 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
215 for (Int_t i=0; i<nentr; i++) {
216 if (i%10==0) cerr<<nentr-i<<" \r";
217 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
218 if (t==0) continue; //this track has been already tracked
219 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
224 if (TMath::Abs(tpcLabel)!=LAB) continue;
225 cout<<tpcLabel<<" *****************\n";
228 ResetTrackToFollow(*t);
231 for (FollowProlongation(); fI<kMaxLayer; fI++) {
232 while (TakeNextProlongation()) FollowProlongation();
236 cout<<fBestTrack.GetNumberOfClusters()<<" number of clusters\n\n";
239 if (fBestTrack.GetNumberOfClusters() < kMaxLayer-kLayersToSkip) continue;
241 if (fConstraint[fPass]) {
242 Int_t index[kMaxLayer];
244 for (k=0; k<kMaxLayer; k++) index[k]=-1;
245 Int_t nc=fBestTrack.GetNumberOfClusters();
246 for (k=0; k<nc; k++) {
247 Int_t idx=fBestTrack.GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
250 fBestTrack.~AliITStrackV2(); new(&fBestTrack) AliITStrackV2(*t);
251 if (!RefitAt(3.7, &fBestTrack, index)) continue;
254 fBestTrack.SetLabel(tpcLabel);
255 fBestTrack.CookdEdx();
256 CookLabel(&fBestTrack,0.); //For comparison only
258 UseClusters(&fBestTrack);
259 delete itsTracks.RemoveAt(i);
268 cerr<<"Number of TPC tracks: "<<nentr<<endl;
269 cerr<<"Number of prolonged tracks: "<<itsTree.GetEntries()<<endl;
274 Int_t AliITStrackerV2::PropagateBack(const TFile *inp, TFile *out) {
275 //--------------------------------------------------------------------
276 //This functions propagates reconstructed ITS tracks back
277 //--------------------------------------------------------------------
278 TFile *in=(TFile*)inp;
279 TDirectory *savedir=gDirectory;
282 cerr<<"AliITStrackerV2::PropagateBack(): ";
283 cerr<<"file with ITS tracks is not open !\n";
287 if (!out->IsOpen()) {
288 cerr<<"AliITStrackerV2::PropagateBack(): ";
289 cerr<<"file for back propagated ITS tracks is not open !\n";
294 TTree *itsTree=(TTree*)in->Get("TreeT_ITS_0");
296 cerr<<"AliITStrackerV2::PropagateBack() ";
297 cerr<<"can't get a tree with ITS tracks !\n";
300 AliITStrackV2 *itrack=new AliITStrackV2;
301 itsTree->SetBranchAddress("tracks",&itrack);
304 TTree backTree("TreeT_ITSb_0","Tree with back propagated ITS tracks");
305 AliTPCtrack *otrack=0;
306 backTree.Branch("tracks","AliTPCtrack",&otrack,32000,0);
310 Int_t nentr=(Int_t)itsTree->GetEntries();
311 for (Int_t i=0; i<nentr; i++) {
312 itsTree->GetEvent(i);
313 ResetTrackToFollow(*itrack);
314 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
315 Int_t itsLabel=fTrackToFollow.GetLabel(); //save the ITS track label
318 if (TMath::Abs(itsLabel)!=LAB) continue;
319 cout<<itsLabel<<" *****************\n";
323 Int_t nc=itrack->GetNumberOfClusters();
325 for (Int_t k=0; k<nc; k++) {
326 Int_t index=itrack->GetClusterIndex(k);
327 AliITSclusterV2 *c=(AliITSclusterV2*)GetCluster(index);
328 cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
332 const AliITSclusterV2 *c=0;
334 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
335 c=(AliITSclusterV2*)GetCluster(idx);
337 for (fI=0; fI<kMaxLayer; fI++) {
338 AliITSlayer &layer=fLayers[fI];
339 Double_t r=layer.GetR();
340 if (fI==2 || fI==4) {
341 Double_t rs=0.5*(fLayers[fI-1].GetR() + r);
342 Double_t d=0.011; if (fI==4) d=0.0053;
343 if (!fTrackToFollow.PropagateTo(rs,-d)) throw "";
347 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z))
348 throw "AliITStrackerV2::PropagateBack: failed to estimate track !";
349 Double_t phi=TMath::ATan2(y,x);
350 Int_t idet=layer.FindDetectorIndex(phi,z);
352 throw "AliITStrackerV2::PropagateBack: failed to find a detector !\n";
353 const AliITSdetector &det=layer.GetDetector(idet);
354 r=det.GetR(); phi=det.GetPhi();
355 if (!fTrackToFollow.Propagate(phi,r)) throw "";
356 fTrackToFollow.SetDetectorIndex(idet);
358 const AliITSclusterV2 *cl=0;
360 Double_t maxchi2=kMaxChi2;
363 idet=c->GetDetectorIndex();
364 if (idet != fTrackToFollow.GetDetectorIndex()) {
365 const AliITSdetector &det=layer.GetDetector(idet);
366 r=det.GetR(); phi=det.GetPhi();
367 if (!fTrackToFollow.Propagate(phi,r)) throw "";
368 fTrackToFollow.SetDetectorIndex(idet);
370 Double_t chi2=fTrackToFollow.GetPredictedChi2(c);
372 cl=c; maxchi2=chi2; index=idx;
375 idx=itrack->GetClusterIndex(nc); l=(idx&0xf0000000)>>28;
376 c=(AliITSclusterV2*)GetCluster(idx);
380 if (fTrackToFollow.GetNumberOfClusters()>2) {
381 Double_t dz=3*TMath::Sqrt(fTrackToFollow.GetSigmaZ2()+kSigmaZ2[fI]);
382 Double_t dy=3*TMath::Sqrt(fTrackToFollow.GetSigmaY2()+kSigmaY2[fI]);
383 Double_t zmin=fTrackToFollow.GetZ() - dz;
384 Double_t zmax=fTrackToFollow.GetZ() + dz;
385 Double_t ymin=fTrackToFollow.GetY() + phi*r - dy;
386 Double_t ymax=fTrackToFollow.GetY() + phi*r + dy;
387 layer.SelectClusters(zmin,zmax,ymin,ymax);
389 const AliITSclusterV2 *cc=0; Int_t ci;
390 while ((cc=layer.GetNextCluster(ci))!=0) {
391 idet=cc->GetDetectorIndex();
392 if (idet != fTrackToFollow.GetDetectorIndex()) continue;
393 Double_t chi2=fTrackToFollow.GetPredictedChi2(cc);
395 cl=cc; index=(fI<<28)+ci; maxchi2=chi2;
401 if (!fTrackToFollow.Update(cl,maxchi2,index))
402 cerr<<"AliITStrackerV2::PropagateBack: filtering failed !\n";
405 x=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ());
406 fTrackToFollow.CorrectForMaterial(-x);
410 Double_t xk=52.,x,y,z; fTrackToFollow.GetGlobalXYZat(xk,x,y,z);
411 if (TMath::Abs(y)<7.77)
412 fTrackToFollow.PropagateTo(xk,-0.19,24.);
413 fTrackToFollow.PropagateTo(61,-0.0110);
414 fTrackToFollow.PropagateTo(80.,-0.0053);
416 fTrackToFollow.SetLabel(itsLabel);
417 otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha());
418 backTree.Fill(); delete otrack;
419 UseClusters(&fTrackToFollow);
420 cerr << ++ntrk << " \r";
422 catch (const Char_t *msg) {
429 cerr<<"Number of ITS tracks: "<<nentr<<endl;
430 cerr<<"Number of back propagated ITS tracks: "<<ntrk<<endl;
434 delete itsTree; //Thanks to Mariana Bondila
440 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
441 //--------------------------------------------------------------------
442 // Return pointer to a given cluster
443 //--------------------------------------------------------------------
444 Int_t l=(index & 0xf0000000) >> 28;
445 Int_t c=(index & 0x0fffffff) >> 00;
446 return fLayers[l].GetCluster(c);
450 void AliITStrackerV2::FollowProlongation() {
451 //--------------------------------------------------------------------
452 //This function finds a track prolongation
453 //--------------------------------------------------------------------
454 Int_t tryAgain=kLayersToSkip;
461 AliITSlayer &layer=fLayers[i];
462 AliITStrackV2 &track=fTracks[i];
464 Double_t r=layer.GetR();
466 Double_t rs=0.5*(fLayers[i+1].GetR() + r);
467 Double_t d=0.011; if (i==3) d=0.0053;
468 if (!fTrackToFollow.PropagateTo(rs,d)) {
469 //cerr<<"AliITStrackerV2::FollowProlongation: "
470 //"propagation failed !\n";
477 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
478 //cerr<<"AliITStrackerV2::FollowProlongation: "
479 //"failed to estimate track !\n";
482 Double_t phi=TMath::ATan2(y,x);
483 Int_t idet=layer.FindDetectorIndex(phi,z);
485 //cerr<<"AliITStrackerV2::FollowProlongation: "
486 //"failed to find a detector !\n";
490 //propagate to the intersection
491 const AliITSdetector &det=layer.GetDetector(idet);
493 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
494 //cerr<<"AliITStrackerV2::FollowProlongation: propagation failed !\n";
497 fTrackToFollow.SetDetectorIndex(idet);
499 //Select possible prolongations and store the current track estimation
500 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
501 Double_t dz=3*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
502 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
504 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n";
508 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) break;
510 Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
511 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
513 //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n";
517 Double_t zmin=track.GetZ() - dz;
518 Double_t zmax=track.GetZ() + dz;
519 Double_t ymin=track.GetY() + r*phi - dy;
520 Double_t ymax=track.GetY() + r*phi + dy;
521 layer.SelectClusters(zmin,zmax,ymin,ymax);
524 //take another prolongation
525 if (!TakeNextProlongation()) if (!tryAgain--) break;
526 tryAgain=kLayersToSkip;
530 //deal with the best track
531 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
532 Int_t nclb=fBestTrack.GetNumberOfClusters();
535 Double_t chi2=fTrackToFollow.GetChi2();
536 if (chi2/ncl < kChi2PerCluster) {
537 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
545 Int_t AliITStrackerV2::TakeNextProlongation() {
546 //--------------------------------------------------------------------
547 // This function takes another track prolongation
549 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
550 //--------------------------------------------------------------------
551 AliITSlayer &layer=fLayers[fI];
552 ResetTrackToFollow(fTracks[fI]);
554 Double_t dz=4*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
555 Double_t dy=4*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
557 const AliITSclusterV2 *c=0; Int_t ci=-1;
558 Double_t chi2=12345.;
559 while ((c=layer.GetNextCluster(ci))!=0) {
560 //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue;
561 Int_t idet=c->GetDetectorIndex();
563 if (fTrackToFollow.GetDetectorIndex()!=idet) {
564 const AliITSdetector &det=layer.GetDetector(idet);
565 ResetTrackToFollow(fTracks[fI]);
566 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
567 //cerr<<"AliITStrackerV2::TakeNextProlongation: "
568 //"propagation failed !\n";
571 fTrackToFollow.SetDetectorIndex(idet);
572 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
575 cout<<fI<<" change detector !\n";
580 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
581 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
583 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
587 cout<<fI<<" chi2="<<chi2<<' '
588 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
592 if (chi2>=kMaxChi2) return 0;
595 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
596 //cerr<<"AliITStrackerV2::TakeNextProlongation: filtering failed !\n";
600 if (fTrackToFollow.GetNumberOfClusters()>1)
601 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
604 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
607 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ());
608 fTrackToFollow.CorrectForMaterial(d);
611 if (fConstraint[fPass]) {
612 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
613 fTrackToFollow.Improve(d,GetY(),GetZ());
617 cout<<"accepted lab="<<c->GetLabel(0)<<' '
618 <<fTrackToFollow.GetNumberOfClusters()<<' '
619 <<fTrackToFollow.GetY()<<' '<<fTrackToFollow.GetZ()<<' '
620 <<fTrackToFollow.Get1Pt()<<endl<<endl;
627 AliITStrackerV2::AliITSlayer::AliITSlayer() {
628 //--------------------------------------------------------------------
629 //default AliITSlayer constructor
630 //--------------------------------------------------------------------
635 AliITStrackerV2::AliITSlayer::
636 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
637 //--------------------------------------------------------------------
638 //main AliITSlayer constructor
639 //--------------------------------------------------------------------
640 fR=r; fPhiOffset=p; fZOffset=z;
641 fNladders=nl; fNdetectors=nd;
642 fDetectors=new AliITSdetector[fNladders*fNdetectors];
648 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
649 //--------------------------------------------------------------------
650 // AliITSlayer destructor
651 //--------------------------------------------------------------------
653 for (Int_t i=0; i<fN; i++) delete fClusters[i];
656 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
657 //--------------------------------------------------------------------
658 //This function adds a cluster to this layer
659 //--------------------------------------------------------------------
660 if (fN==kMaxClusterPerLayer) {
661 cerr<<"AliITStrackerV2::AliITSlayer::InsertCluster(): "
662 "Too many clusters !\n";
666 if (fN==0) {fClusters[fN++]=c; return 0;}
667 Int_t i=FindClusterIndex(c->GetZ());
668 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
669 fClusters[i]=c; fN++;
674 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
675 //--------------------------------------------------------------------
676 // This function returns the index of the nearest cluster
677 //--------------------------------------------------------------------
679 if (z <= fClusters[0]->GetZ()) return 0;
680 if (z > fClusters[fN-1]->GetZ()) return fN;
681 Int_t b=0, e=fN-1, m=(b+e)/2;
682 for (; b<e; m=(b+e)/2) {
683 if (z > fClusters[m]->GetZ()) b=m+1;
689 void AliITStrackerV2::AliITSlayer::
690 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
691 //--------------------------------------------------------------------
692 // This function sets the "window"
693 //--------------------------------------------------------------------
694 fI=FindClusterIndex(zmin); fZmax=zmax;
695 Double_t circle=2*TMath::Pi()*fR;
696 if (ymax>circle) { ymax-=circle; ymin-=circle; }
697 fYmin=ymin; fYmax=ymax;
700 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
701 //--------------------------------------------------------------------
702 // This function returns clusters within the "window"
703 //--------------------------------------------------------------------
704 const AliITSclusterV2 *cluster=0;
705 for (Int_t i=fI; i<fN; i++) {
706 const AliITSclusterV2 *c=fClusters[i];
707 if (c->GetZ() > fZmax) break;
708 if (c->IsUsed()) continue;
709 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
710 Double_t y=fR*det.GetPhi() + c->GetY();
712 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
713 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
715 if (y<fYmin) continue;
716 if (y>fYmax) continue;
725 Int_t AliITStrackerV2::AliITSlayer::
726 FindDetectorIndex(Double_t phi, Double_t z) const {
727 //--------------------------------------------------------------------
728 //This function finds the detector crossed by the track
729 //--------------------------------------------------------------------
730 Double_t dphi=phi-fPhiOffset;
731 if (dphi < 0) dphi += 2*TMath::Pi();
732 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
733 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
734 if (np>=fNladders) np-=fNladders;
735 if (np<0) np+=fNladders;
737 Double_t dz=fZOffset-z;
738 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
739 if (nz>=fNdetectors) return -1;
743 cout<<np<<' '<<nz<<endl;
746 return np*fNdetectors + nz;
750 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y, Double_t z) const {
751 //--------------------------------------------------------------------
752 //This function returns the layer thickness at this point (units X0)
753 //--------------------------------------------------------------------
756 if (43<fR&&fR<45) { //SSD2
758 if (TMath::Abs(y-0.00)>3.40) d+=0.0036;
759 if (TMath::Abs(y-2.50)<0.10) d+=(0.02-0.0036);
760 if (TMath::Abs(y+1.90)<0.10) d+=(0.02-0.0036);
761 for (Int_t i=0; i<12; i++) {
762 if (TMath::Abs(z-3.6*(i+0.5))<0.20) {d+=0.0036; break;}
763 if (TMath::Abs(z+3.6*(i+0.5))<0.20) {d+=0.0036; break;}
764 if (TMath::Abs(z-3.6*(i+0.929))<0.50) {d+=(0.02-0.0036); break;}
765 if (TMath::Abs(z+3.6*(i+0.104))<0.50) {d+=(0.02-0.0036); break;}
768 if (37<fR&&fR<41) { //SSD1
770 if (TMath::Abs(y-0.00)>3.40) d+=0.0036;
771 if (TMath::Abs(y-2.20)<0.10) d+=(0.02-0.0036);
772 if (TMath::Abs(y+2.20)<0.10) d+=(0.02-0.0036);
773 for (Int_t i=0; i<11; i++) {
774 if (TMath::Abs(z-3.6*i)<0.20) {d+=0.0036; break;}
775 if (TMath::Abs(z+3.6*i)<0.20) {d+=0.0036; break;}
776 if (TMath::Abs(z-3.6*(i+0.54))<0.50) {d+=(0.02-0.0036); break;}
777 if (TMath::Abs(z+3.6*(i+0.58))<0.50) {d+=(0.02-0.0036); break;}
780 if (13<fR&&fR<26) { //SDD
782 if (TMath::Abs(y-0.00)>3.30) d+=0.0034;
783 if (TMath::Abs(y-2.10)<0.20) d+=0.0034*3;
784 if (TMath::Abs(y+2.10)<0.20) d+=0.0034*3;
785 for (Int_t i=0; i<4; i++) {
786 if (TMath::Abs(z-7.3*i)<0.60) {d+=0.0034; break;}
787 if (TMath::Abs(z+7.3*i)<0.60) {d+=0.0034; break;}
790 if (6<fR&&fR<8) { //SPD2
792 if (TMath::Abs(y-3.08)>0.45) d+=0.0064;
793 if (TMath::Abs(y-3.03)<0.10) d+=0.0192;
795 if (3<fR&&fR<5) { //SPD1
797 if (TMath::Abs(y+0.21)>0.55) d+=0.0064;
798 if (TMath::Abs(y+0.10)<0.10) d+=0.0192;
806 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
808 //--------------------------------------------------------------------
809 //Returns the thickness between the current layer and the vertex (units X0)
810 //--------------------------------------------------------------------
811 Double_t d=0.1/65.19*1.848;
813 Double_t xn=fLayers[fI].GetR();
814 for (Int_t i=0; i<fI; i++) {
815 Double_t xi=fLayers[i].GetR();
816 d+=fLayers[i].GetThickness(y,z)*xi*xi;
820 Double_t xi=0.5*(fLayers[1].GetR()+fLayers[2].GetR());
825 Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR());
833 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
834 //--------------------------------------------------------------------
835 // This function returns number of clusters within the "window"
836 //--------------------------------------------------------------------
838 for (Int_t i=fI; i<fN; i++) {
839 const AliITSclusterV2 *c=fClusters[i];
840 if (c->GetZ() > fZmax) break;
841 if (c->IsUsed()) continue;
842 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
843 Double_t y=fR*det.GetPhi() + c->GetY();
845 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
846 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
848 if (y<fYmin) continue;
849 if (y>fYmax) continue;
855 Bool_t AliITStrackerV2::RefitAt(Double_t x, AliITStrackV2 *t, Int_t *index) {
856 //--------------------------------------------------------------------
857 // This function refits a track at a given position
858 //--------------------------------------------------------------------
859 Int_t from, to, step;
861 from=0; to=kMaxLayer;
864 from=kMaxLayer-1; to=-1;
868 for (Int_t i=from; i != to; i += step) {
869 AliITSlayer &layer=fLayers[i];
870 Double_t r=layer.GetR();
873 Double_t hI=i-0.5*step;
874 if (hI==1.5 || hI==3.5) {
875 Double_t rs=0.5*(fLayers[i-step].GetR() + r);
876 Double_t ds=0.011; if (hI==3.5) ds=0.0053;
877 if (!t->PropagateTo(rs,ds)) {
884 if (!t->GetGlobalXYZat(r,x,y,z)) {
887 Double_t phi=TMath::ATan2(y,x);
888 Int_t idet=layer.FindDetectorIndex(phi,z);
892 const AliITSdetector &det=layer.GetDetector(idet);
894 if (!t->Propagate(phi,det.GetR())) {
897 t->SetDetectorIndex(idet);
899 const AliITSclusterV2 *cl=0;
900 Double_t maxchi2=kMaxChi2;
904 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
905 if (idet != c->GetDetectorIndex()) {
906 idet=c->GetDetectorIndex();
907 const AliITSdetector &det=layer.GetDetector(idet);
908 if (!t->Propagate(det.GetPhi(),det.GetR())) {
911 t->SetDetectorIndex(idet);
913 Double_t chi2=t->GetPredictedChi2(c);
914 if (chi2<maxchi2) { cl=c; maxchi2=chi2; }
920 if (t->GetNumberOfClusters()>2) {
921 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
922 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
923 Double_t zmin=t->GetZ() - dz;
924 Double_t zmax=t->GetZ() + dz;
925 Double_t ymin=t->GetY() + phi*r - dy;
926 Double_t ymax=t->GetY() + phi*r + dy;
927 layer.SelectClusters(zmin,zmax,ymin,ymax);
929 const AliITSclusterV2 *c=0; Int_t ci=-1;
930 while ((c=layer.GetNextCluster(ci))!=0) {
931 if (idet != c->GetDetectorIndex()) continue;
932 Double_t chi2=t->GetPredictedChi2(c);
933 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
939 if (!t->Update(cl,maxchi2,idx)) {
945 Double_t d=layer.GetThickness(t->GetY(),t->GetZ());
946 t->CorrectForMaterial(-step*d);
951 if (!t->PropagateTo(x,0.,0.)) return kFALSE;