/************************************************************************** * 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. * **************************************************************************/ /* $Log$ Revision 1.10 2001/05/16 14:57:25 alibrary New files for folders and Stack Revision 1.9 2001/05/11 07:16:56 hristov Fix needed on Sun and Alpha Revision 1.8 2001/05/08 15:00:15 hristov Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov) Revision 1.5 2000/12/20 07:51:59 kowal2 Changes suggested by Alessandra and Paolo to avoid overlapped data fields in encapsulated classes. Revision 1.4 2000/11/02 07:27:16 kowal2 code corrections Revision 1.2 2000/06/30 12:07:50 kowal2 Updated from the TPC-PreRelease branch Revision 1.1.2.1 2000/06/25 08:53:55 kowal2 Splitted from AliTPCtracking */ //------------------------------------------------------- // Implementation of the TPC tracker // // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch //------------------------------------------------------- #include #include #include #include #include "AliTPCtracker.h" #include "AliTPCcluster.h" #include "AliTPCParam.h" #include "AliTPCClustersRow.h" //_____________________________________________________________________________ AliTPCtracker::AliTPCtracker(const AliTPCParam *par): fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2) { //--------------------------------------------------------------------- // The main TPC tracker constructor //--------------------------------------------------------------------- fInnerSec=new AliTPCSector[fkNIS]; fOuterSec=new AliTPCSector[fkNOS]; Int_t i; for (i=0; iGetEntries()); for (Int_t i=0; iAdjustSectorRow(s->GetID(),sec,row); if (secGetArray()->GetEntriesFast(); while (ncl--) { AliTPCcluster *c=(AliTPCcluster*)(*clrow)[ncl]; index=(((sec<<8)+row)<<16)+ncl; fOuterSec[(sec-fkNIS*2)%fkNOS][row].InsertCluster(c,index); } } fN=fkNOS; fSectors=fOuterSec; } //_____________________________________________________________________________ void AliTPCtracker::UnloadOuterSectors() { //----------------------------------------------------------------- // This function clears outer TPC sectors. //----------------------------------------------------------------- Int_t nup=fOuterSec->GetNRows(); for (Int_t i=0; iGetEntries()); for (Int_t i=0; iAdjustSectorRow(s->GetID(),sec,row); if (sec>=fkNIS*2) continue; AliTPCClustersRow *clrow=fClustersArray.GetRow(sec,row); Int_t ncl=clrow->GetArray()->GetEntriesFast(); while (ncl--) { AliTPCcluster *c=(AliTPCcluster*)(*clrow)[ncl]; index=(((sec<<8)+row)<<16)+ncl; fInnerSec[sec%fkNIS][row].InsertCluster(c,index); } } fN=fkNIS; fSectors=fInnerSec; } //_____________________________________________________________________________ void AliTPCtracker::UnloadInnerSectors() { //----------------------------------------------------------------- // This function clears inner TPC sectors. //----------------------------------------------------------------- Int_t nlow=fInnerSec->GetNRows(); for (Int_t i=0; iGetNRows()); Int_t tryAgain=kSKIP; Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift(); if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi(); if (alpha < 0. ) alpha += 2.*TMath::Pi(); Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN; for (Int_t nr=fSectors->GetRowNumber(xt)-1; nr>=rf; nr--) { Double_t x=fSectors->GetX(nr), ymax=fSectors->GetMaxY(nr); if (!t.PropagateTo(x)) return 0; AliTPCcluster *cl=0; UInt_t index=0; Double_t maxchi2=kMaxCHI2; const AliTPCRow &krow=fSectors[s][nr]; Double_t pt=t.GetConvConst()/(100/0.299792458/0.2)/t.Get1Pt(); Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),pt); Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl()); Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ(); if (road>kMaxROAD) { if (t.GetNumberOfClusters()>4) cerr<GetY() > y+road) break; if (c->IsUsed()) continue; if ((c->GetZ()-z)*(c->GetZ()-z) > 16.*(t.GetSigmaZ2()+sz2)) continue; Double_t chi2=t.GetPredictedChi2(c); if (chi2 > maxchi2) continue; maxchi2=chi2; cl=c; index=krow.GetIndex(i); } } if (cl) { Float_t l=fSectors->GetPadPitchWidth(); t.SetSampledEdx(cl->GetQ()/l,t.GetNumberOfClusters()); if (!t.Update(cl,maxchi2,index)) { if (!tryAgain--) return 0; } else tryAgain=kSKIP; } else { if (tryAgain==0) break; if (y > ymax) { s = (s+1) % fN; if (!t.Rotate(fSectors->GetAlpha())) return 0; } else if (y <-ymax) { s = (s-1+fN) % fN; if (!t.Rotate(-fSectors->GetAlpha())) return 0; } tryAgain--; } } return 1; } //_____________________________________________________________________________ Int_t AliTPCtracker::FollowBackProlongation (AliTPCseed& seed, const AliTPCtrack &track) { //----------------------------------------------------------------- // This function propagates tracks back through the TPC //----------------------------------------------------------------- Double_t alpha=seed.GetAlpha() - fSectors->GetAlphaShift(); if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi(); if (alpha < 0. ) alpha += 2.*TMath::Pi(); Int_t s=Int_t(alpha/fSectors->GetAlpha())%fN; Int_t idx=-1, sec=-1, row=-1; Int_t nc=seed.GetLabel(); //index of the cluster to start with if (nc--) { idx=track.GetClusterIndex(nc); sec=(idx&0xff000000)>>24; row=(idx&0x00ff0000)>>16; } if (fSectors==fInnerSec) { if (sec >= 2*fkNIS) row=-1; } else { if (sec < 2*fkNIS) row=-1; } Int_t nr=fSectors->GetNRows(); for (Int_t i=0; iGetX(i), ymax=fSectors->GetMaxY(i); if (!seed.PropagateTo(x)) return 0; Double_t y=seed.GetY(); if (y > ymax) { s = (s+1) % fN; if (!seed.Rotate(fSectors->GetAlpha())) return 0; } else if (y <-ymax) { s = (s-1+fN) % fN; if (!seed.Rotate(-fSectors->GetAlpha())) return 0; } AliTPCcluster *cl=0; Int_t index=0; Double_t maxchi2=kMaxCHI2; Double_t pt=seed.GetConvConst()/(100/0.299792458/0.2)/seed.Get1Pt(); Double_t sy2=SigmaY2(seed.GetX(),seed.GetTgl(),pt); Double_t sz2=SigmaZ2(seed.GetX(),seed.GetTgl()); Double_t road=4.*sqrt(seed.GetSigmaY2() + sy2), z=seed.GetZ(); if (road>kMaxROAD) { cerr<>24; row=(idx&0x00ff0000)>>16; } if (fSectors==fInnerSec) { if (sec >= 2*fkNIS) row=-1; } else { if (sec < 2*fkNIS) row=-1; } } if (!cl) { //try to fill the gap const AliTPCRow &krow=fSectors[s][i]; if (accepted>27) if (krow) { for (Int_t i=krow.Find(y-road); iGetY() > y+road) break; if (c->IsUsed()) continue; if ((c->GetZ()-z)*(c->GetZ()-z)>16.*(seed.GetSigmaZ2()+sz2)) continue; Double_t chi2=seed.GetPredictedChi2(c); if (chi2 > maxchi2) continue; maxchi2=chi2; cl=c; index=krow.GetIndex(i); } } } if (cl) { Float_t l=fSectors->GetPadPitchWidth(); seed.SetSampledEdx(cl->GetQ()/l,seed.GetNumberOfClusters()); seed.Update(cl,maxchi2,index); } } seed.SetLabel(nc); return 1; } //_____________________________________________________________________________ void AliTPCtracker::MakeSeeds(Int_t i1, Int_t i2) { //----------------------------------------------------------------- // This function creates track seeds. //----------------------------------------------------------------- if (fSeeds==0) fSeeds=new TObjArray(15000); Double_t x[5], c[15]; Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift(); Double_t cs=cos(alpha), sn=sin(alpha); Double_t x1 =fOuterSec->GetX(i1); Double_t xx2=fOuterSec->GetX(i2); for (Int_t ns=0; nsGetY(), z1=kr1[is]->GetZ(); for (Int_t js=0; js < nl+nm+nu; js++) { const AliTPCcluster *kcl; Double_t x2, y2, z2; Double_t x3=0.,y3=0.; if (jsGetY(); z2=kcl->GetZ(); x2= xx2*cs+y2*sn; y2=-xx2*sn+y2*cs; } else if (jsGetY(); z2=kcl->GetZ(); } else { const AliTPCRow& kr2=fOuterSec[(ns+1)%fkNOS][i2]; kcl=kr2[js-nl-nm]; y2=kcl->GetY(); z2=kcl->GetZ(); x2=xx2*cs-y2*sn; y2=xx2*sn+y2*cs; } Double_t zz=z1 - z1/x1*(x1-x2); if (TMath::Abs(zz-z2)>5.) continue; Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1); if (d==0.) {cerr<<"MakeSeeds warning: Straight seed !\n"; continue;} x[0]=y1; x[1]=z1; x[4]=f1(x1,y1,x2,y2,x3,y3); if (TMath::Abs(x[4]) >= 0.0066) continue; x[2]=f2(x1,y1,x2,y2,x3,y3); //if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue; x[3]=f3(x1,y1,x2,y2,z1,z2); if (TMath::Abs(x[3]) > 1.2) continue; Double_t a=asin(x[2]); Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2])); if (TMath::Abs(zv)>10.) continue; Double_t sy1=kr1[is]->GetSigmaY2(), sz1=kr1[is]->GetSigmaZ2(); Double_t sy2=kcl->GetSigmaY2(), sz2=kcl->GetSigmaZ2(); Double_t sy3=100*0.025, sy=0.1, sz=0.1; Double_t f40=(f1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy; Double_t f42=(f1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy; Double_t f43=(f1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy; Double_t f20=(f2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy; Double_t f22=(f2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy; Double_t f23=(f2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy; Double_t f30=(f3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy; Double_t f31=(f3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz; Double_t f32=(f3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy; Double_t f34=(f3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz; c[0]=sy1; c[1]=0.; c[2]=sz1; c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23; c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22; c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34; c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23; c[13]=f30*sy1*f40+f32*sy2*f42; c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43; UInt_t index=kr1.GetIndex(is); AliTPCseed *track=new AliTPCseed(index, x, c, x1, ns*alpha+shift); Float_t l=fOuterSec->GetPadPitchWidth(); track->SetSampledEdx(kr1[is]->GetQ()/l,0); Int_t rc=FollowProlongation(*track, i2); if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track; else fSeeds->AddLast(track); } } } } //_____________________________________________________________________________ Int_t AliTPCtracker::ReadSeeds(const TFile *inp) { //----------------------------------------------------------------- // This function reades track seeds. //----------------------------------------------------------------- TDirectory *savedir=gDirectory; TFile *in=(TFile*)inp; if (!in->IsOpen()) { cerr<<"AliTPCtracker::ReadSeeds(): input file is not open !\n"; return 1; } in->cd(); TTree *seedTree=(TTree*)in->Get("Seeds"); if (!seedTree) { cerr<<"AliTPCtracker::ReadSeeds(): "; cerr<<"can't get a tree with track seeds !\n"; return 2; } AliTPCtrack *seed=new AliTPCtrack; seedTree->SetBranchAddress("tracks",&seed); if (fSeeds==0) fSeeds=new TObjArray(15000); Int_t n=(Int_t)seedTree->GetEntries(); for (Int_t i=0; iGetEvent(i); fSeeds->AddLast(new AliTPCseed(*seed,seed->GetAlpha())); } delete seed; savedir->cd(); return 0; } //_____________________________________________________________________________ Int_t AliTPCtracker::Clusters2Tracks(const TFile *inp, TFile *out) { //----------------------------------------------------------------- // This is a track finder. //----------------------------------------------------------------- TDirectory *savedir=gDirectory; if (inp) { TFile *in=(TFile*)inp; if (!in->IsOpen()) { cerr<<"AliTPCtracker::Clusters2Tracks(): input file is not open !\n"; return 1; } } if (!out->IsOpen()) { cerr<<"AliTPCtracker::Clusters2Tracks(): output file is not open !\n"; return 2; } out->cd(); TTree tracktree("TPCf","Tree with TPC tracks"); AliTPCtrack *iotrack=0; tracktree.Branch("tracks","AliTPCtrack",&iotrack,32000,0); LoadOuterSectors(); //find track seeds Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows(); Int_t nrows=nlow+nup; if (fSeeds==0) { Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap); MakeSeeds(nup-1, nup-1-gap); MakeSeeds(nup-1-shift, nup-1-shift-gap); } fSeeds->Sort(); //tracking in outer sectors Int_t nseed=fSeeds->GetEntriesFast(); Int_t i; for (i=0; iUncheckedAt(i), &t=*pt; if (FollowProlongation(t)) { UseClusters(&t); continue; } delete fSeeds->RemoveAt(i); } UnloadOuterSectors(); //tracking in inner sectors LoadInnerSectors(); Int_t found=0; for (i=0; iUncheckedAt(i), &t=*pt; if (!pt) continue; Int_t nc=t.GetNumberOfClusters(); Double_t alpha=t.GetAlpha() - fInnerSec->GetAlphaShift(); if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi(); if (alpha < 0. ) alpha += 2.*TMath::Pi(); Int_t ns=Int_t(alpha/fInnerSec->GetAlpha())%fkNIS; alpha=ns*fInnerSec->GetAlpha() + fInnerSec->GetAlphaShift() - t.GetAlpha(); if (t.Rotate(alpha)) { if (FollowProlongation(t)) { if (t.GetNumberOfClusters() >= Int_t(0.4*nrows)) { t.CookdEdx(); iotrack=pt; tracktree.Fill(); UseClusters(&t,nc); cerr<RemoveAt(i); } UnloadInnerSectors(); tracktree.Write(); cerr<<"Number of found tracks : "<cd(); return 0; } //_____________________________________________________________________________ Int_t AliTPCtracker::PropagateBack(const TFile *inp, TFile *out) { //----------------------------------------------------------------- // This function propagates tracks back through the TPC. //----------------------------------------------------------------- fSeeds=new TObjArray(15000); TFile *in=(TFile*)inp; TDirectory *savedir=gDirectory; if (!in->IsOpen()) { cerr<<"AliTPCtracker::PropagateBack(): "; cerr<<"file with back propagated ITS tracks is not open !\n"; return 1; } if (!out->IsOpen()) { cerr<<"AliTPCtracker::PropagateBack(): "; cerr<<"file for back propagated TPC tracks is not open !\n"; return 2; } in->cd(); TTree *bckTree=(TTree*)in->Get("ITSb"); if (!bckTree) { cerr<<"AliTPCtracker::PropagateBack() "; cerr<<"can't get a tree with back propagated ITS tracks !\n"; return 3; } AliTPCtrack *bckTrack=new AliTPCtrack; bckTree->SetBranchAddress("tracks",&bckTrack); TTree *tpcTree=(TTree*)in->Get("TPCf"); if (!tpcTree) { cerr<<"AliTPCtracker::PropagateBack() "; cerr<<"can't get a tree with TPC tracks !\n"; return 4; } AliTPCtrack *tpcTrack=new AliTPCtrack; tpcTree->SetBranchAddress("tracks",&tpcTrack); //*** Prepare an array of tracks to be back propagated Int_t nup=fOuterSec->GetNRows(), nlow=fInnerSec->GetNRows(); Int_t nrows=nlow+nup; TObjArray tracks(15000); Int_t i=0,j=0; Int_t tpcN=(Int_t)tpcTree->GetEntries(); Int_t bckN=(Int_t)bckTree->GetEntries(); if (jGetEvent(j++); for (i=0; iGetEvent(i); Double_t alpha=tpcTrack->GetAlpha(); if (jGetLabel())==TMath::Abs(bckTrack->GetLabel())) { if (!bckTrack->Rotate(alpha-bckTrack->GetAlpha())) continue; fSeeds->AddLast(new AliTPCseed(*bckTrack,bckTrack->GetAlpha())); bckTree->GetEvent(j++); } else { tpcTrack->ResetCovariance(); fSeeds->AddLast(new AliTPCseed(*tpcTrack,alpha)); } tracks.AddLast(new AliTPCtrack(*tpcTrack)); } out->cd(); TTree backTree("TPCb","Tree with back propagated TPC tracks"); AliTPCtrack *otrack=0; backTree.Branch("tracks","AliTPCtrack",&otrack,32000,0); //*** Back propagation through inner sectors LoadInnerSectors(); Int_t nseed=fSeeds->GetEntriesFast(); for (i=0; iUncheckedAt(i), &s=*ps; const AliTPCtrack *pt=(AliTPCtrack*)tracks.UncheckedAt(i), &t=*pt; Int_t nc=t.GetNumberOfClusters(); s.SetLabel(nc-1); //set number of the cluster to start with if (FollowBackProlongation(s,t)) { UseClusters(&s); continue; } delete fSeeds->RemoveAt(i); } UnloadInnerSectors(); //*** Back propagation through outer sectors LoadOuterSectors(); Int_t found=0; for (i=0; iUncheckedAt(i), &s=*ps; if (!ps) continue; Int_t nc=s.GetNumberOfClusters(); const AliTPCtrack *pt=(AliTPCtrack*)tracks.UncheckedAt(i), &t=*pt; Double_t alpha=s.GetAlpha() - fSectors->GetAlphaShift(); if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi(); if (alpha < 0. ) alpha += 2.*TMath::Pi(); Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN; alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift(); alpha-=s.GetAlpha(); if (s.Rotate(alpha)) { if (FollowBackProlongation(s,t)) { if (s.GetNumberOfClusters() >= Int_t(0.4*nrows)) { s.CookdEdx(); s.SetLabel(t.GetLabel()); UseClusters(&s,nc); otrack=ps; backTree.Fill(); cerr<RemoveAt(i); } UnloadOuterSectors(); backTree.Write(); savedir->cd(); cerr<<"Number of seeds: "<>24; Int_t row=(index&0x00ff0000)>>16; Int_t ncl=(index&0x0000ffff)>>00; AliTPCClustersRow *clrow=((AliTPCtracker *) this)->fClustersArray.GetRow(sec,row); return (AliCluster*)(*clrow)[ncl]; } //__________________________________________________________________________ void AliTPCtracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const { //-------------------------------------------------------------------- //This function "cooks" a track label. If label<0, this track is fake. //-------------------------------------------------------------------- Int_t noc=t->GetNumberOfClusters(); Int_t *lb=new Int_t[noc]; Int_t *mx=new Int_t[noc]; AliCluster **clusters=new AliCluster*[noc]; Int_t i; for (i=0; iGetClusterIndex(i); clusters[i]=GetCluster(index); } Int_t lab=123456789; for (i=0; iGetLabel(0)); Int_t j; for (j=0; jmax) {max=mx[i]; lab=lb[i];} for (i=0; iGetLabel(1)) == lab || TMath::Abs(c->GetLabel(2)) == lab ) max++; } if ((1.- Float_t(max)/noc) > wrong) lab=-lab; else { Int_t tail=Int_t(0.10*noc); max=0; for (i=1; i<=tail; i++) { AliCluster *c=clusters[noc-i]; if (lab == TMath::Abs(c->GetLabel(0)) || lab == TMath::Abs(c->GetLabel(1)) || lab == TMath::Abs(c->GetLabel(2))) max++; } if (max < Int_t(0.5*tail)) lab=-lab; } t->SetLabel(lab); delete[] lb; delete[] mx; delete[] clusters; } //_________________________________________________________________________ void AliTPCtracker::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) { //----------------------------------------------------------------------- // Setup inner sector //----------------------------------------------------------------------- if (f==0) { fAlpha=par->GetInnerAngle(); fAlphaShift=par->GetInnerAngleShift(); fPadPitchWidth=par->GetInnerPadPitchWidth(); fPadPitchLength=par->GetInnerPadPitchLength(); fN=par->GetNRowLow(); fRow=new AliTPCRow[fN]; for (Int_t i=0; iGetPadRowRadiiLow(i)); } else { fAlpha=par->GetOuterAngle(); fAlphaShift=par->GetOuterAngleShift(); fPadPitchWidth=par->GetOuterPadPitchWidth(); fPadPitchLength=par->GetOuterPadPitchLength(); fN=par->GetNRowUp(); fRow=new AliTPCRow[fN]; for (Int_t i=0; iGetPadRowRadiiUp(i)); } } //_________________________________________________________________________ void AliTPCtracker::AliTPCRow::InsertCluster(const AliTPCcluster* c, UInt_t index) { //----------------------------------------------------------------------- // Insert a cluster into this pad row in accordence with its y-coordinate //----------------------------------------------------------------------- if (fN==kMaxClusterPerRow) { cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return; } if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;} Int_t i=Find(c->GetY()); memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCcluster*)); memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t)); fIndex[i]=index; fClusters[i]=c; fN++; } //___________________________________________________________________ Int_t AliTPCtracker::AliTPCRow::Find(Double_t y) const { //----------------------------------------------------------------------- // Return the index of the nearest cluster //----------------------------------------------------------------------- if (y <= fClusters[0]->GetY()) return 0; if (y > fClusters[fN-1]->GetY()) return fN; Int_t b=0, e=fN-1, m=(b+e)/2; for (; b fClusters[m]->GetY()) b=m+1; else e=m; } return m; } //_____________________________________________________________________________ void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) { //----------------------------------------------------------------- // This funtion calculates dE/dX within the "low" and "up" cuts. //----------------------------------------------------------------- Int_t i; Int_t nc=GetNumberOfClusters(); Int_t swap;//stupid sorting do { swap=0; for (i=0; i