ROOT logo
/**************************************************************************
 * 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.                  *
 **************************************************************************/

/* $Id$ */

//-----------------------------------------------------------------
//           Implementation of the TPC track class
//        This class is used by the AliTPCtracker class
//      Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
//-----------------------------------------------------------------

#include <Riostream.h>

#include "AliTPCtrack.h"
#include "AliCluster.h"
#include "AliTracker.h"
#include "AliESDtrack.h"
#include "AliESDVertex.h"
#include "TTreeStream.h"
#include  "AliTPCRecoParam.h"
#include  "AliTPCReconstructor.h"
ClassImp(AliTPCtrack)

//_________________________________________________________________________
AliTPCtrack::AliTPCtrack(): 
  AliKalmanTrack(),
  fdEdx(0),
  fSdEdx(1e10),
  fNFoundable(0),
  fBConstrain(kFALSE),
  fLastPoint(-1),
  fFirstPoint(-1),
  fRemoval(0),
  fTrackType(0),
  fLab2(-1),
  fNShared(0),
  fReference()
{
  //-------------------------------------------------
  // default constructor
  //-------------------------------------------------
  for (Int_t i=0; i<kMaxRow;i++) fIndex[i]=-2;
  for (Int_t i=0; i<4;i++) fPoints[i]=0.;
  for (Int_t i=0; i<12;i++) fKinkPoint[i]=0.;
  for (Int_t i=0; i<3;i++) fKinkIndexes[i]=0;
  for (Int_t i=0; i<3;i++) fV0Indexes[i]=0;
}

//_________________________________________________________________________



AliTPCtrack::AliTPCtrack(Double_t x, Double_t alpha, const Double_t p[5],
		         const Double_t cov[15], Int_t index) :
  AliKalmanTrack(),
  fdEdx(0),
  fSdEdx(1e10),
  fNFoundable(0),
  fBConstrain(kFALSE),
  fLastPoint(0),
  fFirstPoint(0),
  fRemoval(0),
  fTrackType(0),
  fLab2(0),
  fNShared(0),
  fReference()
{
  //-----------------------------------------------------------------
  // This is the main track constructor.
  //-----------------------------------------------------------------
  Double_t cnv=1./(AliTracker::GetBz()*kB2C);

  Double_t pp[5]={
    p[0],
    p[1],
    x*p[4] - p[2],
    p[3],
    p[4]*cnv
  };

  Double_t c22 = x*x*cov[14] - 2*x*cov[12] + cov[5];
  Double_t c32 = x*cov[13] - cov[8];
  Double_t c20 = x*cov[10] - cov[3], 
           c21 = x*cov[11] - cov[4], c42 = x*cov[14] - cov[12];

  Double_t cc[15]={
    cov[0 ],
    cov[1 ],     cov[2 ],
    c20,         c21,         c22,
    cov[6 ],     cov[7 ],     c32,     cov[9 ],
    cov[10]*cnv, cov[11]*cnv, c42*cnv, cov[13]*cnv, cov[14]*cnv*cnv
  };

  Double_t mostProbablePt=AliExternalTrackParam::GetMostProbablePt();
  Double_t p0=TMath::Sign(1/mostProbablePt,pp[4]);
  Double_t w0=cc[14]/(cc[14] + p0*p0), w1=p0*p0/(cc[14] + p0*p0);
  pp[4] = w0*p0 + w1*pp[4]; 
  cc[10]*=w1; cc[11]*=w1; cc[12]*=w1; cc[13]*=w1; cc[14]*=w1;

  Set(x,alpha,pp,cc);

  SetNumberOfClusters(1);
  
  fIndex[0]=index;
  for (Int_t i=1; i<kMaxRow;i++) fIndex[i]=-2;
  for (Int_t i=0; i<4;i++) fPoints[i]=0.;
  for (Int_t i=0; i<12;i++) fKinkPoint[i]=0.;
  for (Int_t i=0; i<3;i++) fKinkIndexes[i]=0;
  for (Int_t i=0; i<3;i++) fV0Indexes[i]=0;
}

