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 ITS track class
//
//          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
//     dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
///////////////////////////////////////////////////////////////////////////
#include <TMath.h>

#include "AliCluster.h"
#include "AliESDVertex.h"
#include "AliITSReconstructor.h"
#include "AliITStrackV2.h"
#include "AliTracker.h"
#include "AliLog.h"
#include "AliPID.h"

const Int_t AliITStrackV2::fgkWARN = 5;

ClassImp(AliITStrackV2)


//____________________________________________________________________________
AliITStrackV2::AliITStrackV2() : AliKalmanTrack(),
  fCheckInvariant(kTRUE),
  fdEdx(0),
  fESDtrack(0)
{
  for(Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) {fIndex[i]=-1; fModule[i]=-1;}
  for(Int_t i=0; i<AliITSgeomTGeo::kNLayers; i++) {fSharedWeight[i]=0;}
  for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
}


//____________________________________________________________________________
AliITStrackV2::AliITStrackV2(AliESDtrack& t,Bool_t c):
  AliKalmanTrack(),
  fCheckInvariant(kTRUE),
  fdEdx(t.GetITSsignal()),
  fESDtrack(&t)
{
  //------------------------------------------------------------------
  // Conversion ESD track -> ITS track.
  // If c==kTRUE, create the ITS track out of the constrained params.
  //------------------------------------------------------------------
  const AliExternalTrackParam *par=&t;
  if (c) {
    par=t.GetConstrainedParam();
    if (!par) AliError("AliITStrackV2: conversion failed !\n");
  }
  Set(par->GetX(),par->GetAlpha(),par->GetParameter(),par->GetCovariance());

  SetLabel(t.GetITSLabel());
  SetMass(t.GetMassForTracking());
  SetNumberOfClusters(t.GetITSclusters(fIndex));

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

  for(Int_t i=0; i<AliITSgeomTGeo::kNLayers; i++) {fSharedWeight[i]=0;}
  for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
}

//____________________________________________________________________________
void AliITStrackV2::ResetClusters() {
  //------------------------------------------------------------------
  // Reset the array of attached clusters.
  //------------------------------------------------------------------
  for (Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) fIndex[i]=-1;
  for (Int_t i=0; i<AliITSgeomTGeo::kNLayers; i++) {fSharedWeight[i]=0;}
  SetChi2(0.); 
  SetNumberOfClusters(0);
} 

//____________________________________________________________________________
void AliITStrackV2::UpdateESDtrack(ULong_t flags) const {
  // Update track params
  fESDtrack->UpdateTrackParams(this,flags);
  //
  // set correctly the global label
  if (fESDtrack->IsOn(AliESDtrack::kTPCin)) { 
    // for global track the GetLabel should be negative if
    // 1) GetTPCLabel<0
    // 2) this->GetLabel()<0
    // 3) GetTPCLabel() != this->GetLabel()
    int label = fESDtrack->GetTPCLabel();
    int itsLabel = GetLabel();
    if (label<0 || itsLabel<0 || itsLabel!=label) label = -TMath::Abs(label);
    fESDtrack->SetLabel(label);
  }
  //
  // copy the module indices
  Int_t i;
  for(i=0;i<2*AliITSgeomTGeo::kNLayers;i++) {
    //   printf("     %d\n",GetModuleIndex(i));
    fESDtrack->SetITSModuleIndex(i,GetModuleIndex(i));
  }
  // copy the map of shared clusters
  if(flags==AliESDtrack::kITSin) {
    UChar_t itsSharedMap=0;
    for(i=0;i<AliITSgeomTGeo::kNLayers;i++) {
      if(fSharedWeight[i]>0) SETBIT(itsSharedMap,i);
      
    }
    fESDtrack->SetITSSharedMap(itsSharedMap);
  }

  // copy the 4 dedx samples
  Double_t sdedx[4]={0.,0.,0.,0.};
  for(i=0; i<4; i++) sdedx[i]=fdEdxSample[i];
  fESDtrack->SetITSdEdxSamples(sdedx);
}

//____________________________________________________________________________
AliITStrackV2::AliITStrackV2(const AliITStrackV2& t) : 
  AliKalmanTrack(t),
  fCheckInvariant(t.fCheckInvariant),
  fdEdx(t.fdEdx),
  fESDtrack(t.fESDtrack) 
{
  //------------------------------------------------------------------
  //Copy constructor
  //------------------------------------------------------------------
  Int_t i;
  for (i=0; i<4; i++) fdEdxSample[i]=t.fdEdxSample[i];
  for (i=0; i<2*AliITSgeomTGeo::GetNLayers(); i++) {
    fIndex[i]=t.fIndex[i];
    fModule[i]=t.fModule[i];
  }
  for (i=0; i<AliITSgeomTGeo::kNLayers; i++) {fSharedWeight[i]=t.fSharedWeight[i];}
}

//_____________________________________________________________________________
Int_t AliITStrackV2::Compare(const TObject *o) const {
  //-----------------------------------------------------------------
  // This function compares tracks according to the their curvature
  //-----------------------------------------------------------------
  AliITStrackV2 *t=(AliITStrackV2*)o;
  //Double_t co=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;
}

