/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ //------------------------------------------------------------------------- // Implementation of the ITS tracker class // // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch //------------------------------------------------------------------------- #include #include #include #include #include "AliITSgeom.h" #include "AliITSRecPoint.h" #include "../TPC/AliTPCtrack.h" #include "AliITSclusterV2.h" #include "AliITStrackerV2.h" //#define DEBUG #ifdef DEBUG Int_t LAB=236; #endif extern TRandom *gRandom; AliITStrackerV2::AliITSlayer AliITStrackerV2::fLayers[kMaxLayer]; // ITS layers AliITStrackerV2:: AliITStrackerV2(const AliITSgeom *geom) throw (const Char_t *) { //-------------------------------------------------------------------- //This is the AliITStracker constructor //It also reads clusters from a root file //-------------------------------------------------------------------- fYV=fZV=0.; AliITSgeom *g=(AliITSgeom*)geom; Float_t x,y,z; Int_t i; for (i=1; iGetNladders(i); Int_t ndet=g->GetNdetectors(i); g->GetTrans(i,1,1,x,y,z); Double_t r=TMath::Sqrt(x*x + y*y); Double_t poff=TMath::ATan2(y,x); Double_t zoff=z; g->GetTrans(i,1,2,x,y,z); r += TMath::Sqrt(x*x + y*y); g->GetTrans(i,2,1,x,y,z); r += TMath::Sqrt(x*x + y*y); g->GetTrans(i,2,2,x,y,z); r += TMath::Sqrt(x*x + y*y); r*=0.25; new (fLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet); for (Int_t j=1; jGetTrans(i,j,k,x,y,zshift); Double_t rot[9]; g->GetRotMatrix(i,j,k,rot); Double_t r =-x*rot[1] + y*rot[0]; if (i==1) r=-r; Double_t phi=TMath::ATan2(rot[1],rot[0]); if (i==1) phi-=3.1415927; phi+=0.5*TMath::Pi(); AliITSdetector &det=fLayers[i-1].GetDetector((j-1)*ndet + k-1); if (phi<0) phi += 2*TMath::Pi(); new(&det) AliITSdetector(r,phi); } } } try { //Read clusters TTree *cTree=(TTree*)gDirectory->Get("cTree"); if (!cTree) throw ("AliITStrackerV2::AliITStrackerV2 can't get cTree !\n"); TBranch *branch=cTree->GetBranch("Clusters"); if (!branch) throw ("AliITStrackerV2::AliITStrackerV2 can't get Clusters branch !\n"); TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy; branch->SetAddress(&clusters); Int_t nentr=(Int_t)cTree->GetEntries(); for (i=0; iGetEvent(i)) continue; Int_t lay,lad,det; g->GetModuleId(i-1,lay,lad,det); Int_t ncl=clusters->GetEntriesFast(); while (ncl--) { AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl); #ifdef DEBUG if (c->GetLabel(2)!=LAB) if (c->GetLabel(1)!=LAB) if (c->GetLabel(0)!=LAB) continue; cout<GetY()<<' '<GetZ()<Delete(); } } catch (const Char_t *msg) { cerr<IsOpen()) { cerr<<"AliITStrackerV2::Clusters2Tracks(): "; cerr<<"file with TPC tracks is not open !\n"; return 1; } if (!out->IsOpen()) { cerr<<"AliITStrackerV2::Clusters2Tracks(): "; cerr<<"file for ITS tracks is not open !\n"; return 2; } in->cd(); TTree *tpcTree=(TTree*)in->Get("TreeT"); if (!tpcTree) { cerr<<"AliITStrackerV2::Clusters2Tracks() "; cerr<<"can't get a tree with TPC tracks !\n"; return 3; } AliTPCtrack *itrack=new AliTPCtrack; tpcTree->SetBranchAddress("tracks",&itrack); out->cd(); TTree itsTree("TreeT","Tree with ITS tracks"); AliITStrackV2 *otrack=&fBestTrack; itsTree.Branch("tracks","AliITStrackV2",&otrack,32000,0); Int_t ntrk=0; Int_t nentr=(Int_t)tpcTree->GetEntries(); for (Int_t i=0; iGetEvent(i)) continue; Int_t tpcLabel=itrack->GetLabel(); //save the TPC track label lbl=tpcLabel; #ifdef DEBUG if (TMath::Abs(tpcLabel)!=LAB) continue; cout<Gaus(0.,kSigmaYV); fZV=gRandom->Gaus(0.,kSigmaZV); Double_t r2=fTrackToFollow.GetX()*fTrackToFollow.GetX(); Double_t x0=0.082/21.82*2.33*(45*45+40*40)/r2+2.000/41.28*0.03*75*75/r2; fTrackToFollow.Improve(x0,fYV,fZV); //Double_t xk=77.2; Double_t xk=88.; fTrackToFollow.PropagateTo(xk,0.,0.); //Ne if it's still there xk-=0.005; fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar xk-=0.02; fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar xk-=2.0; fTrackToFollow.PropagateTo(xk, 2.0/41.28*0.029, 2.0*0.029);//Nomex xk-=0.02; fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar xk-=0.005; fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar xk-=14.5; fTrackToFollow.PropagateTo(xk,0.,0.); //C02 xk -=0.005; fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7); //Al xk -=0.005; fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar xk -=0.02; fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar xk -=0.5; fTrackToFollow.PropagateTo(xk, 0.5/41.28*.029, 0.5*0.029); //Nomex xk -=0.02; fTrackToFollow.PropagateTo(xk, 0.02/44.86*1.45, 0.02*1.45); //Kevlar xk -=0.005; fTrackToFollow.PropagateTo(xk, 0.005/44.77*1.71, 0.005*1.71); //Tedlar xk -=0.005; fTrackToFollow.PropagateTo(xk, 0.005/24.01*2.7, 0.005*2.7); //Al for (FollowProlongation(); fI= kMaxLayer-kLayersToSkip) { cerr << ++ntrk << " \r"; fBestTrack.SetLabel(tpcLabel); UseClusters(&fBestTrack); itsTree.Fill(); } } itsTree.Write(); savedir->cd(); cerr<<"Number of TPC tracks: "<> 28; Int_t c=(index & 0x0fffffff) >> 00; return fLayers[l].GetCluster(c); } void AliITStrackerV2::FollowProlongation() { //-------------------------------------------------------------------- //This function finds a track prolongation //-------------------------------------------------------------------- Int_t tryAgain=kLayersToSkip; while (fI) { fI--; #ifdef DEBUG cout< kMaxRoad/4) { //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Z !\n"; break; } Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[fI]); if (dy<0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp()); if (dy > kMaxRoad/4) { //cerr<<"AliITStrackerV2::FollowProlongation: too broad road in Y !\n"; break; } Double_t zmin=track.GetZ() - dz; Double_t zmax=track.GetZ() + dz; Double_t ymin=track.GetY() + layer.GetR()*det.GetPhi() - dy; Double_t ymax=track.GetY() + layer.GetR()*det.GetPhi() + dy; if (ymax>layer.GetR()*2*TMath::Pi()) { ymax-=layer.GetR()*2*TMath::Pi(); ymin-=layer.GetR()*2*TMath::Pi(); } layer.SelectClusters(zmin,zmax,ymin,ymax); //if (1/TMath::Abs(track.Get1Pt())<0.2) //cout<= nclb) { Double_t chi2=fTrackToFollow.GetChi2(); if (chi2/ncl < kChi2PerCluster) { if (ncl > nclb || chi2 < fBestTrack.GetChi2()) { ResetBestTrack(); } } } if (fI) fI++; } Int_t AliITStrackerV2::TakeNextProlongation() { //-------------------------------------------------------------------- //This function takes another track prolongation //-------------------------------------------------------------------- //Double_t m[20]; Double_t d=GetEffectiveThickness(0,0); //Think of this !!!! AliITSlayer &layer=fLayers[fI]; AliITStrackV2 &t=fTracks[fI]; Double_t dz=4*TMath::Sqrt(t.GetSigmaZ2() + kSigmaZ2[fI]); Double_t dy=4*TMath::Sqrt(t.GetSigmaY2() + kSigmaY2[fI]); const AliITSclusterV2 *c=0; Int_t ci=-1; Double_t chi2=12345.; while ((c=layer.GetNextCluster(ci))!=0) { //if (fI==0) //if (c->GetLabel(0)!=TMath::Abs(lbl)) continue; //88888888888888888888 Int_t idet=c->GetDetectorIndex(); if (t.GetDetectorIndex()!=idet) { const AliITSdetector &det=layer.GetDetector(idet); if (!t.Propagate(det.GetPhi(),det.GetR(),0.,0.)) { cerr<<"AliITStrackerV2::TakeNextProlongation: propagation failed !\n"; continue; } t.SetDetectorIndex(idet); #ifdef DEBUG cout<GetZ()) > dz) continue; if (TMath::Abs(t.GetY() - c->GetY()) > dy) continue; //m[0]=fYV; m[1]=fZV; //chi2=t.GetPredictedChi2(c,m,d/21.82*2.33); chi2=t.GetPredictedChi2(c); if (chi2GetZ()); memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*)); fClusters[i]=c; fN++; return 0; } Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const { //-------------------------------------------------------------------- // This function returns the index of the nearest cluster //-------------------------------------------------------------------- if (fN==0) return 0; if (z <= fClusters[0]->GetZ()) return 0; if (z > fClusters[fN-1]->GetZ()) return fN; Int_t b=0, e=fN-1, m=(b+e)/2; for (; b fClusters[m]->GetZ()) b=m+1; else e=m; } return m; } void AliITStrackerV2::AliITSlayer:: SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) { //-------------------------------------------------------------------- // This function sets the "window" //-------------------------------------------------------------------- fI=FindClusterIndex(zmin); fZmax=zmax; fYmin=ymin; fYmax=ymax; } const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){ //-------------------------------------------------------------------- // This function returns clusters within the "window" //-------------------------------------------------------------------- const AliITSclusterV2 *cluster=0; for (Int_t i=fI; iGetZ() > fZmax) break; if (c->IsUsed()) continue; const AliITSdetector &det=GetDetector(c->GetDetectorIndex()); Double_t y=fR*det.GetPhi() + c->GetY(); if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi(); if (y>1.*fR*TMath::Pi() && fYmaxfYmax) continue; cluster=c; ci=i; fI=i+1; break; } return cluster; } Int_t AliITStrackerV2::AliITSlayer:: FindDetectorIndex(Double_t phi, Double_t z) const { //-------------------------------------------------------------------- //This function finds the detector crossed by the track //-------------------------------------------------------------------- Double_t dphi=phi-fPhiOffset; if (dphi < 0) dphi += 2*TMath::Pi(); else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi(); Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5); Double_t dz=fZOffset-z; Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5); if (np>=fNladders) np-=fNladders; if (np<0) np+=fNladders; #ifdef DEBUG cout<1) { Double_t xi=0.5*(fLayers[1].GetR()+fLayers[2].GetR()); d+=0.034*xi*xi; } if (fI>3) { Double_t xi=0.5*(fLayers[3].GetR()+fLayers[4].GetR()); d+=0.039*xi*xi; } return d/(xn*xn); } Int_t AliITStrackerV2::AliITSlayer::InRoad() const { //-------------------------------------------------------------------- // This function returns number of clusters within the "window" //-------------------------------------------------------------------- Int_t ncl=0; for (Int_t i=fI; iGetZ() > fZmax) break; //if (c->IsUsed()) continue; const AliITSdetector &det=GetDetector(c->GetDetectorIndex()); Double_t y=fR*det.GetPhi() + c->GetY(); if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi(); if (y>1.*fR*TMath::Pi() && fYmaxfYmax) continue; ncl++; } return ncl; }