//_____________________________________________________________________________
AliTPCtrack::AliTPCtrack(const AliESDtrack& t, TTreeSRedirector *pcstream) :
  AliKalmanTrack(),
  fdEdx(t.GetTPCsignal()),
  fSdEdx(1e10),
  fNFoundable(0),
  fBConstrain(kFALSE),
  fLastPoint(0),
  fFirstPoint(0),
  fRemoval(0),
  fTrackType(0),
  fLab2(0),
  fNShared(0),
  fReference()
{
  //-----------------------------------------------------------------
  // Conversion AliESDtrack -> AliTPCtrack.
  //-----------------------------------------------------------------
  const Double_t kmaxC[4]={10,10,0.1,0.1};  // cuts on the rms /fP0,fP1,fP2,fP3
  SetNumberOfClusters(t.GetTPCclusters(fIndex));
  SetLabel(t.GetLabel());
  SetMass(t.GetMassForTracking());
  for (Int_t i=0; i<4;i++) fPoints[i]=0.;
  for (Int_t i=0; i<12;i++) fKinkPoint[i]=0.;
  for (Int_t i=0; i<3;i++) fKinkIndexes[i]=0;
  for (Int_t i=0; i<3;i++) fV0Indexes[i]=0;
  //
  // choose parameters to start
  //
  const AliTPCRecoParam * recoParam = AliTPCReconstructor::GetRecoParam();
  Int_t reject=0;
  AliExternalTrackParam param(t);

  const AliExternalTrackParam  *tpcout=(t.GetFriendTrack())? ((AliESDfriendTrack*)(t.GetFriendTrack()))->GetTPCOut():0;
  const AliExternalTrackParam  *tpcin = t.GetInnerParam();
  const AliExternalTrackParam  *tpc=(tpcout)?tpcout:tpcin;
  Bool_t isBackProp = tpcout==0; // is this backpropagation?
  if (!tpc) tpc=&param;

  Bool_t isOK = (recoParam->GetUseOuterDetectors() && t.IsOn(AliESDtrack::kTRDrefit)) || isBackProp;
  if (param.GetCovariance()[0]>kmaxC[0]*kmaxC[0] ||
      param.GetCovariance()[2]>kmaxC[1]*kmaxC[1] ||
      param.GetCovariance()[5]>kmaxC[2]*kmaxC[2] ||
      param.GetCovariance()[9]>kmaxC[3]*kmaxC[3]) isOK=kFALSE;
  //
  if (isOK) isOK &= param.Rotate(tpc->GetAlpha()); // using external seed
  Double_t oldX=param.GetX(),  oldY=param.GetY(),  oldZ=param.GetZ();
  if (!isOK ){
    param=*tpc;
    isOK=kTRUE;
    reject=1;
  }
  else { // using external seed
    //  param.Rotate(tpc->GetAlpha()); // not needed
    if (!AliTracker::PropagateTrackToBxByBz(&param,tpc->GetX(),GetMass(),2.,kFALSE) ||
	param.GetCovariance()[0]>kmaxC[0]*kmaxC[0] ||
	param.GetCovariance()[2]>kmaxC[1]*kmaxC[1] ||
	param.GetCovariance()[5]>kmaxC[2]*kmaxC[2] ||
	param.GetCovariance()[9]>kmaxC[3]*kmaxC[3]) isOK=kFALSE;
  }
  if (isOK) {
    Double_t chi2= param.GetPredictedChi2(tpc);
    if (isBackProp) {
      if (chi2>recoParam->GetMaxChi2TPCITS()) isOK=kFALSE; // protection against outliers in the ITS
    }
    else if (chi2>recoParam->GetMaxChi2TPCTRD()) isOK=kFALSE; // protection against outliers in the TRD
  }

  if (!isOK){
    param=*tpc;
    isOK=kTRUE;
    reject=2;
  }
  if (reject>0){
    param.ResetCovariance(4.);  // reset covariance if start from backup param
  }
  //
  //
  if (pcstream){
    AliExternalTrackParam dummy;
    AliExternalTrackParam *ptpc=(AliExternalTrackParam *)tpc;
    //    if (!ptpc) ptpc=&dummy;
    AliESDtrack *esd= (AliESDtrack *)&t;
    (*pcstream)<<"trackP"<<
      "reject="<<reject<<   // flag - rejection of current esd track parameters
      "esd.="<<esd<<        // original esd track
      "tr.="<<&param<<      // starting track parameters
      "out.="<<ptpc<<       // backup tpc parameters
      "\n";
  }

  Set(param.GetX(),param.GetAlpha(),param.GetParameter(),param.GetCovariance());

  if ((t.GetStatus()&AliESDtrack::kTIME) == 0) return;
  StartTimeIntegral();
  Double_t times[AliPID::kSPECIESC]; 
  t.GetIntegratedTimes(times,AliPID::kSPECIESC); 
  SetIntegratedTimes(times);
  SetIntegratedLength(t.GetIntegratedLength());

  if (GetX()>oldX) {
     Double_t dX=GetX()-oldX, dY=GetY()-oldY, dZ=GetZ()-oldZ;
     Double_t d=TMath::Sqrt(dX*dX + dY*dY + dZ*dZ);
     AddTimeStep(d);
  }
}

//_____________________________________________________________________________
AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) :
  AliKalmanTrack(t),
  fdEdx(t.fdEdx),
  fSdEdx(t.fSdEdx),
  fNFoundable(t.fNFoundable),
  fBConstrain(t.fBConstrain),
  fLastPoint(t.fLastPoint),
  fFirstPoint(t.fFirstPoint),
  fRemoval(t.fRemoval),
  fTrackType(t.fTrackType),
  fLab2(t.fLab2),
  fNShared(t.fNShared),
  fReference(t.fReference)