//____________________________________________________________________________
Bool_t 
AliITStrackV2::PropagateToVertex(const AliESDVertex *v,Double_t d,Double_t x0) 
{
  //------------------------------------------------------------------
  //This function propagates a track to the minimal distance from the origin
  //------------------------------------------------------------------  
  Double_t bz=GetBz();
  if (PropagateToDCA(v,bz,kVeryBig)) {
    Double_t xOverX0,xTimesRho; 
    xOverX0 = d; xTimesRho = d*x0;
    if (CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kTRUE;
  }
  return kFALSE;
}

//____________________________________________________________________________
Bool_t AliITStrackV2::
GetGlobalXYZat(Double_t xloc, Double_t &x, Double_t &y, Double_t &z) const {
  //------------------------------------------------------------------
  //This function returns a track position in the global system
  //------------------------------------------------------------------
  Double_t r[3];
  Bool_t rc=GetXYZAt(xloc, GetBz(), r);
  x=r[0]; y=r[1]; z=r[2]; 
  return rc;
}

//_____________________________________________________________________________
Double_t AliITStrackV2::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 AliITStrackV2::PropagateTo(Double_t xk, Double_t d, Double_t x0) {
  //------------------------------------------------------------------
  //This function propagates a track
  //------------------------------------------------------------------

  Double_t oldX=GetX(), oldY=GetY(), oldZ=GetZ();
  
  //Double_t bz=GetBz();
  //if (!AliExternalTrackParam::PropagateTo(xk,bz)) return kFALSE;
  Double_t b[3]; GetBxByBz(b);
  if (!AliExternalTrackParam::PropagateToBxByBz(xk,b)) return kFALSE;
  Double_t xOverX0,xTimesRho; 
  xOverX0 = d; xTimesRho = d*x0;
  if (!CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kFALSE;

  Double_t x=GetX(), y=GetY(), z=GetZ();
  if (IsStartedTimeIntegral() && x>oldX) {
    Double_t l2 = (x-oldX)*(x-oldX) + (y-oldY)*(y-oldY) + (z-oldZ)*(z-oldZ);
    AddTimeStep(TMath::Sqrt(l2));
  }

  return kTRUE;
}

//____________________________________________________________________________
Bool_t AliITStrackV2::PropagateToTGeo(Double_t xToGo, Int_t nstep, Double_t &xOverX0, Double_t &xTimesRho, Bool_t addTime) {
  //-------------------------------------------------------------------
  //  Propagates the track to a reference plane x=xToGo in n steps.
  //  These n steps are only used to take into account the curvature.
  //  The material is calculated with TGeo. (L.Gaudichet)
  //-------------------------------------------------------------------
  
  Double_t startx = GetX(), starty = GetY(), startz = GetZ();
  Double_t sign = (startx<xToGo) ? -1.:1.;
  Double_t step = (xToGo-startx)/TMath::Abs(nstep);

  Double_t start[3], end[3], mparam[7];
  //Double_t bz = GetBz();
  Double_t b[3]; GetBxByBz(b);
  Double_t bz = b[2];

  Double_t x = startx;
  
  for (Int_t i=0; i<nstep; i++) {
    
    GetXYZ(start);   //starting global position
    x += step;
    if (!GetXYZAt(x, bz, end)) return kFALSE;
    //if (!AliExternalTrackParam::PropagateTo(x, bz)) return kFALSE;
    if (!AliExternalTrackParam::PropagateToBxByBz(x, b)) return kFALSE;
    AliTracker::MeanMaterialBudget(start, end, mparam);
    xTimesRho = sign*mparam[4]*mparam[0];
    xOverX0   = mparam[1];
    if (mparam[1]<900000) {
      if (!AliExternalTrackParam::CorrectForMeanMaterial(xOverX0,
			   xTimesRho,GetMass())) return kFALSE;
    } else { // this happens when MeanMaterialBudget cannot cross a boundary
      return kFALSE;
    }
  }

  if (addTime && IsStartedTimeIntegral() && GetX()>startx) {
    Double_t l2 = ( (GetX()-startx)*(GetX()-startx) +
		    (GetY()-starty)*(GetY()-starty) +
		    (GetZ()-startz)*(GetZ()-startz) );
    AddTimeStep(TMath::Sqrt(l2));
  }

  return kTRUE;
}

//____________________________________________________________________________
Bool_t AliITStrackV2::Update(const AliCluster* c, Double_t chi2, Int_t index) 
{
  //------------------------------------------------------------------
  //This function updates track parameters
  //------------------------------------------------------------------
  Double_t p[2]={c->GetY(), c->GetZ()};
  Double_t cov[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};

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

  Int_t n=GetNumberOfClusters();
  if (!Invariant()) {
    if (n>fgkWARN) AliDebug(1,"Wrong invariant !");
     return kFALSE;
  }

  if (chi2<0) return kTRUE;

  // fill residuals for ITS+TPC tracks 
  if (fESDtrack) {
    if (fESDtrack->GetStatus()&AliESDtrack::kTPCin) {
      AliTracker::FillResiduals(this,p,cov,c->GetVolumeId());
    }
  }

  fIndex[n]=index;
  SetNumberOfClusters(n+1);
  SetChi2(GetChi2()+chi2);

  return kTRUE;
}

Bool_t AliITStrackV2::Invariant() const {
  //------------------------------------------------------------------
  // This function is for debugging purpose only
  //------------------------------------------------------------------
  if(!fCheckInvariant) return kTRUE;

  Int_t n=GetNumberOfClusters();
  static Float_t bz = GetBz();
  // take into account the misalignment error
  Float_t maxMisalErrY2=0,maxMisalErrZ2=0;
  //RS
  const AliITSRecoParam* recopar = AliITSReconstructor::GetRecoParam();
  if (!recopar) recopar = AliITSRecoParam::GetHighFluxParam();

  for (Int_t lay=0; lay<AliITSgeomTGeo::kNLayers; lay++) {
    maxMisalErrY2 = TMath::Max(maxMisalErrY2,recopar->GetClusterMisalErrorY(lay,bz));
    maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,recopar->GetClusterMisalErrorZ(lay,bz));
  }
  maxMisalErrY2 *= maxMisalErrY2;
  maxMisalErrZ2 *= maxMisalErrZ2;
  // this is because when we reset before refitting, we multiply the
  // matrix by 10
  maxMisalErrY2 *= 10.; 
  maxMisalErrZ2 *= 10.;

  Double_t sP2=GetParameter()[2];
  if (TMath::Abs(sP2) >= kAlmost1){
    if (n>fgkWARN) AliDebug(1,Form("fP2=%f\n",sP2));
     return kFALSE;
  }
  Double_t sC00=GetCovariance()[0];
  if (sC00<=0 || sC00>(9.+maxMisalErrY2)) {
    if (n>fgkWARN) AliDebug(1,Form("fC00=%f\n",sC00)); 
     return kFALSE;
  }
  Double_t sC11=GetCovariance()[2];
  if (sC11<=0 || sC11>(9.+maxMisalErrZ2)) {
    if (n>fgkWARN) AliDebug(1,Form("fC11=%f\n",sC11)); 
     return kFALSE;
  }
  Double_t sC22=GetCovariance()[5];
  if (sC22<=0 || sC22>1.) {
    if (n>fgkWARN) AliDebug(1,Form("fC22=%f\n",sC22)); 
     return kFALSE;
  }
  Double_t sC33=GetCovariance()[9];
  if (sC33<=0 || sC33>1.) {
    if (n>fgkWARN) AliDebug(1,Form("fC33=%f\n",sC33)); 
     return kFALSE;
  }
  Double_t sC44=GetCovariance()[14];
  if (sC44<=0 /*|| sC44>6e-5*/) {
    if (n>fgkWARN) AliDebug(1,Form("fC44=%f\n",sC44));
     return kFALSE;
  }

  return kTRUE;
}

//____________________________________________________________________________
Bool_t AliITStrackV2::Propagate(Double_t alp,Double_t xk) {
  //------------------------------------------------------------------
  //This function propagates a track
  //------------------------------------------------------------------
  //Double_t bz=GetBz();
  //if (!AliExternalTrackParam::Propagate(alp,xk,bz)) return kFALSE;
  Double_t b[3]; GetBxByBz(b);
  if (!AliExternalTrackParam::PropagateBxByBz(alp,xk,b)) return kFALSE;

  if (!Invariant()) {
    Int_t n=GetNumberOfClusters();
    if (n>fgkWARN) AliDebug(1,"Wrong invariant !");
    return kFALSE;
  }

  return kTRUE;
}

Bool_t AliITStrackV2::MeanBudgetToPrimVertex(Double_t xyz[3], Double_t step, Double_t &d) const {

  //-------------------------------------------------------------------
  //  Get the mean material budget between the actual point and the
  //  primary vertex. (L.Gaudichet)
  //-------------------------------------------------------------------

  Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
  Double_t vertexX = xyz[0]*cs + xyz[1]*sn;

  Int_t nstep = Int_t((GetX()-vertexX)/step);
  if (nstep<1) nstep = 1;
  step = (GetX()-vertexX)/nstep;

  //  Double_t mparam[7], densMean=0, radLength=0, length=0;
  Double_t mparam[7];
  Double_t p1[3], p2[3], x = GetX(), bz = GetBz();
  GetXYZ(p1);

  d=0.;

  for (Int_t i=0; i<nstep; i++) {
    x  += step;
    if (!GetXYZAt(x, bz, p2)) return kFALSE;
    AliTracker::MeanMaterialBudget(p1, p2, mparam);
    if (mparam[1]>900000) return kFALSE;
    d  += mparam[1];

    p1[0] = p2[0];
    p1[1] = p2[1];
    p1[2] = p2[2];
  }

  return kTRUE;
}

Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) {
  //------------------------------------------------------------------
  //This function improves angular track parameters
  //------------------------------------------------------------------
  //Store the initail track parameters

  Double_t x = GetX();
  Double_t alpha = GetAlpha();
  Double_t par[] = {GetY(),GetZ(),GetSnp(),GetTgl(),GetSigned1Pt()};
  Double_t cov[] = {
    GetSigmaY2(),
    GetSigmaZY(),
    GetSigmaZ2(),
    GetSigmaSnpY(),
    GetSigmaSnpZ(),
    GetSigmaSnp2(),
    GetSigmaTglY(),
    GetSigmaTglZ(),
    GetSigmaTglSnp(),
    GetSigmaTgl2(),
    GetSigma1PtY(),
    GetSigma1PtZ(),
    GetSigma1PtSnp(),
    GetSigma1PtTgl(),
    GetSigma1Pt2()
  }; 


  Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
  Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex
  Double_t yv =-xyz[0]*sn + xyz[1]*cs; // in the
  Double_t zv = xyz[2];                // local frame

  Double_t dx = x - xv, dy = par[0] - yv, dz = par[1] - zv;
  Double_t r2=dx*dx + dy*dy;
  Double_t p2=(1.+ GetTgl()*GetTgl())/(GetSigned1Pt()*GetSigned1Pt());
  if (GetMass()<0) p2 *= 4; // q=2
  Double_t beta2=p2/(p2 + GetMass()*GetMass());
  x0*=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(1.- GetSnp()*GetSnp()));
  Double_t theta2=14.1*14.1/(beta2*p2*1e6)*x0;
  //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*x0*9.36*2.33;

  Double_t bz=GetBz();
  Double_t cnv=bz*kB2C;
  Double_t curv=GetC(bz);
  {
    Double_t dummy = 4/r2 - curv*curv;
    if (dummy < 0) return kFALSE;
    Double_t parp = 0.5*(curv*dx + dy*TMath::Sqrt(dummy));
    Double_t sigma2p = theta2*(1.-GetSnp())*(1.+GetSnp())*(1. + GetTgl()*GetTgl());
    Double_t ovSqr2 = 1./TMath::Sqrt(r2);
    Double_t tfact = ovSqr2*(1.-dy*ovSqr2)*(1.+dy*ovSqr2);
    sigma2p += cov[0]*tfact*tfact;
    sigma2p += ers[1]*ers[1]/r2;
    sigma2p += 0.25*cov[14]*cnv*cnv*dx*dx;
    Double_t eps2p=sigma2p/(cov[5] + sigma2p);
    par[0] += cov[3]/(cov[5] + sigma2p)*(parp - GetSnp());
    par[2]  = eps2p*GetSnp() + (1 - eps2p)*parp;
    cov[5] *= eps2p;
    cov[3] *= eps2p;
  }
  {
    Double_t parl=0.5*curv*dz/TMath::ASin(0.5*curv*TMath::Sqrt(r2));
    Double_t sigma2l=theta2;
    sigma2l += cov[2]/r2 + cov[0]*dy*dy*dz*dz/(r2*r2*r2);
    sigma2l += ers[2]*ers[2]/r2;
    Double_t eps2l = sigma2l/(cov[9] + sigma2l);
    par[1] += cov[7 ]/(cov[9] + sigma2l)*(parl - par[3]);
    par[4] += cov[13]/(cov[9] + sigma2l)*(parl - par[3]);
    par[3]  = eps2l*par[3] + (1-eps2l)*parl;
    cov[9] *= eps2l; 
    cov[13]*= eps2l; 
    cov[7] *= eps2l; 
  }

  Set(x,alpha,par,cov);

  if (!Invariant()) return kFALSE;

  return kTRUE;
}