{
  //-----------------------------------------------------------------
  // This is a track copy constructor.
  //-----------------------------------------------------------------
  Set(t.GetX(),t.GetAlpha(),t.GetParameter(),t.GetCovariance());

  for (Int_t i=0; i<kMaxRow; i++) fIndex[i]=t.fIndex[i];
  for (Int_t i=0; i<4;i++) fPoints[i]=t.fPoints[i];
  for (Int_t i=0; i<12;i++) fKinkPoint[i]=t.fKinkPoint[i];
  for (Int_t i=0; i<3;i++) fKinkIndexes[i]=t.fKinkIndexes[i];
  for (Int_t i=0; i<3;i++) fV0Indexes[i]=t.fV0Indexes[i];
}

AliTPCtrack& AliTPCtrack::operator=(const AliTPCtrack& o){
  if(this!=&o){
    AliKalmanTrack::operator=(o);
    fdEdx = o.fdEdx;
    for(Int_t i = 0;i<kMaxRow;++i)fIndex[i] = o.fIndex[i];
    for(Int_t i = 0;i<4;++i)fPoints[i] = o.fPoints[i];
    fSdEdx = o.fSdEdx;
    fNFoundable = o.fNFoundable;
    fBConstrain = o.fBConstrain;
    fLastPoint  = o.fLastPoint;
    fFirstPoint = o.fFirstPoint;
    fTrackType  = o.fTrackType;
    fLab2       = o.fLab2;
    fNShared    = o.fNShared;
    fReference  = o.fReference;
    for(Int_t i = 0;i<12;++i) fKinkPoint[i] = o.fKinkPoint[i];

    for(Int_t i = 0;i<3;++i){
      fKinkIndexes[i] = o.fKinkIndexes[i];
      fV0Indexes[i] = o.fV0Indexes[i];
    }
  }
  return *this;

}


//_____________________________________________________________________________
Int_t AliTPCtrack::Compare(const TObject *o) const {
  //-----------------------------------------------------------------
  // This function compares tracks according to the their curvature
  //-----------------------------------------------------------------
  AliTPCtrack *t=(AliTPCtrack*)o;
  //Double_t co=t->OneOverPt();
  //Double_t c = OneOverPt();
  Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
  Double_t c =GetSigmaY2()*GetSigmaZ2();
  if (c>co) return 1;
  else if (c<co) return -1;
  return 0;
}

Double_t AliTPCtrack::GetPredictedChi2(const AliCluster *c) const {
  //-----------------------------------------------------------------
  // This function calculates a predicted chi2 increment.
  //-----------------------------------------------------------------
  Double_t p[2]={c->GetY(), c->GetZ()};
  Double_t cov[3]={c->GetSigmaY2(), 0., c->GetSigmaZ2()};
  return AliExternalTrackParam::GetPredictedChi2(p,cov);
}

//_____________________________________________________________________________
Bool_t AliTPCtrack::PropagateTo(Double_t xk,Double_t rho,Double_t x0) {
  //-----------------------------------------------------------------
  //  This function propagates a track to a reference plane x=xk.
  //  rho - density of the crossed matrial (g/cm^3)
  //  x0  - radiation length of the crossed material (g/cm^2) 
  //-----------------------------------------------------------------
  //
  Double_t bz=GetBz();
  Double_t zat=0;
  if (!GetZAt(xk, bz,zat)) return kFALSE;
  if (TMath::Abs(zat)>250.){
    // Don't propagate track outside of the fiducial volume - material budget not proper one
    //
    //AliWarning("Propagate outside of fiducial volume");
    return kFALSE;
  }

  Double_t oldX=GetX(), oldY=GetY(), oldZ=GetZ();
  //if (!AliExternalTrackParam::PropagateTo(xk,bz)) return kFALSE;
  Double_t b[3]; GetBxByBz(b);
  if (!AliExternalTrackParam::PropagateToBxByBz(xk,b)) return kFALSE;

  Double_t d = TMath::Sqrt((GetX()-oldX)*(GetX()-oldX) + 
                           (GetY()-oldY)*(GetY()-oldY) + 
                           (GetZ()-oldZ)*(GetZ()-oldZ));
  if (IsStartedTimeIntegral() && GetX()>oldX) AddTimeStep(d);

  if (oldX < xk) d = -d;
  if (!AliExternalTrackParam::CorrectForMeanMaterial(d*rho/x0,d*rho,GetMass(),
      kFALSE,AliExternalTrackParam::BetheBlochGas)) return kFALSE;

  return kTRUE;
}

//_____________________________________________________________________________
Bool_t 
AliTPCtrack::PropagateToVertex(const AliESDVertex *v,Double_t rho,Double_t x0) 
{
  //-----------------------------------------------------------------
  // This function propagates tracks to the vertex
  // rho - density of the crossed matrial (g/cm3)
  // x0  - radiation length of the crossed material (g/cm2) 
  //-----------------------------------------------------------------
  Double_t oldX=GetX(), oldY=GetY(), oldZ=GetZ();

  //Double_t bz=GetBz();
  //if (!PropagateToDCA(v,bz,kVeryBig)) return kFALSE;
  Double_t b[3]; GetBxByBz(b);
  if (!PropagateToDCABxByBz(v,b,kVeryBig)) return kFALSE;

  Double_t d = TMath::Sqrt((GetX()-oldX)*(GetX()-oldX) + 
                           (GetY()-oldY)*(GetY()-oldY) + 
                           (GetZ()-oldZ)*(GetZ()-oldZ));

  if (oldX < GetX()) d = -d;
  if (!AliExternalTrackParam::CorrectForMeanMaterial(d*rho/x0,d*rho,GetMass(),
      kFALSE,AliExternalTrackParam::BetheBlochGas)) return kFALSE;

  return kTRUE;
}

//_____________________________________________________________________________
Bool_t AliTPCtrack::Update(const AliCluster *c, Double_t chisq, Int_t index) {
  //-----------------------------------------------------------------
  // This function associates a cluster with this track.
  //-----------------------------------------------------------------
  Double_t p[2]={c->GetY(), c->GetZ()};
  Double_t cov[3]={c->GetSigmaY2(), 0., c->GetSigmaZ2()};

  if (!AliExternalTrackParam::Update(p,cov)) return kFALSE;

  AliTracker::FillResiduals(this,p,cov,c->GetVolumeId());

  Int_t n=GetNumberOfClusters();
  fIndex[n]=index;
  SetNumberOfClusters(n+1);
  SetChi2(GetChi2()+chisq);

  return kTRUE;
}

////////////////////////////////////////////////////////////////////////
// MI ADDITION

Float_t AliTPCtrack::Density(Int_t row0, Int_t row1)
{
  //
  // calculate cluster density
  Int_t good  = 0;
  Int_t found = 0;
  //if (row0<fFirstPoint) row0 = fFirstPoint;
  if (row1>fLastPoint) row1 = fLastPoint;

  
  for (Int_t i=row0;i<=row1;i++){ 
    //    Int_t index = fClusterIndex[i];
    Int_t index = fIndex[i];
    if (index!=-1)  good++;
    if (index>0)    found++;
  }
  Float_t density=0;
  if (good>0) density = Float_t(found)/Float_t(good);
  return density;
}


Float_t AliTPCtrack::Density2(Int_t row0, Int_t row1)
{
  //
  // calculate cluster density
  Int_t good  = 0;
  Int_t found = 0;
  //  
  for (Int_t i=row0;i<=row1;i++){     
    Int_t index = fIndex[i];
    if (index!=-1)  good++;
    if (index>0)    found++;
  }
  Float_t density=0;
  if (good>0) density = Float_t(found)/Float_t(good);
  return density;
}