void AliITStrackV2::CookdEdx(Double_t /*low*/, Double_t /*up*/) {
  //-----------------------------------------------------------------
  // This function calculates dE/dX within the "low" and "up" cuts.
  // Origin: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
  // Updated: F. Prino 8-June-2009
  //-----------------------------------------------------------------
  // The cluster order is: SDD-1, SDD-2, SSD-1, SSD-2

  Int_t nc=0;
  Float_t dedx[4];
  for (Int_t il=0; il<4; il++) { // count good (>0) dE/dx values
    if(fdEdxSample[il]>0.){
      dedx[nc]= fdEdxSample[il];
      nc++;
    }
  }
  if(nc<1){
    SetdEdx(0.);
    return;
  }

  Int_t swap; // sort in ascending order
  do {
    swap=0;
    for (Int_t i=0; i<nc-1; i++) {
      if (dedx[i]<=dedx[i+1]) continue;
      Float_t tmp=dedx[i];
      dedx[i]=dedx[i+1]; 
      dedx[i+1]=tmp;
      swap++;
    }
  } while (swap);


  Double_t sumamp=0,sumweight=0;
  Double_t weight[4]={1.,1.,0.,0.};
  if(nc==3) weight[1]=0.5;
  else if(nc<3) weight[1]=0.;
  for (Int_t i=0; i<nc; i++) {
    sumamp+= dedx[i]*weight[i];
    sumweight+=weight[i];
  }
  SetdEdx(sumamp/sumweight);
}

//____________________________________________________________________________
Bool_t AliITStrackV2::
GetPhiZat(Double_t r, Double_t &phi, Double_t &z) const {
  //------------------------------------------------------------------
  // This function returns the global cylindrical (phi,z) of the track 
  // position estimated at the radius r. 
  // The track curvature is neglected.
  //------------------------------------------------------------------
  Double_t d=GetD(0.,0.);
  if (TMath::Abs(d) > r) {
    if (r>1e-1) return kFALSE;
    r = TMath::Abs(d);
  }

  Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
  if (TMath::Abs(d) > rcurr) return kFALSE;
  Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
  Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);

  if (GetX()>=0.) {
    phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
  } else {
    phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
  }

  // return a phi in [0,2pi[ 
  if (phi<0.) phi+=2.*TMath::Pi();
  else if (phi>=2.*TMath::Pi()) phi-=2.*TMath::Pi();
  z=GetZ()+GetTgl()*(TMath::Sqrt((r-d)*(r+d))-TMath::Sqrt((rcurr-d)*(rcurr+d)));
  return kTRUE;
}
//____________________________________________________________________________
Bool_t AliITStrackV2::
GetLocalXat(Double_t r,Double_t &xloc) const {
  //------------------------------------------------------------------
  // This function returns the local x of the track 
  // position estimated at the radius r. 
  // The track curvature is neglected.
  //------------------------------------------------------------------
  Double_t d=GetD(0.,0.);
  if (TMath::Abs(d) > r) { 
    if (r>1e-1) return kFALSE; 
    r = TMath::Abs(d); 
  } 

  Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
  Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
  Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
  Double_t phi;
  if (GetX()>=0.) {
    phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
  } else {
    phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
  }

  xloc=r*(TMath::Cos(phi)*TMath::Cos(GetAlpha())
         +TMath::Sin(phi)*TMath::Sin(GetAlpha())); 

  return kTRUE;
}