void  AliTPCtrack::UpdatePoints()
{
  //--------------------------------------------------
  //calculates first ,amx dens and last points
  //--------------------------------------------------
  Float_t density[160];
  for (Int_t i=0;i<160;i++) density[i]=-1.;
  fPoints[0]= 160;
  fPoints[1] = -1;
  //
  Int_t ngood=0;
  Int_t undeff=0;
  Int_t nall =0;
  Int_t range=20;
  for (Int_t i=0;i<160;i++){
    Int_t last = i-range;
    if (nall<range) nall++;
    if (last>=0){
      if (fIndex[last]>0&& (fIndex[last]&0x8000)==0) ngood--;
      if (fIndex[last]==-1) undeff--;
    }
    if (fIndex[i]>0&& (fIndex[i]&0x8000)==0)   ngood++;
    if (fIndex[i]==-1) undeff++;
    if (nall==range &&undeff<range/2) density[i-range/2] = Float_t(ngood)/Float_t(nall-undeff);
  }
  Float_t maxdens=0;
  Int_t indexmax =0;
  for (Int_t i=0;i<160;i++){
    if (density[i]<0) continue;
    if (density[i]>maxdens){
      maxdens=density[i];
      indexmax=i;
    }
  }
  //
  //max dens point
  fPoints[3] = maxdens;
  fPoints[1] = indexmax;
  //
  // last point
  for (Int_t i=indexmax;i<160;i++){
    if (density[i]<0) continue;
    if (density[i]<maxdens/2.) {
      break;
    }
    fPoints[2]=i;
  }
  //
  // first point
  for (Int_t i=indexmax;i>0;i--){
    if (density[i]<0) continue;
    if (density[i]<maxdens/2.) {
      break;
    }
    fPoints[0]=i;
  }
  //
}

 AliTPCtrack.cxx:1
 AliTPCtrack.cxx:2
 AliTPCtrack.cxx:3
 AliTPCtrack.cxx:4
 AliTPCtrack.cxx:5
 AliTPCtrack.cxx:6
 AliTPCtrack.cxx:7
 AliTPCtrack.cxx:8
 AliTPCtrack.cxx:9
 AliTPCtrack.cxx:10
 AliTPCtrack.cxx:11
 AliTPCtrack.cxx:12
 AliTPCtrack.cxx:13
 AliTPCtrack.cxx:14
 AliTPCtrack.cxx:15
 AliTPCtrack.cxx:16
 AliTPCtrack.cxx:17
 AliTPCtrack.cxx:18
 AliTPCtrack.cxx:19
 AliTPCtrack.cxx:20
 AliTPCtrack.cxx:21
 AliTPCtrack.cxx:22
 AliTPCtrack.cxx:23
 AliTPCtrack.cxx:24
 AliTPCtrack.cxx:25
 AliTPCtrack.cxx:26
 AliTPCtrack.cxx:27
 AliTPCtrack.cxx:28
 AliTPCtrack.cxx:29
 AliTPCtrack.cxx:30
 AliTPCtrack.cxx:31
 AliTPCtrack.cxx:32
 AliTPCtrack.cxx:33
 AliTPCtrack.cxx:34
 AliTPCtrack.cxx:35
 AliTPCtrack.cxx:36
 AliTPCtrack.cxx:37
 AliTPCtrack.cxx:38
 AliTPCtrack.cxx:39
 AliTPCtrack.cxx:40
 AliTPCtrack.cxx:41
 AliTPCtrack.cxx:42
 AliTPCtrack.cxx:43
 AliTPCtrack.cxx:44
 AliTPCtrack.cxx:45
 AliTPCtrack.cxx:46
 AliTPCtrack.cxx:47
 AliTPCtrack.cxx:48
 AliTPCtrack.cxx:49
 AliTPCtrack.cxx:50
 AliTPCtrack.cxx:51
 AliTPCtrack.cxx:52
 AliTPCtrack.cxx:53
 AliTPCtrack.cxx:54
 AliTPCtrack.cxx:55
 AliTPCtrack.cxx:56
 AliTPCtrack.cxx:57
 AliTPCtrack.cxx:58
 AliTPCtrack.cxx:59
 AliTPCtrack.cxx:60
 AliTPCtrack.cxx:61
 AliTPCtrack.cxx:62
 AliTPCtrack.cxx:63
 AliTPCtrack.cxx:64
 AliTPCtrack.cxx:65
 AliTPCtrack.cxx:66
 AliTPCtrack.cxx:67
 AliTPCtrack.cxx:68
 AliTPCtrack.cxx:69
 AliTPCtrack.cxx:70
 AliTPCtrack.cxx:71
 AliTPCtrack.cxx:72
 AliTPCtrack.cxx:73
 AliTPCtrack.cxx:74
 AliTPCtrack.cxx:75
 AliTPCtrack.cxx:76
 AliTPCtrack.cxx:77
 AliTPCtrack.cxx:78
 AliTPCtrack.cxx:79
 AliTPCtrack.cxx:80
 AliTPCtrack.cxx:81
 AliTPCtrack.cxx:82
 AliTPCtrack.cxx:83
 AliTPCtrack.cxx:84
 AliTPCtrack.cxx:85
 AliTPCtrack.cxx:86
 AliTPCtrack.cxx:87
 AliTPCtrack.cxx:88
 AliTPCtrack.cxx:89
 AliTPCtrack.cxx:90
 AliTPCtrack.cxx:91
 AliTPCtrack.cxx:92
 AliTPCtrack.cxx:93
 AliTPCtrack.cxx:94
 AliTPCtrack.cxx:95
 AliTPCtrack.cxx:96
 AliTPCtrack.cxx:97
 AliTPCtrack.cxx:98
 AliTPCtrack.cxx:99
 AliTPCtrack.cxx:100
 AliTPCtrack.cxx:101
 AliTPCtrack.cxx:102
 AliTPCtrack.cxx:103
 AliTPCtrack.cxx:104
 AliTPCtrack.cxx:105
 AliTPCtrack.cxx:106
 AliTPCtrack.cxx:107
 AliTPCtrack.cxx:108
 AliTPCtrack.cxx:109
 AliTPCtrack.cxx:110
 AliTPCtrack.cxx:111
 AliTPCtrack.cxx:112
 AliTPCtrack.cxx:113
 AliTPCtrack.cxx:114
 AliTPCtrack.cxx:115
 AliTPCtrack.cxx:116
 AliTPCtrack.cxx:117
 AliTPCtrack.cxx:118
 AliTPCtrack.cxx:119
 AliTPCtrack.cxx:120
 AliTPCtrack.cxx:121
 AliTPCtrack.cxx:122
 AliTPCtrack.cxx:123
 AliTPCtrack.cxx:124
 AliTPCtrack.cxx:125
 AliTPCtrack.cxx:126
 AliTPCtrack.cxx:127
 AliTPCtrack.cxx:128
 AliTPCtrack.cxx:129
 AliTPCtrack.cxx:130
 AliTPCtrack.cxx:131
 AliTPCtrack.cxx:132
 AliTPCtrack.cxx:133
 AliTPCtrack.cxx:134
 AliTPCtrack.cxx:135
 AliTPCtrack.cxx:136
 AliTPCtrack.cxx:137
 AliTPCtrack.cxx:138
 AliTPCtrack.cxx:139
 AliTPCtrack.cxx:140
 AliTPCtrack.cxx:141
 AliTPCtrack.cxx:142
 AliTPCtrack.cxx:143
 AliTPCtrack.cxx:144
 AliTPCtrack.cxx:145
 AliTPCtrack.cxx:146
 AliTPCtrack.cxx:147
 AliTPCtrack.cxx:148
 AliTPCtrack.cxx:149
 AliTPCtrack.cxx:150
 AliTPCtrack.cxx:151
 AliTPCtrack.cxx:152
 AliTPCtrack.cxx:153
 AliTPCtrack.cxx:154
 AliTPCtrack.cxx:155
 AliTPCtrack.cxx:156
 AliTPCtrack.cxx:157
 AliTPCtrack.cxx:158
 AliTPCtrack.cxx:159
 AliTPCtrack.cxx:160
 AliTPCtrack.cxx:161
 AliTPCtrack.cxx:162
 AliTPCtrack.cxx:163
 AliTPCtrack.cxx:164
 AliTPCtrack.cxx:165
 AliTPCtrack.cxx:166
 AliTPCtrack.cxx:167
 AliTPCtrack.cxx:168
 AliTPCtrack.cxx:169
 AliTPCtrack.cxx:170
 AliTPCtrack.cxx:171
 AliTPCtrack.cxx:172
 AliTPCtrack.cxx:173
 AliTPCtrack.cxx:174
 AliTPCtrack.cxx:175
 AliTPCtrack.cxx:176
 AliTPCtrack.cxx:177
 AliTPCtrack.cxx:178
 AliTPCtrack.cxx:179
 AliTPCtrack.cxx:180
 AliTPCtrack.cxx:181
 AliTPCtrack.cxx:182
 AliTPCtrack.cxx:183
 AliTPCtrack.cxx:184
 AliTPCtrack.cxx:185
 AliTPCtrack.cxx:186
 AliTPCtrack.cxx:187
 AliTPCtrack.cxx:188
 AliTPCtrack.cxx:189
 AliTPCtrack.cxx:190
 AliTPCtrack.cxx:191
 AliTPCtrack.cxx:192
 AliTPCtrack.cxx:193
 AliTPCtrack.cxx:194
 AliTPCtrack.cxx:195
 AliTPCtrack.cxx:196
 AliTPCtrack.cxx:197
 AliTPCtrack.cxx:198
 AliTPCtrack.cxx:199
 AliTPCtrack.cxx:200
 AliTPCtrack.cxx:201
 AliTPCtrack.cxx:202
 AliTPCtrack.cxx:203
 AliTPCtrack.cxx:204
 AliTPCtrack.cxx:205
 AliTPCtrack.cxx:206
 AliTPCtrack.cxx:207
 AliTPCtrack.cxx:208
 AliTPCtrack.cxx:209
 AliTPCtrack.cxx:210
 AliTPCtrack.cxx:211
 AliTPCtrack.cxx:212
 AliTPCtrack.cxx:213
 AliTPCtrack.cxx:214
 AliTPCtrack.cxx:215
 AliTPCtrack.cxx:216
 AliTPCtrack.cxx:217
 AliTPCtrack.cxx:218
 AliTPCtrack.cxx:219
 AliTPCtrack.cxx:220
 AliTPCtrack.cxx:221
 AliTPCtrack.cxx:222
 AliTPCtrack.cxx:223
 AliTPCtrack.cxx:224
 AliTPCtrack.cxx:225
 AliTPCtrack.cxx:226
 AliTPCtrack.cxx:227
 AliTPCtrack.cxx:228
 AliTPCtrack.cxx:229
 AliTPCtrack.cxx:230
 AliTPCtrack.cxx:231
 AliTPCtrack.cxx:232
 AliTPCtrack.cxx:233
 AliTPCtrack.cxx:234
 AliTPCtrack.cxx:235
 AliTPCtrack.cxx:236
 AliTPCtrack.cxx:237
 AliTPCtrack.cxx:238
 AliTPCtrack.cxx:239
 AliTPCtrack.cxx:240
 AliTPCtrack.cxx:241
 AliTPCtrack.cxx:242
 AliTPCtrack.cxx:243
 AliTPCtrack.cxx:244
 AliTPCtrack.cxx:245
 AliTPCtrack.cxx:246
 AliTPCtrack.cxx:247
 AliTPCtrack.cxx:248
 AliTPCtrack.cxx:249
 AliTPCtrack.cxx:250
 AliTPCtrack.cxx:251
 AliTPCtrack.cxx:252
 AliTPCtrack.cxx:253
 AliTPCtrack.cxx:254
 AliTPCtrack.cxx:255
 AliTPCtrack.cxx:256
 AliTPCtrack.cxx:257
 AliTPCtrack.cxx:258
 AliTPCtrack.cxx:259
 AliTPCtrack.cxx:260
 AliTPCtrack.cxx:261
 AliTPCtrack.cxx:262
 AliTPCtrack.cxx:263
 AliTPCtrack.cxx:264
 AliTPCtrack.cxx:265
 AliTPCtrack.cxx:266
 AliTPCtrack.cxx:267
 AliTPCtrack.cxx:268
 AliTPCtrack.cxx:269
 AliTPCtrack.cxx:270
 AliTPCtrack.cxx:271
 AliTPCtrack.cxx:272
 AliTPCtrack.cxx:273
 AliTPCtrack.cxx:274
 AliTPCtrack.cxx:275
 AliTPCtrack.cxx:276
 AliTPCtrack.cxx:277
 AliTPCtrack.cxx:278
 AliTPCtrack.cxx:279
 AliTPCtrack.cxx:280
 AliTPCtrack.cxx:281
 AliTPCtrack.cxx:282
 AliTPCtrack.cxx:283
 AliTPCtrack.cxx:284
 AliTPCtrack.cxx:285
 AliTPCtrack.cxx:286
 AliTPCtrack.cxx:287
 AliTPCtrack.cxx:288
 AliTPCtrack.cxx:289
 AliTPCtrack.cxx:290
 AliTPCtrack.cxx:291
 AliTPCtrack.cxx:292
 AliTPCtrack.cxx:293
 AliTPCtrack.cxx:294
 AliTPCtrack.cxx:295
 AliTPCtrack.cxx:296
 AliTPCtrack.cxx:297
 AliTPCtrack.cxx:298
 AliTPCtrack.cxx:299
 AliTPCtrack.cxx:300
 AliTPCtrack.cxx:301
 AliTPCtrack.cxx:302
 AliTPCtrack.cxx:303
 AliTPCtrack.cxx:304
 AliTPCtrack.cxx:305
 AliTPCtrack.cxx:306
 AliTPCtrack.cxx:307
 AliTPCtrack.cxx:308
 AliTPCtrack.cxx:309
 AliTPCtrack.cxx:310
 AliTPCtrack.cxx:311
 AliTPCtrack.cxx:312
 AliTPCtrack.cxx:313
 AliTPCtrack.cxx:314
 AliTPCtrack.cxx:315
 AliTPCtrack.cxx:316
 AliTPCtrack.cxx:317
 AliTPCtrack.cxx:318
 AliTPCtrack.cxx:319
 AliTPCtrack.cxx:320
 AliTPCtrack.cxx:321
 AliTPCtrack.cxx:322
 AliTPCtrack.cxx:323
 AliTPCtrack.cxx:324
 AliTPCtrack.cxx:325
 AliTPCtrack.cxx:326
 AliTPCtrack.cxx:327
 AliTPCtrack.cxx:328
 AliTPCtrack.cxx:329
 AliTPCtrack.cxx:330
 AliTPCtrack.cxx:331
 AliTPCtrack.cxx:332
 AliTPCtrack.cxx:333
 AliTPCtrack.cxx:334
 AliTPCtrack.cxx:335
 AliTPCtrack.cxx:336
 AliTPCtrack.cxx:337
 AliTPCtrack.cxx:338
 AliTPCtrack.cxx:339
 AliTPCtrack.cxx:340
 AliTPCtrack.cxx:341
 AliTPCtrack.cxx:342
 AliTPCtrack.cxx:343
 AliTPCtrack.cxx:344
 AliTPCtrack.cxx:345
 AliTPCtrack.cxx:346
 AliTPCtrack.cxx:347
 AliTPCtrack.cxx:348
 AliTPCtrack.cxx:349
 AliTPCtrack.cxx:350
 AliTPCtrack.cxx:351
 AliTPCtrack.cxx:352
 AliTPCtrack.cxx:353
 AliTPCtrack.cxx:354
 AliTPCtrack.cxx:355
 AliTPCtrack.cxx:356
 AliTPCtrack.cxx:357
 AliTPCtrack.cxx:358
 AliTPCtrack.cxx:359
 AliTPCtrack.cxx:360
 AliTPCtrack.cxx:361
 AliTPCtrack.cxx:362
 AliTPCtrack.cxx:363
 AliTPCtrack.cxx:364
 AliTPCtrack.cxx:365
 AliTPCtrack.cxx:366
 AliTPCtrack.cxx:367
 AliTPCtrack.cxx:368
 AliTPCtrack.cxx:369
 AliTPCtrack.cxx:370
 AliTPCtrack.cxx:371
 AliTPCtrack.cxx:372
 AliTPCtrack.cxx:373
 AliTPCtrack.cxx:374
 AliTPCtrack.cxx:375
 AliTPCtrack.cxx:376
 AliTPCtrack.cxx:377
 AliTPCtrack.cxx:378
 AliTPCtrack.cxx:379
 AliTPCtrack.cxx:380
 AliTPCtrack.cxx:381
 AliTPCtrack.cxx:382
 AliTPCtrack.cxx:383
 AliTPCtrack.cxx:384
 AliTPCtrack.cxx:385
 AliTPCtrack.cxx:386
 AliTPCtrack.cxx:387
 AliTPCtrack.cxx:388
 AliTPCtrack.cxx:389
 AliTPCtrack.cxx:390
 AliTPCtrack.cxx:391
 AliTPCtrack.cxx:392
 AliTPCtrack.cxx:393
 AliTPCtrack.cxx:394
 AliTPCtrack.cxx:395
 AliTPCtrack.cxx:396
 AliTPCtrack.cxx:397
 AliTPCtrack.cxx:398
 AliTPCtrack.cxx:399
 AliTPCtrack.cxx:400
 AliTPCtrack.cxx:401
 AliTPCtrack.cxx:402
 AliTPCtrack.cxx:403
 AliTPCtrack.cxx:404
 AliTPCtrack.cxx:405
 AliTPCtrack.cxx:406
 AliTPCtrack.cxx:407
 AliTPCtrack.cxx:408
 AliTPCtrack.cxx:409
 AliTPCtrack.cxx:410
 AliTPCtrack.cxx:411
 AliTPCtrack.cxx:412
 AliTPCtrack.cxx:413
 AliTPCtrack.cxx:414
 AliTPCtrack.cxx:415
 AliTPCtrack.cxx:416
 AliTPCtrack.cxx:417
 AliTPCtrack.cxx:418
 AliTPCtrack.cxx:419
 AliTPCtrack.cxx:420
 AliTPCtrack.cxx:421
 AliTPCtrack.cxx:422
 AliTPCtrack.cxx:423
 AliTPCtrack.cxx:424
 AliTPCtrack.cxx:425
 AliTPCtrack.cxx:426
 AliTPCtrack.cxx:427
 AliTPCtrack.cxx:428
 AliTPCtrack.cxx:429
 AliTPCtrack.cxx:430
 AliTPCtrack.cxx:431
 AliTPCtrack.cxx:432
 AliTPCtrack.cxx:433
 AliTPCtrack.cxx:434
 AliTPCtrack.cxx:435
 AliTPCtrack.cxx:436
 AliTPCtrack.cxx:437
 AliTPCtrack.cxx:438
 AliTPCtrack.cxx:439
 AliTPCtrack.cxx:440
 AliTPCtrack.cxx:441
 AliTPCtrack.cxx:442
 AliTPCtrack.cxx:443
 AliTPCtrack.cxx:444
 AliTPCtrack.cxx:445
 AliTPCtrack.cxx:446
 AliTPCtrack.cxx:447
 AliTPCtrack.cxx:448
 AliTPCtrack.cxx:449
 AliTPCtrack.cxx:450
 AliTPCtrack.cxx:451
 AliTPCtrack.cxx:452
 AliTPCtrack.cxx:453
 AliTPCtrack.cxx:454
 AliTPCtrack.cxx:455
 AliTPCtrack.cxx:456
 AliTPCtrack.cxx:457
 AliTPCtrack.cxx:458
 AliTPCtrack.cxx:459
 AliTPCtrack.cxx:460
 AliTPCtrack.cxx:461
 AliTPCtrack.cxx:462
 AliTPCtrack.cxx:463
 AliTPCtrack.cxx:464
 AliTPCtrack.cxx:465
 AliTPCtrack.cxx:466
 AliTPCtrack.cxx:467
 AliTPCtrack.cxx:468
 AliTPCtrack.cxx:469
 AliTPCtrack.cxx:470
 AliTPCtrack.cxx:471
 AliTPCtrack.cxx:472
 AliTPCtrack.cxx:473
 AliTPCtrack.cxx:474
 AliTPCtrack.cxx:475
 AliTPCtrack.cxx:476
 AliTPCtrack.cxx:477
 AliTPCtrack.cxx:478
 AliTPCtrack.cxx:479
 AliTPCtrack.cxx:480
 AliTPCtrack.cxx:481
 AliTPCtrack.cxx:482
 AliTPCtrack.cxx:483
 AliTPCtrack.cxx:484
 AliTPCtrack.cxx:485
 AliTPCtrack.cxx:486
 AliTPCtrack.cxx:487
 AliTPCtrack.cxx:488
 AliTPCtrack.cxx:489
 AliTPCtrack.cxx:490
 AliTPCtrack.cxx:491
 AliTPCtrack.cxx:492