//____________________________________________________________________________
Bool_t AliITStrackV2::ImproveKalman(Double_t xyz[3],Double_t ers[3], const Double_t* xlMS, const Double_t* x2X0MS, Int_t nMS)
{
  // Substitute the state of the track (p_{k|k},C_{k|k}) at the k-th measumerent by its
  // smoothed value from the k-th measurement + measurement at the vertex.
  // Account for the MS on nMS layers at x-postions xlMS with x/x0 = x2X0MS
  // p_{k|kv} = p_{k|k} + C_{k|k}*D^Tr_{k+1} B^{-1}_{k+1} ( vtx - D_{k+1}*p_{k|k})
  // C_{k|kv} = C_{k|k}*( I - D^Tr_{k+1} B^{-1}_{k+1} D_{k+1} C^Tr_{k|k})
  // 
  // where D_{k} = H_{k} F_{k} with H being the matrix converting the tracks parameters
  // to measurements m_{k} = H_{k} p_{k} and F_{k} the matrix propagating the track between the
  // the point k-1 and k:  p_{k|k-1} = F_{k} p_{k-1|k-1}
  //
  // B_{k+1} = V_{k+1} + H_{k+1} C_{k+1|k} H^Tr_{k+1} with V_{k+1} being the error of the measurment
  // at point k+1 (i.e. vertex), and C_{k+1|k} - error matrix extrapolated from k-th measurement to
  // k+1 (vtx) and accounting for the MS inbetween
  //
  // H = {{1,0,0,0,0},{0,1,0,0,0}}
  //
  double covc[15], *cori = (double*) GetCovariance(),par[5] = {GetY(),GetZ(),GetSnp(),GetTgl(),GetSigned1Pt()},
    &c00=cori[0],
    &c01=cori[1],&c11=cori[2],
    &c02=cori[3],&c12=cori[4],&c22=cori[5],
    &c03=cori[6],&c13=cori[7],&c23=cori[8],&c33=cori[9],
    &c04=cori[10],&c14=cori[11],&c24=cori[12],&c34=cori[13],&c44=cori[14],
    // for smoothed cov matrix 
    &cov00=covc[0],
    &cov01=covc[1],&cov11=covc[2],
    &cov02=covc[3],&cov12=covc[4],&cov22=covc[5],
    &cov03=covc[6],&cov13=covc[7],&cov23=covc[8],&cov33=covc[9],
    &cov04=covc[10],&cov14=covc[11],&cov24=covc[12],&cov34=covc[13],&cov44=covc[14];
  //
  double x = GetX(), alpha = GetAlpha();
  // vertex in the track frame
  double cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
  double xv = xyz[0]*cs + xyz[1]*sn, yv =-xyz[0]*sn + xyz[1]*cs, zv = xyz[2];
  double dx = xv - GetX();
  if (TMath::Abs(dx)<=kAlmost0)  return kTRUE;
  //
  double cnv=GetBz()*kB2C, x2r=cnv*par[4]*dx, f1=par[2], f2=f1+x2r;
  if (TMath::Abs(f1) >= kAlmost1 || TMath::Abs(f2) >= kAlmost1) {
    AliInfo(Form("Fail: %+e %+e",f1,f2));
    return kFALSE;
  }
  double r1=TMath::Sqrt((1.-f1)*(1.+f1)), r2=TMath::Sqrt((1.-f2)*(1.+f2)), dx2r=dx/(r1+r2);
  // elements of matrix F_{k+1} (1s on diagonal)
  double f02 = 2*dx2r, f04 = cnv*dx*dx2r, f13/*, f24 = cnv*dx*/;
  if (TMath::Abs(x2r)<0.05) f13 = dx*r2+f2*(f1+f2)*dx2r; // see AliExternalTrackParam::PropagateTo
  else {
    double dy2dx = (f1+f2)/(r1+r2);
    f13 = 2*TMath::ASin(0.5*TMath::Sqrt(1+dy2dx*dy2dx)*x2r)/(cnv*par[4]);
  }  
  // elements of matrix D_{k+1} = H_{k+1} * F_{k+1}
  // double d00 = 1., d11 = 1.;
  double &d02 = f02,  &d04 = f04, &d13 = f13;
  //
  // elements of matrix DC = D_{k+1}*C_{kk}^T
  double dc00 = c00+c02*d02+c04*d04,   dc10 = c01+c03*d13;
  double dc01 = c01+c12*d02+c14*d04,   dc11 = c11+c13*d13;
  double dc02 = c02+c22*d02+c24*d04,   dc12 = c12+c23*d13;
  double dc03 = c03+c23*d02+c34*d04,   dc13 = c13+c33*d13;
  double dc04 = c04+c24*d02+c44*d04,   dc14 = c14+c34*d13;
  //
  // difference between the vertex and the the track extrapolated to vertex
  yv -= par[0] + par[2]*d02 + par[4]*d04;
  zv -= par[1] + par[3]*d13;
  //
  // y,z part of the cov.matrix extrapolated to vtx (w/o MS contribution)
  // C_{k+1,k} = H F_{k+1} C_{k,k} F^Tr_{k+1} H^Tr = D C D^Tr
  double cv00 = dc00+dc02*d02+dc04*d04, cv01 = dc01+dc03*d13, cv11 = dc11+dc13*d13;
  //
  // add MS contribution layer by layer
  double xCurr = x;
  double p2Curr = par[2];
  //
  // precalculated factors of MS contribution matrix:
  double ms22t = (1. + par[3]*par[3]);
  double ms33t = ms22t*ms22t;
  double p34 = par[3]*par[4];
  double ms34t = p34*ms22t;
  double ms44t = p34*p34;
  //
  double p2=(1.+ par[3]*par[3])/(par[4]*par[4]);
  if (GetMass()<0) p2 *= 4; // q=2
  double beta2 = p2/(p2+GetMass()*GetMass());
  double theta2t = 14.1*14.1/(beta2*p2*1e6) * (1. + par[3]*par[3]);
  //
  // account for the MS in the layers between the last measurement and the vertex
  for (int il=0;il<nMS;il++) {
    double dfx = xlMS[il] - xCurr;
    xCurr = xlMS[il];
    p2Curr += dfx*cnv*par[4];   // p2 at the scattering layer
    double dxL=xv - xCurr;    // distance from scatering layer to vtx
    double x2rL=cnv*par[4]*dxL, f1L=p2Curr, f2L=f1L+x2rL;
    if (TMath::Abs(f1L) >= kAlmost1 || TMath::Abs(f2L) >= kAlmost1) {
      AliInfo(Form("FailMS at step %d of %d: dfx:%e dxL:%e %e %e",il,nMS,dfx,dxL,f1L,f2L));
      return kFALSE;
    }
    double r1L=TMath::Sqrt((1.-f1L)*(1.+f1L)), r2L=TMath::Sqrt((1.-f2L)*(1.+f2L)), dx2rL=dxL/(r1L+r2L);
    // elements of matrix for propagation from scatering layer to vertex
    double f02L = 2*dx2rL, f04L = cnv*dxL*dx2rL, f13L/*, f24L = cnv*dxL*/;
    if (TMath::Abs(x2rL)<0.05) f13L = dxL*r2L+f2L*(f1L+f2L)*dx2rL; // see AliExternalTrackParam::PropagateTo
    else {
      double dy2dxL = (f1L+f2L)/(r1L+r2L);
      f13L = 2*TMath::ASin(0.5*TMath::Sqrt(1+dy2dxL*dy2dxL)*x2rL)/(cnv*par[4]);
    }
    // MS contribution matrix:
    double theta2 = theta2t*TMath::Abs(x2X0MS[il]);
    double ms22 = theta2*(1.-p2Curr)*(1.+p2Curr)*ms22t;
    double ms33 = theta2*ms33t;
    double ms34 = theta2*ms34t;
    double ms44 = theta2*ms44t;
    //
    // add  H F MS F^Tr H^Tr to cv
    cv00 += f02L*f02L*ms22 + f04L*f04L*ms44;
    cv01 += f04L*f13L*ms34;
    cv11 += f13L*f13L*ms33;
  }
  //
  // inverse of matrix B
  double b11 = ers[1]*ers[1] + cv00;
  double b00 = ers[2]*ers[2] + cv11;
  double det = b11*b00 - cv01*cv01;
  if (TMath::Abs(det)<kAlmost0) {
    AliInfo(Form("Fail on det %e: %e %e %e",det,cv00,cv11,cv01));
    return kFALSE;
  }
  det = 1./det;
  b00 *= det; b11 *= det; 
  double b01 = -cv01*det;
  //
  // elements of matrix DC^Tr * B^-1
  double dcb00 = b00*dc00+b01*dc10, dcb01 = b01*dc00+b11*dc10;
  double dcb10 = b00*dc01+b01*dc11, dcb11 = b01*dc01+b11*dc11;
  double dcb20 = b00*dc02+b01*dc12, dcb21 = b01*dc02+b11*dc12;
  double dcb30 = b00*dc03+b01*dc13, dcb31 = b01*dc03+b11*dc13;
  double dcb40 = b00*dc04+b01*dc14, dcb41 = b01*dc04+b11*dc14;
  //
  // p_{k|k+1} = p_{k|k} +  C_{k|k}*D^Tr_{k+1} B^{-1}_{k+1} ( vtx - D_{k+1}*p_{k|k})
  par[0] += dcb00*yv + dcb01*zv;
  par[1] += dcb10*yv + dcb11*zv;
  par[2] += dcb20*yv + dcb21*zv;
  par[3] += dcb30*yv + dcb31*zv;
  par[4] += dcb40*yv + dcb41*zv;
  //
  // C_{k|kv} = C_{k|k} - C_{k|k} D^Tr_{k+1} B^{-1}_{k+1} D_{k+1} C^Tr_{k|k})
  //
  cov00 = c00 - (dc00*dcb00 + dc10*dcb01);
  cov01 = c01 - (dc01*dcb00 + dc11*dcb01);
  cov02 = c02 - (dc02*dcb00 + dc12*dcb01);
  cov03 = c03 - (dc03*dcb00 + dc13*dcb01);
  cov04 = c04 - (dc04*dcb00 + dc14*dcb01);
  //
  cov11 = c11 - (dc01*dcb10 + dc11*dcb11);
  cov12 = c12 - (dc02*dcb10 + dc12*dcb11);
  cov13 = c13 - (dc03*dcb10 + dc13*dcb11);
  cov14 = c14 - (dc04*dcb10 + dc14*dcb11);
  //
  cov22 = c22 - (dc02*dcb20 + dc12*dcb21);
  cov23 = c23 - (dc03*dcb20 + dc13*dcb21);
  cov24 = c24 - (dc04*dcb20 + dc14*dcb21);
  //
  cov33 = c33 - (dc03*dcb30 + dc13*dcb31);
  cov34 = c34 - (dc04*dcb30 + dc14*dcb31);
  //
  cov44 = c44 - (dc04*dcb40 + dc14*dcb41);
  //
  Set(x,alpha,par,covc);
  if (!Invariant()) {
    AliInfo(Form("Fail on Invariant, X=%e",GetX()));
    return kFALSE;
  }
  return kTRUE;
  //
}
 AliITStrackV2.cxx:1
 AliITStrackV2.cxx:2
 AliITStrackV2.cxx:3
 AliITStrackV2.cxx:4
 AliITStrackV2.cxx:5
 AliITStrackV2.cxx:6
 AliITStrackV2.cxx:7
 AliITStrackV2.cxx:8
 AliITStrackV2.cxx:9
 AliITStrackV2.cxx:10
 AliITStrackV2.cxx:11
 AliITStrackV2.cxx:12
 AliITStrackV2.cxx:13
 AliITStrackV2.cxx:14
 AliITStrackV2.cxx:15
 AliITStrackV2.cxx:16
 AliITStrackV2.cxx:17
 AliITStrackV2.cxx:18
 AliITStrackV2.cxx:19
 AliITStrackV2.cxx:20
 AliITStrackV2.cxx:21
 AliITStrackV2.cxx:22
 AliITStrackV2.cxx:23
 AliITStrackV2.cxx:24
 AliITStrackV2.cxx:25
 AliITStrackV2.cxx:26
 AliITStrackV2.cxx:27
 AliITStrackV2.cxx:28
 AliITStrackV2.cxx:29
 AliITStrackV2.cxx:30
 AliITStrackV2.cxx:31
 AliITStrackV2.cxx:32
 AliITStrackV2.cxx:33
 AliITStrackV2.cxx:34
 AliITStrackV2.cxx:35
 AliITStrackV2.cxx:36
 AliITStrackV2.cxx:37
 AliITStrackV2.cxx:38
 AliITStrackV2.cxx:39
 AliITStrackV2.cxx:40
 AliITStrackV2.cxx:41
 AliITStrackV2.cxx:42
 AliITStrackV2.cxx:43
 AliITStrackV2.cxx:44
 AliITStrackV2.cxx:45
 AliITStrackV2.cxx:46
 AliITStrackV2.cxx:47
 AliITStrackV2.cxx:48
 AliITStrackV2.cxx:49
 AliITStrackV2.cxx:50
 AliITStrackV2.cxx:51
 AliITStrackV2.cxx:52
 AliITStrackV2.cxx:53
 AliITStrackV2.cxx:54
 AliITStrackV2.cxx:55
 AliITStrackV2.cxx:56
 AliITStrackV2.cxx:57
 AliITStrackV2.cxx:58
 AliITStrackV2.cxx:59
 AliITStrackV2.cxx:60
 AliITStrackV2.cxx:61
 AliITStrackV2.cxx:62
 AliITStrackV2.cxx:63
 AliITStrackV2.cxx:64
 AliITStrackV2.cxx:65
 AliITStrackV2.cxx:66
 AliITStrackV2.cxx:67
 AliITStrackV2.cxx:68
 AliITStrackV2.cxx:69
 AliITStrackV2.cxx:70
 AliITStrackV2.cxx:71
 AliITStrackV2.cxx:72
 AliITStrackV2.cxx:73
 AliITStrackV2.cxx:74
 AliITStrackV2.cxx:75
 AliITStrackV2.cxx:76
 AliITStrackV2.cxx:77
 AliITStrackV2.cxx:78
 AliITStrackV2.cxx:79
 AliITStrackV2.cxx:80
 AliITStrackV2.cxx:81
 AliITStrackV2.cxx:82
 AliITStrackV2.cxx:83
 AliITStrackV2.cxx:84
 AliITStrackV2.cxx:85
 AliITStrackV2.cxx:86
 AliITStrackV2.cxx:87
 AliITStrackV2.cxx:88
 AliITStrackV2.cxx:89
 AliITStrackV2.cxx:90
 AliITStrackV2.cxx:91
 AliITStrackV2.cxx:92
 AliITStrackV2.cxx:93
 AliITStrackV2.cxx:94
 AliITStrackV2.cxx:95
 AliITStrackV2.cxx:96
 AliITStrackV2.cxx:97
 AliITStrackV2.cxx:98
 AliITStrackV2.cxx:99
 AliITStrackV2.cxx:100
 AliITStrackV2.cxx:101
 AliITStrackV2.cxx:102
 AliITStrackV2.cxx:103
 AliITStrackV2.cxx:104
 AliITStrackV2.cxx:105
 AliITStrackV2.cxx:106
 AliITStrackV2.cxx:107
 AliITStrackV2.cxx:108
 AliITStrackV2.cxx:109
 AliITStrackV2.cxx:110
 AliITStrackV2.cxx:111
 AliITStrackV2.cxx:112
 AliITStrackV2.cxx:113
 AliITStrackV2.cxx:114
 AliITStrackV2.cxx:115
 AliITStrackV2.cxx:116
 AliITStrackV2.cxx:117
 AliITStrackV2.cxx:118
 AliITStrackV2.cxx:119
 AliITStrackV2.cxx:120
 AliITStrackV2.cxx:121
 AliITStrackV2.cxx:122
 AliITStrackV2.cxx:123
 AliITStrackV2.cxx:124
 AliITStrackV2.cxx:125
 AliITStrackV2.cxx:126
 AliITStrackV2.cxx:127
 AliITStrackV2.cxx:128
 AliITStrackV2.cxx:129
 AliITStrackV2.cxx:130
 AliITStrackV2.cxx:131
 AliITStrackV2.cxx:132
 AliITStrackV2.cxx:133
 AliITStrackV2.cxx:134
 AliITStrackV2.cxx:135
 AliITStrackV2.cxx:136
 AliITStrackV2.cxx:137
 AliITStrackV2.cxx:138
 AliITStrackV2.cxx:139
 AliITStrackV2.cxx:140
 AliITStrackV2.cxx:141
 AliITStrackV2.cxx:142
 AliITStrackV2.cxx:143
 AliITStrackV2.cxx:144
 AliITStrackV2.cxx:145
 AliITStrackV2.cxx:146
 AliITStrackV2.cxx:147
 AliITStrackV2.cxx:148
 AliITStrackV2.cxx:149
 AliITStrackV2.cxx:150
 AliITStrackV2.cxx:151
 AliITStrackV2.cxx:152
 AliITStrackV2.cxx:153
 AliITStrackV2.cxx:154
 AliITStrackV2.cxx:155
 AliITStrackV2.cxx:156
 AliITStrackV2.cxx:157
 AliITStrackV2.cxx:158
 AliITStrackV2.cxx:159
 AliITStrackV2.cxx:160
 AliITStrackV2.cxx:161
 AliITStrackV2.cxx:162
 AliITStrackV2.cxx:163
 AliITStrackV2.cxx:164
 AliITStrackV2.cxx:165
 AliITStrackV2.cxx:166
 AliITStrackV2.cxx:167
 AliITStrackV2.cxx:168
 AliITStrackV2.cxx:169
 AliITStrackV2.cxx:170
 AliITStrackV2.cxx:171
 AliITStrackV2.cxx:172
 AliITStrackV2.cxx:173
 AliITStrackV2.cxx:174
 AliITStrackV2.cxx:175
 AliITStrackV2.cxx:176
 AliITStrackV2.cxx:177
 AliITStrackV2.cxx:178
 AliITStrackV2.cxx:179
 AliITStrackV2.cxx:180
 AliITStrackV2.cxx:181
 AliITStrackV2.cxx:182
 AliITStrackV2.cxx:183
 AliITStrackV2.cxx:184
 AliITStrackV2.cxx:185
 AliITStrackV2.cxx:186
 AliITStrackV2.cxx:187
 AliITStrackV2.cxx:188
 AliITStrackV2.cxx:189
 AliITStrackV2.cxx:190
 AliITStrackV2.cxx:191
 AliITStrackV2.cxx:192
 AliITStrackV2.cxx:193
 AliITStrackV2.cxx:194
 AliITStrackV2.cxx:195
 AliITStrackV2.cxx:196
 AliITStrackV2.cxx:197
 AliITStrackV2.cxx:198
 AliITStrackV2.cxx:199
 AliITStrackV2.cxx:200
 AliITStrackV2.cxx:201
 AliITStrackV2.cxx:202
 AliITStrackV2.cxx:203
 AliITStrackV2.cxx:204
 AliITStrackV2.cxx:205
 AliITStrackV2.cxx:206
 AliITStrackV2.cxx:207
 AliITStrackV2.cxx:208
 AliITStrackV2.cxx:209
 AliITStrackV2.cxx:210
 AliITStrackV2.cxx:211
 AliITStrackV2.cxx:212
 AliITStrackV2.cxx:213
 AliITStrackV2.cxx:214
 AliITStrackV2.cxx:215
 AliITStrackV2.cxx:216
 AliITStrackV2.cxx:217
 AliITStrackV2.cxx:218
 AliITStrackV2.cxx:219
 AliITStrackV2.cxx:220
 AliITStrackV2.cxx:221
 AliITStrackV2.cxx:222
 AliITStrackV2.cxx:223
 AliITStrackV2.cxx:224
 AliITStrackV2.cxx:225
 AliITStrackV2.cxx:226
 AliITStrackV2.cxx:227
 AliITStrackV2.cxx:228
 AliITStrackV2.cxx:229
 AliITStrackV2.cxx:230
 AliITStrackV2.cxx:231
 AliITStrackV2.cxx:232
 AliITStrackV2.cxx:233
 AliITStrackV2.cxx:234
 AliITStrackV2.cxx:235
 AliITStrackV2.cxx:236
 AliITStrackV2.cxx:237
 AliITStrackV2.cxx:238
 AliITStrackV2.cxx:239
 AliITStrackV2.cxx:240
 AliITStrackV2.cxx:241
 AliITStrackV2.cxx:242
 AliITStrackV2.cxx:243
 AliITStrackV2.cxx:244
 AliITStrackV2.cxx:245
 AliITStrackV2.cxx:246
 AliITStrackV2.cxx:247
 AliITStrackV2.cxx:248
 AliITStrackV2.cxx:249
 AliITStrackV2.cxx:250
 AliITStrackV2.cxx:251
 AliITStrackV2.cxx:252
 AliITStrackV2.cxx:253
 AliITStrackV2.cxx:254
 AliITStrackV2.cxx:255
 AliITStrackV2.cxx:256
 AliITStrackV2.cxx:257
 AliITStrackV2.cxx:258
 AliITStrackV2.cxx:259
 AliITStrackV2.cxx:260
 AliITStrackV2.cxx:261
 AliITStrackV2.cxx:262
 AliITStrackV2.cxx:263
 AliITStrackV2.cxx:264
 AliITStrackV2.cxx:265
 AliITStrackV2.cxx:266
 AliITStrackV2.cxx:267
 AliITStrackV2.cxx:268
 AliITStrackV2.cxx:269
 AliITStrackV2.cxx:270
 AliITStrackV2.cxx:271
 AliITStrackV2.cxx:272
 AliITStrackV2.cxx:273
 AliITStrackV2.cxx:274
 AliITStrackV2.cxx:275
 AliITStrackV2.cxx:276
 AliITStrackV2.cxx:277
 AliITStrackV2.cxx:278
 AliITStrackV2.cxx:279
 AliITStrackV2.cxx:280
 AliITStrackV2.cxx:281
 AliITStrackV2.cxx:282
 AliITStrackV2.cxx:283
 AliITStrackV2.cxx:284
 AliITStrackV2.cxx:285
 AliITStrackV2.cxx:286
 AliITStrackV2.cxx:287
 AliITStrackV2.cxx:288
 AliITStrackV2.cxx:289
 AliITStrackV2.cxx:290
 AliITStrackV2.cxx:291
 AliITStrackV2.cxx:292
 AliITStrackV2.cxx:293
 AliITStrackV2.cxx:294
 AliITStrackV2.cxx:295
 AliITStrackV2.cxx:296
 AliITStrackV2.cxx:297
 AliITStrackV2.cxx:298
 AliITStrackV2.cxx:299
 AliITStrackV2.cxx:300
 AliITStrackV2.cxx:301
 AliITStrackV2.cxx:302
 AliITStrackV2.cxx:303
 AliITStrackV2.cxx:304
 AliITStrackV2.cxx:305
 AliITStrackV2.cxx:306
 AliITStrackV2.cxx:307
 AliITStrackV2.cxx:308
 AliITStrackV2.cxx:309
 AliITStrackV2.cxx:310
 AliITStrackV2.cxx:311
 AliITStrackV2.cxx:312
 AliITStrackV2.cxx:313
 AliITStrackV2.cxx:314
 AliITStrackV2.cxx:315
 AliITStrackV2.cxx:316
 AliITStrackV2.cxx:317
 AliITStrackV2.cxx:318
 AliITStrackV2.cxx:319
 AliITStrackV2.cxx:320
 AliITStrackV2.cxx:321
 AliITStrackV2.cxx:322
 AliITStrackV2.cxx:323
 AliITStrackV2.cxx:324
 AliITStrackV2.cxx:325
 AliITStrackV2.cxx:326
 AliITStrackV2.cxx:327
 AliITStrackV2.cxx:328
 AliITStrackV2.cxx:329
 AliITStrackV2.cxx:330
 AliITStrackV2.cxx:331
 AliITStrackV2.cxx:332
 AliITStrackV2.cxx:333
 AliITStrackV2.cxx:334
 AliITStrackV2.cxx:335
 AliITStrackV2.cxx:336
 AliITStrackV2.cxx:337
 AliITStrackV2.cxx:338
 AliITStrackV2.cxx:339
 AliITStrackV2.cxx:340
 AliITStrackV2.cxx:341
 AliITStrackV2.cxx:342
 AliITStrackV2.cxx:343
 AliITStrackV2.cxx:344
 AliITStrackV2.cxx:345
 AliITStrackV2.cxx:346
 AliITStrackV2.cxx:347
 AliITStrackV2.cxx:348
 AliITStrackV2.cxx:349
 AliITStrackV2.cxx:350
 AliITStrackV2.cxx:351
 AliITStrackV2.cxx:352
 AliITStrackV2.cxx:353
 AliITStrackV2.cxx:354
 AliITStrackV2.cxx:355
 AliITStrackV2.cxx:356
 AliITStrackV2.cxx:357
 AliITStrackV2.cxx:358
 AliITStrackV2.cxx:359
 AliITStrackV2.cxx:360
 AliITStrackV2.cxx:361
 AliITStrackV2.cxx:362
 AliITStrackV2.cxx:363
 AliITStrackV2.cxx:364
 AliITStrackV2.cxx:365
 AliITStrackV2.cxx:366
 AliITStrackV2.cxx:367
 AliITStrackV2.cxx:368
 AliITStrackV2.cxx:369
 AliITStrackV2.cxx:370
 AliITStrackV2.cxx:371
 AliITStrackV2.cxx:372
 AliITStrackV2.cxx:373
 AliITStrackV2.cxx:374
 AliITStrackV2.cxx:375
 AliITStrackV2.cxx:376
 AliITStrackV2.cxx:377
 AliITStrackV2.cxx:378
 AliITStrackV2.cxx:379
 AliITStrackV2.cxx:380
 AliITStrackV2.cxx:381
 AliITStrackV2.cxx:382
 AliITStrackV2.cxx:383
 AliITStrackV2.cxx:384
 AliITStrackV2.cxx:385
 AliITStrackV2.cxx:386
 AliITStrackV2.cxx:387
 AliITStrackV2.cxx:388
 AliITStrackV2.cxx:389
 AliITStrackV2.cxx:390
 AliITStrackV2.cxx:391
 AliITStrackV2.cxx:392
 AliITStrackV2.cxx:393
 AliITStrackV2.cxx:394
 AliITStrackV2.cxx:395
 AliITStrackV2.cxx:396
 AliITStrackV2.cxx:397
 AliITStrackV2.cxx:398
 AliITStrackV2.cxx:399
 AliITStrackV2.cxx:400
 AliITStrackV2.cxx:401
 AliITStrackV2.cxx:402
 AliITStrackV2.cxx:403
 AliITStrackV2.cxx:404
 AliITStrackV2.cxx:405
 AliITStrackV2.cxx:406
 AliITStrackV2.cxx:407
 AliITStrackV2.cxx:408
 AliITStrackV2.cxx:409
 AliITStrackV2.cxx:410
 AliITStrackV2.cxx:411
 AliITStrackV2.cxx:412
 AliITStrackV2.cxx:413
 AliITStrackV2.cxx:414
 AliITStrackV2.cxx:415
 AliITStrackV2.cxx:416
 AliITStrackV2.cxx:417
 AliITStrackV2.cxx:418
 AliITStrackV2.cxx:419
 AliITStrackV2.cxx:420
 AliITStrackV2.cxx:421
 AliITStrackV2.cxx:422
 AliITStrackV2.cxx:423
 AliITStrackV2.cxx:424
 AliITStrackV2.cxx:425
 AliITStrackV2.cxx:426
 AliITStrackV2.cxx:427
 AliITStrackV2.cxx:428
 AliITStrackV2.cxx:429
 AliITStrackV2.cxx:430
 AliITStrackV2.cxx:431
 AliITStrackV2.cxx:432
 AliITStrackV2.cxx:433
 AliITStrackV2.cxx:434
 AliITStrackV2.cxx:435
 AliITStrackV2.cxx:436
 AliITStrackV2.cxx:437
 AliITStrackV2.cxx:438
 AliITStrackV2.cxx:439
 AliITStrackV2.cxx:440
 AliITStrackV2.cxx:441
 AliITStrackV2.cxx:442
 AliITStrackV2.cxx:443
 AliITStrackV2.cxx:444
 AliITStrackV2.cxx:445
 AliITStrackV2.cxx:446
 AliITStrackV2.cxx:447
 AliITStrackV2.cxx:448
 AliITStrackV2.cxx:449
 AliITStrackV2.cxx:450
 AliITStrackV2.cxx:451
 AliITStrackV2.cxx:452
 AliITStrackV2.cxx:453
 AliITStrackV2.cxx:454
 AliITStrackV2.cxx:455
 AliITStrackV2.cxx:456
 AliITStrackV2.cxx:457
 AliITStrackV2.cxx:458
 AliITStrackV2.cxx:459
 AliITStrackV2.cxx:460
 AliITStrackV2.cxx:461
 AliITStrackV2.cxx:462
 AliITStrackV2.cxx:463
 AliITStrackV2.cxx:464
 AliITStrackV2.cxx:465
 AliITStrackV2.cxx:466
 AliITStrackV2.cxx:467
 AliITStrackV2.cxx:468
 AliITStrackV2.cxx:469
 AliITStrackV2.cxx:470
 AliITStrackV2.cxx:471
 AliITStrackV2.cxx:472
 AliITStrackV2.cxx:473
 AliITStrackV2.cxx:474
 AliITStrackV2.cxx:475
 AliITStrackV2.cxx:476
 AliITStrackV2.cxx:477
 AliITStrackV2.cxx:478
 AliITStrackV2.cxx:479
 AliITStrackV2.cxx:480
 AliITStrackV2.cxx:481
 AliITStrackV2.cxx:482
 AliITStrackV2.cxx:483
 AliITStrackV2.cxx:484
 AliITStrackV2.cxx:485
 AliITStrackV2.cxx:486
 AliITStrackV2.cxx:487
 AliITStrackV2.cxx:488
 AliITStrackV2.cxx:489
 AliITStrackV2.cxx:490
 AliITStrackV2.cxx:491
 AliITStrackV2.cxx:492
 AliITStrackV2.cxx:493
 AliITStrackV2.cxx:494
 AliITStrackV2.cxx:495
 AliITStrackV2.cxx:496
 AliITStrackV2.cxx:497
 AliITStrackV2.cxx:498
 AliITStrackV2.cxx:499
 AliITStrackV2.cxx:500
 AliITStrackV2.cxx:501
 AliITStrackV2.cxx:502
 AliITStrackV2.cxx:503
 AliITStrackV2.cxx:504
 AliITStrackV2.cxx:505
 AliITStrackV2.cxx:506
 AliITStrackV2.cxx:507
 AliITStrackV2.cxx:508
 AliITStrackV2.cxx:509
 AliITStrackV2.cxx:510
 AliITStrackV2.cxx:511
 AliITStrackV2.cxx:512
 AliITStrackV2.cxx:513
 AliITStrackV2.cxx:514
 AliITStrackV2.cxx:515
 AliITStrackV2.cxx:516
 AliITStrackV2.cxx:517
 AliITStrackV2.cxx:518
 AliITStrackV2.cxx:519
 AliITStrackV2.cxx:520
 AliITStrackV2.cxx:521
 AliITStrackV2.cxx:522
 AliITStrackV2.cxx:523
 AliITStrackV2.cxx:524
 AliITStrackV2.cxx:525
 AliITStrackV2.cxx:526
 AliITStrackV2.cxx:527
 AliITStrackV2.cxx:528
 AliITStrackV2.cxx:529
 AliITStrackV2.cxx:530
 AliITStrackV2.cxx:531
 AliITStrackV2.cxx:532
 AliITStrackV2.cxx:533
 AliITStrackV2.cxx:534
 AliITStrackV2.cxx:535
 AliITStrackV2.cxx:536
 AliITStrackV2.cxx:537
 AliITStrackV2.cxx:538
 AliITStrackV2.cxx:539
 AliITStrackV2.cxx:540
 AliITStrackV2.cxx:541
 AliITStrackV2.cxx:542
 AliITStrackV2.cxx:543
 AliITStrackV2.cxx:544
 AliITStrackV2.cxx:545
 AliITStrackV2.cxx:546
 AliITStrackV2.cxx:547
 AliITStrackV2.cxx:548
 AliITStrackV2.cxx:549
 AliITStrackV2.cxx:550
 AliITStrackV2.cxx:551
 AliITStrackV2.cxx:552
 AliITStrackV2.cxx:553
 AliITStrackV2.cxx:554
 AliITStrackV2.cxx:555
 AliITStrackV2.cxx:556
 AliITStrackV2.cxx:557
 AliITStrackV2.cxx:558
 AliITStrackV2.cxx:559
 AliITStrackV2.cxx:560
 AliITStrackV2.cxx:561
 AliITStrackV2.cxx:562
 AliITStrackV2.cxx:563
 AliITStrackV2.cxx:564
 AliITStrackV2.cxx:565
 AliITStrackV2.cxx:566
 AliITStrackV2.cxx:567
 AliITStrackV2.cxx:568
 AliITStrackV2.cxx:569
 AliITStrackV2.cxx:570
 AliITStrackV2.cxx:571
 AliITStrackV2.cxx:572
 AliITStrackV2.cxx:573
 AliITStrackV2.cxx:574
 AliITStrackV2.cxx:575
 AliITStrackV2.cxx:576
 AliITStrackV2.cxx:577
 AliITStrackV2.cxx:578
 AliITStrackV2.cxx:579
 AliITStrackV2.cxx:580
 AliITStrackV2.cxx:581
 AliITStrackV2.cxx:582
 AliITStrackV2.cxx:583
 AliITStrackV2.cxx:584
 AliITStrackV2.cxx:585
 AliITStrackV2.cxx:586
 AliITStrackV2.cxx:587
 AliITStrackV2.cxx:588
 AliITStrackV2.cxx:589
 AliITStrackV2.cxx:590
 AliITStrackV2.cxx:591
 AliITStrackV2.cxx:592
 AliITStrackV2.cxx:593
 AliITStrackV2.cxx:594
 AliITStrackV2.cxx:595
 AliITStrackV2.cxx:596
 AliITStrackV2.cxx:597
 AliITStrackV2.cxx:598
 AliITStrackV2.cxx:599
 AliITStrackV2.cxx:600
 AliITStrackV2.cxx:601
 AliITStrackV2.cxx:602
 AliITStrackV2.cxx:603
 AliITStrackV2.cxx:604
 AliITStrackV2.cxx:605
 AliITStrackV2.cxx:606
 AliITStrackV2.cxx:607
 AliITStrackV2.cxx:608
 AliITStrackV2.cxx:609
 AliITStrackV2.cxx:610
 AliITStrackV2.cxx:611
 AliITStrackV2.cxx:612
 AliITStrackV2.cxx:613
 AliITStrackV2.cxx:614
 AliITStrackV2.cxx:615
 AliITStrackV2.cxx:616
 AliITStrackV2.cxx:617
 AliITStrackV2.cxx:618
 AliITStrackV2.cxx:619
 AliITStrackV2.cxx:620
 AliITStrackV2.cxx:621
 AliITStrackV2.cxx:622
 AliITStrackV2.cxx:623
 AliITStrackV2.cxx:624
 AliITStrackV2.cxx:625
 AliITStrackV2.cxx:626
 AliITStrackV2.cxx:627
 AliITStrackV2.cxx:628
 AliITStrackV2.cxx:629
 AliITStrackV2.cxx:630
 AliITStrackV2.cxx:631
 AliITStrackV2.cxx:632
 AliITStrackV2.cxx:633
 AliITStrackV2.cxx:634
 AliITStrackV2.cxx:635
 AliITStrackV2.cxx:636
 AliITStrackV2.cxx:637
 AliITStrackV2.cxx:638
 AliITStrackV2.cxx:639
 AliITStrackV2.cxx:640
 AliITStrackV2.cxx:641
 AliITStrackV2.cxx:642
 AliITStrackV2.cxx:643
 AliITStrackV2.cxx:644
 AliITStrackV2.cxx:645
 AliITStrackV2.cxx:646
 AliITStrackV2.cxx:647
 AliITStrackV2.cxx:648
 AliITStrackV2.cxx:649
 AliITStrackV2.cxx:650
 AliITStrackV2.cxx:651
 AliITStrackV2.cxx:652
 AliITStrackV2.cxx:653
 AliITStrackV2.cxx:654
 AliITStrackV2.cxx:655
 AliITStrackV2.cxx:656
 AliITStrackV2.cxx:657
 AliITStrackV2.cxx:658
 AliITStrackV2.cxx:659
 AliITStrackV2.cxx:660
 AliITStrackV2.cxx:661
 AliITStrackV2.cxx:662
 AliITStrackV2.cxx:663
 AliITStrackV2.cxx:664
 AliITStrackV2.cxx:665
 AliITStrackV2.cxx:666
 AliITStrackV2.cxx:667
 AliITStrackV2.cxx:668
 AliITStrackV2.cxx:669
 AliITStrackV2.cxx:670
 AliITStrackV2.cxx:671
 AliITStrackV2.cxx:672
 AliITStrackV2.cxx:673
 AliITStrackV2.cxx:674
 AliITStrackV2.cxx:675
 AliITStrackV2.cxx:676
 AliITStrackV2.cxx:677
 AliITStrackV2.cxx:678
 AliITStrackV2.cxx:679
 AliITStrackV2.cxx:680
 AliITStrackV2.cxx:681
 AliITStrackV2.cxx:682
 AliITStrackV2.cxx:683
 AliITStrackV2.cxx:684
 AliITStrackV2.cxx:685
 AliITStrackV2.cxx:686
 AliITStrackV2.cxx:687
 AliITStrackV2.cxx:688
 AliITStrackV2.cxx:689
 AliITStrackV2.cxx:690
 AliITStrackV2.cxx:691
 AliITStrackV2.cxx:692
 AliITStrackV2.cxx:693
 AliITStrackV2.cxx:694
 AliITStrackV2.cxx:695
 AliITStrackV2.cxx:696
 AliITStrackV2.cxx:697
 AliITStrackV2.cxx:698
 AliITStrackV2.cxx:699
 AliITStrackV2.cxx:700
 AliITStrackV2.cxx:701
 AliITStrackV2.cxx:702
 AliITStrackV2.cxx:703
 AliITStrackV2.cxx:704
 AliITStrackV2.cxx:705
 AliITStrackV2.cxx:706
 AliITStrackV2.cxx:707
 AliITStrackV2.cxx:708
 AliITStrackV2.cxx:709
 AliITStrackV2.cxx:710
 AliITStrackV2.cxx:711
 AliITStrackV2.cxx:712
 AliITStrackV2.cxx:713
 AliITStrackV2.cxx:714
 AliITStrackV2.cxx:715
 AliITStrackV2.cxx:716
 AliITStrackV2.cxx:717
 AliITStrackV2.cxx:718
 AliITStrackV2.cxx:719
 AliITStrackV2.cxx:720
 AliITStrackV2.cxx:721
 AliITStrackV2.cxx:722
 AliITStrackV2.cxx:723
 AliITStrackV2.cxx:724
 AliITStrackV2.cxx:725
 AliITStrackV2.cxx:726
 AliITStrackV2.cxx:727
 AliITStrackV2.cxx:728
 AliITStrackV2.cxx:729
 AliITStrackV2.cxx:730
 AliITStrackV2.cxx:731
 AliITStrackV2.cxx:732
 AliITStrackV2.cxx:733
 AliITStrackV2.cxx:734
 AliITStrackV2.cxx:735
 AliITStrackV2.cxx:736
 AliITStrackV2.cxx:737
 AliITStrackV2.cxx:738
 AliITStrackV2.cxx:739
 AliITStrackV2.cxx:740
 AliITStrackV2.cxx:741
 AliITStrackV2.cxx:742
 AliITStrackV2.cxx:743
 AliITStrackV2.cxx:744
 AliITStrackV2.cxx:745
 AliITStrackV2.cxx:746
 AliITStrackV2.cxx:747
 AliITStrackV2.cxx:748
 AliITStrackV2.cxx:749
 AliITStrackV2.cxx:750
 AliITStrackV2.cxx:751
 AliITStrackV2.cxx:752
 AliITStrackV2.cxx:753
 AliITStrackV2.cxx:754
 AliITStrackV2.cxx:755
 AliITStrackV2.cxx:756
 AliITStrackV2.cxx:757
 AliITStrackV2.cxx:758
 AliITStrackV2.cxx:759
 AliITStrackV2.cxx:760
 AliITStrackV2.cxx:761
 AliITStrackV2.cxx:762
 AliITStrackV2.cxx:763
 AliITStrackV2.cxx:764
 AliITStrackV2.cxx:765
 AliITStrackV2.cxx:766
 AliITStrackV2.cxx:767
 AliITStrackV2.cxx:768
 AliITStrackV2.cxx:769
 AliITStrackV2.cxx:770
 AliITStrackV2.cxx:771
 AliITStrackV2.cxx:772
 AliITStrackV2.cxx:773
 AliITStrackV2.cxx:774
 AliITStrackV2.cxx:775
 AliITStrackV2.cxx:776
 AliITStrackV2.cxx:777
 AliITStrackV2.cxx:778
 AliITStrackV2.cxx:779
 AliITStrackV2.cxx:780
 AliITStrackV2.cxx:781
 AliITStrackV2.cxx:782
 AliITStrackV2.cxx:783
 AliITStrackV2.cxx:784
 AliITStrackV2.cxx:785
 AliITStrackV2.cxx:786
 AliITStrackV2.cxx:787
 AliITStrackV2.cxx:788