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 derived class for track residuals
//   based on linear chi2 minimization (in approximation of
//   small alignment angles and translations)
//
//-----------------------------------------------------------------

#include <TMath.h>
#include <TGeoMatrix.h>

#include "AliLog.h"
#include "AliAlignObj.h"
#include "AliTrackPointArray.h"
#include "AliTrackResidualsFast.h"

#include <TMatrixDSym.h>
#include <TMatrixDSymEigen.h>

ClassImp(AliTrackResidualsFast)

//______________________________________________________________________________
AliTrackResidualsFast::AliTrackResidualsFast():
  AliTrackResiduals(),
  fSumR(0)
{
  // Default constructor
  for (Int_t i = 0; i < 27; i++) fSum[i] = 0;
}

//______________________________________________________________________________
AliTrackResidualsFast::AliTrackResidualsFast(Int_t ntracks):
  AliTrackResiduals(ntracks),
  fSumR(0)
{
  // Constructor
  for (Int_t i = 0; i < 27; i++) fSum[i] = 0;
}
 
//______________________________________________________________________________
AliTrackResidualsFast::AliTrackResidualsFast(const AliTrackResidualsFast &res):
  AliTrackResiduals(res),
  fSumR(res.fSumR)
{
  // Copy constructor
  memcpy(fSum,res.fSum,27*sizeof(Double_t));
}

//______________________________________________________________________________
AliTrackResidualsFast &AliTrackResidualsFast::operator= (const AliTrackResidualsFast& res)
{
  //
  // Assignment operator
  //
  if(this != &res) {
    AliTrackResiduals::operator=(res);
    memcpy(fSum,res.fSum,27*sizeof(Double_t));
    fSumR = res.fSumR;
  }

 return *this;
}

//______________________________________________________________________________
Bool_t AliTrackResidualsFast::Minimize()
{
  // Implementation of fast linear Chi2
  // based minimization of track residuals sum

  //  if(fBFixed[0]||fBFixed[1]||fBFixed[2]||fBFixed[3]||fBFixed[4]||fBFixed[5])
  //    AliError("Cannot yet fix parameters in this minimizer");


  for (Int_t i = 0; i < 27; i++) fSum[i] = 0;
  fSumR = 0;

  AliTrackPoint p1,p2;

  for (Int_t itrack = 0; itrack < fLast; itrack++) {
    if (!fVolArray[itrack] || !fTrackArray[itrack]) continue;
    for (Int_t ipoint = 0; ipoint < fVolArray[itrack]->GetNPoints(); ipoint++) {
      fVolArray[itrack]->GetPoint(p1,ipoint);
      fTrackArray[itrack]->GetPoint(p2,ipoint);
      AddPoints(p1,p2);
    }
  }

  return Update();

  // debug info
//   Float_t chi2 = 0;
//   for (Int_t itrack = 0; itrack < fLast; itrack++) {
//     if (!fVolArray[itrack] || !fTrackArray[itrack]) continue;
//     for (Int_t ipoint = 0; ipoint < fVolArray[itrack]->GetNPoints(); ipoint++) {
//       fVolArray[itrack]->GetPoint(p1,ipoint);
//       fAlignObj->Transform(p1);
//       fTrackArray[itrack]->GetPoint(p2,ipoint);
//       Float_t residual = p2.GetResidual(p1,kFALSE);
//       chi2 += residual;
//     }
//   }
//   printf("Final chi2 = %f\n",chi2);
}

//______________________________________________________________________________
void AliTrackResidualsFast::AddPoints(AliTrackPoint &p, AliTrackPoint &pprime)
{
  // Update the sums used for
  // the linear chi2 minimization
  Float_t xyz[3],xyzp[3];
  Float_t cov[6],covp[6];
  p.GetXYZ(xyz,cov); pprime.GetXYZ(xyzp,covp);
  TMatrixDSym mcov(3);
  mcov(0,0) = cov[0]; mcov(0,1) = cov[1]; mcov(0,2) = cov[2];
  mcov(1,0) = cov[1]; mcov(1,1) = cov[3]; mcov(1,2) = cov[4];
  mcov(2,0) = cov[2]; mcov(2,1) = cov[4]; mcov(2,2) = cov[5];
  TMatrixDSym mcovp(3);
  mcovp(0,0) = covp[0]; mcovp(0,1) = covp[1]; mcovp(0,2) = covp[2];
  mcovp(1,0) = covp[1]; mcovp(1,1) = covp[3]; mcovp(1,2) = covp[4];
  mcovp(2,0) = covp[2]; mcovp(2,1) = covp[4]; mcovp(2,2) = covp[5];
  TMatrixDSym msum = mcov + mcovp;


  msum.Invert();


  if (!msum.IsValid()) return;

  TMatrixD        sums(3,1);
  sums(0,0) = (xyzp[0]-xyz[0]); 
  sums(1,0) = (xyzp[1]-xyz[1]);
  sums(2,0) = (xyzp[2]-xyz[2]);   
  TMatrixD sumst = sums.T(); sums.T();

  TMatrixD        mf(3,6);
  mf(0,0) = 1;      mf(1,0) = 0;       mf(2,0) = 0;
  mf(0,1) = 0;      mf(1,1) = 1;       mf(2,1) = 0;
  mf(0,2) = 0;      mf(1,2) = 0;       mf(2,2) = 1;
  mf(0,3) = 0;      mf(1,3) = -xyz[2]; mf(2,3) = xyz[1];
  mf(0,4) = xyz[2]; mf(1,4) = 0;       mf(2,4) =-xyz[0];
  mf(0,5) =-xyz[1]; mf(1,5) = xyz[0];  mf(2,5) = 0;

  for(Int_t j=0;j<6;j++){
    if(fBFixed[j]==kTRUE){
      mf(0,j)=0.;mf(1,j)=0.;mf(2,j)=0.;
    }
  }

  TMatrixD        mft = mf.T(); mf.T();
  TMatrixD sums2 = mft * msum * sums;

  TMatrixD smatrix = mft * msum * mf;

  fSum[0] += smatrix(0,0);
  fSum[1] += smatrix(0,1);
  fSum[2] += smatrix(0,2);
  fSum[3] += smatrix(0,3);
  fSum[4] += smatrix(0,4);
  fSum[5] += smatrix(0,5);
  fSum[6] += smatrix(1,1);
  fSum[7] += smatrix(1,2);
  fSum[8] += smatrix(1,3);
  fSum[9] += smatrix(1,4);
  fSum[10]+= smatrix(1,5);
  fSum[11]+= smatrix(2,2);
  fSum[12]+= smatrix(2,3);
  fSum[13]+= smatrix(2,4);
  fSum[14]+= smatrix(2,5);
  fSum[15]+= smatrix(3,3);
  fSum[16]+= smatrix(3,4);
  fSum[17]+= smatrix(3,5);
  fSum[18]+= smatrix(4,4);
  fSum[19]+= smatrix(4,5);
  fSum[20]+= smatrix(5,5);
  fSum[21] += sums2(0,0);
  fSum[22] += sums2(1,0);
  fSum[23] += sums2(2,0);
  fSum[24] += sums2(3,0);
  fSum[25] += sums2(4,0);
  fSum[26] += sums2(5,0);

  TMatrixD  tmp = sumst * msum * sums;
  fSumR += tmp(0,0);

  fNdf += 3;
}

//______________________________________________________________________________
Bool_t AliTrackResidualsFast::Update()
{
  // Find the alignment parameters
  // by using the already accumulated
  // sums
  TMatrixDSym     smatrix(6);
  TMatrixD        sums(1,6);

  smatrix(0,0) = fSum[0];
  smatrix(0,1) = smatrix(1,0) = fSum[1];
  smatrix(0,2) = smatrix(2,0) = fSum[2];
  smatrix(0,3) = smatrix(3,0) = fSum[3];
  smatrix(0,4) = smatrix(4,0) = fSum[4];
  smatrix(0,5) = smatrix(5,0) = fSum[5];
  smatrix(1,1) = fSum[6];
  smatrix(1,2) = smatrix(2,1) = fSum[7];
  smatrix(1,3) = smatrix(3,1) = fSum[8];
  smatrix(1,4) = smatrix(4,1) = fSum[9];
  smatrix(1,5) = smatrix(5,1) = fSum[10];
  smatrix(2,2) = fSum[11];
  smatrix(2,3) = smatrix(3,2) = fSum[12];
  smatrix(2,4) = smatrix(4,2) = fSum[13];
  smatrix(2,5) = smatrix(5,2) = fSum[14];
  smatrix(3,3) = fSum[15];
  smatrix(3,4) = smatrix(4,3) = fSum[16];
  smatrix(3,5) = smatrix(5,3) = fSum[17];
  smatrix(4,4) = fSum[18];
  smatrix(4,5) = smatrix(5,4) = fSum[19];
  smatrix(5,5) = fSum[20];

  sums(0,0)    = fSum[21]; sums(0,1) = fSum[22]; sums(0,2) = fSum[23];
  sums(0,3)    = fSum[24]; sums(0,4) = fSum[25]; sums(0,5) = fSum[26];

  
  Int_t fixedparamat[6]={0,0,0,0,0,0}; 
  const Int_t unfixedparam=GetNFreeParam();
  Int_t position[6]={0,0,0,0,0,0};
  Int_t last=0;//position is of size 6 but only unfiexedparam indeces will be used
  
  if(fBFixed[0]==kTRUE){
    fixedparamat[0]=1;
  }
  else {
    position[0]=0;
    last++;
  }
  
  for(Int_t j=1;j<6;j++){
    if(fBFixed[j]==kTRUE){
      fixedparamat[j]=fixedparamat[j-1]+1;
    }
    else {
      fixedparamat[j]=fixedparamat[j-1];
      position[last]=j;
      last++;
    }
  }

  TMatrixDSym smatrixRedu(unfixedparam);
  for(Int_t i=0;i<unfixedparam;i++){
    for(Int_t j=0;j<unfixedparam;j++){
      smatrixRedu(i,j)=smatrix(position[i],position[j]);
    }
  }
  
  //  smatrixRedu.Print();
  smatrixRedu.Invert();
  
  if (!smatrixRedu.IsValid()) {
    printf("Minimization Failed! \n");
    return kFALSE;
  }

  TMatrixDSym smatrixUp(6);
  for(Int_t i=0;i<6;i++){
    for(Int_t j=0;j<6;j++){
      if(fBFixed[i]==kTRUE||fBFixed[j]==kTRUE)smatrixUp(i,j)=0.;
      else smatrixUp(i,j)=smatrixRedu(i-fixedparamat[i],j-fixedparamat[j]);
    }
  }
  
  Double_t covmatrarray[21];
  
  for(Int_t i=0;i<6;i++){
      for(Int_t j=0;j<=i;j++){
	if(fBFixed[i]==kFALSE&&fBFixed[j]==kFALSE){
	  if(TMath::Abs(smatrixUp(i,j)/TMath::Sqrt(TMath::Abs(smatrixUp(i,i)*smatrixUp(j,j))))>1.01)printf("Too large Correlation number!\n");
	}
	covmatrarray[i*(i+1)/2+j]=smatrixUp(i,j);
      }
  }
    
  TMatrixD res = sums*smatrixUp;
  fAlignObj->SetPars(res(0,0),res(0,1),res(0,2),
		     TMath::RadToDeg()*res(0,3),
		     TMath::RadToDeg()*res(0,4),
		     TMath::RadToDeg()*res(0,5));
  
  fAlignObj->SetCorrMatrix(covmatrarray);
  TMatrixD  tmp = res*sums.T();
  fChi2 = fSumR - tmp(0,0);
  fNdf -= unfixedparam;
  
  return kTRUE;
}



 AliTrackResidualsFast.cxx:1
 AliTrackResidualsFast.cxx:2
 AliTrackResidualsFast.cxx:3
 AliTrackResidualsFast.cxx:4
 AliTrackResidualsFast.cxx:5
 AliTrackResidualsFast.cxx:6
 AliTrackResidualsFast.cxx:7
 AliTrackResidualsFast.cxx:8
 AliTrackResidualsFast.cxx:9
 AliTrackResidualsFast.cxx:10
 AliTrackResidualsFast.cxx:11
 AliTrackResidualsFast.cxx:12
 AliTrackResidualsFast.cxx:13
 AliTrackResidualsFast.cxx:14
 AliTrackResidualsFast.cxx:15
 AliTrackResidualsFast.cxx:16
 AliTrackResidualsFast.cxx:17
 AliTrackResidualsFast.cxx:18
 AliTrackResidualsFast.cxx:19
 AliTrackResidualsFast.cxx:20
 AliTrackResidualsFast.cxx:21
 AliTrackResidualsFast.cxx:22
 AliTrackResidualsFast.cxx:23
 AliTrackResidualsFast.cxx:24
 AliTrackResidualsFast.cxx:25
 AliTrackResidualsFast.cxx:26
 AliTrackResidualsFast.cxx:27
 AliTrackResidualsFast.cxx:28
 AliTrackResidualsFast.cxx:29
 AliTrackResidualsFast.cxx:30
 AliTrackResidualsFast.cxx:31
 AliTrackResidualsFast.cxx:32
 AliTrackResidualsFast.cxx:33
 AliTrackResidualsFast.cxx:34
 AliTrackResidualsFast.cxx:35
 AliTrackResidualsFast.cxx:36
 AliTrackResidualsFast.cxx:37
 AliTrackResidualsFast.cxx:38
 AliTrackResidualsFast.cxx:39
 AliTrackResidualsFast.cxx:40
 AliTrackResidualsFast.cxx:41
 AliTrackResidualsFast.cxx:42
 AliTrackResidualsFast.cxx:43
 AliTrackResidualsFast.cxx:44
 AliTrackResidualsFast.cxx:45
 AliTrackResidualsFast.cxx:46
 AliTrackResidualsFast.cxx:47
 AliTrackResidualsFast.cxx:48
 AliTrackResidualsFast.cxx:49
 AliTrackResidualsFast.cxx:50
 AliTrackResidualsFast.cxx:51
 AliTrackResidualsFast.cxx:52
 AliTrackResidualsFast.cxx:53
 AliTrackResidualsFast.cxx:54
 AliTrackResidualsFast.cxx:55
 AliTrackResidualsFast.cxx:56
 AliTrackResidualsFast.cxx:57
 AliTrackResidualsFast.cxx:58
 AliTrackResidualsFast.cxx:59
 AliTrackResidualsFast.cxx:60
 AliTrackResidualsFast.cxx:61
 AliTrackResidualsFast.cxx:62
 AliTrackResidualsFast.cxx:63
 AliTrackResidualsFast.cxx:64
 AliTrackResidualsFast.cxx:65
 AliTrackResidualsFast.cxx:66
 AliTrackResidualsFast.cxx:67
 AliTrackResidualsFast.cxx:68
 AliTrackResidualsFast.cxx:69
 AliTrackResidualsFast.cxx:70
 AliTrackResidualsFast.cxx:71
 AliTrackResidualsFast.cxx:72
 AliTrackResidualsFast.cxx:73
 AliTrackResidualsFast.cxx:74
 AliTrackResidualsFast.cxx:75
 AliTrackResidualsFast.cxx:76
 AliTrackResidualsFast.cxx:77
 AliTrackResidualsFast.cxx:78
 AliTrackResidualsFast.cxx:79
 AliTrackResidualsFast.cxx:80
 AliTrackResidualsFast.cxx:81
 AliTrackResidualsFast.cxx:82
 AliTrackResidualsFast.cxx:83
 AliTrackResidualsFast.cxx:84
 AliTrackResidualsFast.cxx:85
 AliTrackResidualsFast.cxx:86
 AliTrackResidualsFast.cxx:87
 AliTrackResidualsFast.cxx:88
 AliTrackResidualsFast.cxx:89
 AliTrackResidualsFast.cxx:90
 AliTrackResidualsFast.cxx:91
 AliTrackResidualsFast.cxx:92
 AliTrackResidualsFast.cxx:93
 AliTrackResidualsFast.cxx:94
 AliTrackResidualsFast.cxx:95
 AliTrackResidualsFast.cxx:96
 AliTrackResidualsFast.cxx:97
 AliTrackResidualsFast.cxx:98
 AliTrackResidualsFast.cxx:99
 AliTrackResidualsFast.cxx:100
 AliTrackResidualsFast.cxx:101
 AliTrackResidualsFast.cxx:102
 AliTrackResidualsFast.cxx:103
 AliTrackResidualsFast.cxx:104
 AliTrackResidualsFast.cxx:105
 AliTrackResidualsFast.cxx:106
 AliTrackResidualsFast.cxx:107
 AliTrackResidualsFast.cxx:108
 AliTrackResidualsFast.cxx:109
 AliTrackResidualsFast.cxx:110
 AliTrackResidualsFast.cxx:111
 AliTrackResidualsFast.cxx:112
 AliTrackResidualsFast.cxx:113
 AliTrackResidualsFast.cxx:114
 AliTrackResidualsFast.cxx:115
 AliTrackResidualsFast.cxx:116
 AliTrackResidualsFast.cxx:117
 AliTrackResidualsFast.cxx:118
 AliTrackResidualsFast.cxx:119
 AliTrackResidualsFast.cxx:120
 AliTrackResidualsFast.cxx:121
 AliTrackResidualsFast.cxx:122
 AliTrackResidualsFast.cxx:123
 AliTrackResidualsFast.cxx:124
 AliTrackResidualsFast.cxx:125
 AliTrackResidualsFast.cxx:126
 AliTrackResidualsFast.cxx:127
 AliTrackResidualsFast.cxx:128
 AliTrackResidualsFast.cxx:129
 AliTrackResidualsFast.cxx:130
 AliTrackResidualsFast.cxx:131
 AliTrackResidualsFast.cxx:132
 AliTrackResidualsFast.cxx:133
 AliTrackResidualsFast.cxx:134
 AliTrackResidualsFast.cxx:135
 AliTrackResidualsFast.cxx:136
 AliTrackResidualsFast.cxx:137
 AliTrackResidualsFast.cxx:138
 AliTrackResidualsFast.cxx:139
 AliTrackResidualsFast.cxx:140
 AliTrackResidualsFast.cxx:141
 AliTrackResidualsFast.cxx:142
 AliTrackResidualsFast.cxx:143
 AliTrackResidualsFast.cxx:144
 AliTrackResidualsFast.cxx:145
 AliTrackResidualsFast.cxx:146
 AliTrackResidualsFast.cxx:147
 AliTrackResidualsFast.cxx:148
 AliTrackResidualsFast.cxx:149
 AliTrackResidualsFast.cxx:150
 AliTrackResidualsFast.cxx:151
 AliTrackResidualsFast.cxx:152
 AliTrackResidualsFast.cxx:153
 AliTrackResidualsFast.cxx:154
 AliTrackResidualsFast.cxx:155
 AliTrackResidualsFast.cxx:156
 AliTrackResidualsFast.cxx:157
 AliTrackResidualsFast.cxx:158
 AliTrackResidualsFast.cxx:159
 AliTrackResidualsFast.cxx:160
 AliTrackResidualsFast.cxx:161
 AliTrackResidualsFast.cxx:162
 AliTrackResidualsFast.cxx:163
 AliTrackResidualsFast.cxx:164
 AliTrackResidualsFast.cxx:165
 AliTrackResidualsFast.cxx:166
 AliTrackResidualsFast.cxx:167
 AliTrackResidualsFast.cxx:168
 AliTrackResidualsFast.cxx:169
 AliTrackResidualsFast.cxx:170
 AliTrackResidualsFast.cxx:171
 AliTrackResidualsFast.cxx:172
 AliTrackResidualsFast.cxx:173
 AliTrackResidualsFast.cxx:174
 AliTrackResidualsFast.cxx:175
 AliTrackResidualsFast.cxx:176
 AliTrackResidualsFast.cxx:177
 AliTrackResidualsFast.cxx:178
 AliTrackResidualsFast.cxx:179
 AliTrackResidualsFast.cxx:180
 AliTrackResidualsFast.cxx:181
 AliTrackResidualsFast.cxx:182
 AliTrackResidualsFast.cxx:183
 AliTrackResidualsFast.cxx:184
 AliTrackResidualsFast.cxx:185
 AliTrackResidualsFast.cxx:186
 AliTrackResidualsFast.cxx:187
 AliTrackResidualsFast.cxx:188
 AliTrackResidualsFast.cxx:189
 AliTrackResidualsFast.cxx:190
 AliTrackResidualsFast.cxx:191
 AliTrackResidualsFast.cxx:192
 AliTrackResidualsFast.cxx:193
 AliTrackResidualsFast.cxx:194
 AliTrackResidualsFast.cxx:195
 AliTrackResidualsFast.cxx:196
 AliTrackResidualsFast.cxx:197
 AliTrackResidualsFast.cxx:198
 AliTrackResidualsFast.cxx:199
 AliTrackResidualsFast.cxx:200
 AliTrackResidualsFast.cxx:201
 AliTrackResidualsFast.cxx:202
 AliTrackResidualsFast.cxx:203
 AliTrackResidualsFast.cxx:204
 AliTrackResidualsFast.cxx:205
 AliTrackResidualsFast.cxx:206
 AliTrackResidualsFast.cxx:207
 AliTrackResidualsFast.cxx:208
 AliTrackResidualsFast.cxx:209
 AliTrackResidualsFast.cxx:210
 AliTrackResidualsFast.cxx:211
 AliTrackResidualsFast.cxx:212
 AliTrackResidualsFast.cxx:213
 AliTrackResidualsFast.cxx:214
 AliTrackResidualsFast.cxx:215
 AliTrackResidualsFast.cxx:216
 AliTrackResidualsFast.cxx:217
 AliTrackResidualsFast.cxx:218
 AliTrackResidualsFast.cxx:219
 AliTrackResidualsFast.cxx:220
 AliTrackResidualsFast.cxx:221
 AliTrackResidualsFast.cxx:222
 AliTrackResidualsFast.cxx:223
 AliTrackResidualsFast.cxx:224
 AliTrackResidualsFast.cxx:225
 AliTrackResidualsFast.cxx:226
 AliTrackResidualsFast.cxx:227
 AliTrackResidualsFast.cxx:228
 AliTrackResidualsFast.cxx:229
 AliTrackResidualsFast.cxx:230
 AliTrackResidualsFast.cxx:231
 AliTrackResidualsFast.cxx:232
 AliTrackResidualsFast.cxx:233
 AliTrackResidualsFast.cxx:234
 AliTrackResidualsFast.cxx:235
 AliTrackResidualsFast.cxx:236
 AliTrackResidualsFast.cxx:237
 AliTrackResidualsFast.cxx:238
 AliTrackResidualsFast.cxx:239
 AliTrackResidualsFast.cxx:240
 AliTrackResidualsFast.cxx:241
 AliTrackResidualsFast.cxx:242
 AliTrackResidualsFast.cxx:243
 AliTrackResidualsFast.cxx:244
 AliTrackResidualsFast.cxx:245
 AliTrackResidualsFast.cxx:246
 AliTrackResidualsFast.cxx:247
 AliTrackResidualsFast.cxx:248
 AliTrackResidualsFast.cxx:249
 AliTrackResidualsFast.cxx:250
 AliTrackResidualsFast.cxx:251
 AliTrackResidualsFast.cxx:252
 AliTrackResidualsFast.cxx:253
 AliTrackResidualsFast.cxx:254
 AliTrackResidualsFast.cxx:255
 AliTrackResidualsFast.cxx:256
 AliTrackResidualsFast.cxx:257
 AliTrackResidualsFast.cxx:258
 AliTrackResidualsFast.cxx:259
 AliTrackResidualsFast.cxx:260
 AliTrackResidualsFast.cxx:261
 AliTrackResidualsFast.cxx:262
 AliTrackResidualsFast.cxx:263
 AliTrackResidualsFast.cxx:264
 AliTrackResidualsFast.cxx:265
 AliTrackResidualsFast.cxx:266
 AliTrackResidualsFast.cxx:267
 AliTrackResidualsFast.cxx:268
 AliTrackResidualsFast.cxx:269
 AliTrackResidualsFast.cxx:270
 AliTrackResidualsFast.cxx:271
 AliTrackResidualsFast.cxx:272
 AliTrackResidualsFast.cxx:273
 AliTrackResidualsFast.cxx:274
 AliTrackResidualsFast.cxx:275
 AliTrackResidualsFast.cxx:276
 AliTrackResidualsFast.cxx:277
 AliTrackResidualsFast.cxx:278
 AliTrackResidualsFast.cxx:279
 AliTrackResidualsFast.cxx:280
 AliTrackResidualsFast.cxx:281
 AliTrackResidualsFast.cxx:282
 AliTrackResidualsFast.cxx:283
 AliTrackResidualsFast.cxx:284
 AliTrackResidualsFast.cxx:285
 AliTrackResidualsFast.cxx:286
 AliTrackResidualsFast.cxx:287
 AliTrackResidualsFast.cxx:288
 AliTrackResidualsFast.cxx:289
 AliTrackResidualsFast.cxx:290
 AliTrackResidualsFast.cxx:291
 AliTrackResidualsFast.cxx:292
 AliTrackResidualsFast.cxx:293
 AliTrackResidualsFast.cxx:294
 AliTrackResidualsFast.cxx:295
 AliTrackResidualsFast.cxx:296
 AliTrackResidualsFast.cxx:297
 AliTrackResidualsFast.cxx:298
 AliTrackResidualsFast.cxx:299
 AliTrackResidualsFast.cxx:300
 AliTrackResidualsFast.cxx:301
 AliTrackResidualsFast.cxx:302
 AliTrackResidualsFast.cxx:303
 AliTrackResidualsFast.cxx:304
 AliTrackResidualsFast.cxx:305
 AliTrackResidualsFast.cxx:306
 AliTrackResidualsFast.cxx:307
 AliTrackResidualsFast.cxx:308
 AliTrackResidualsFast.cxx:309
 AliTrackResidualsFast.cxx:310
 AliTrackResidualsFast.cxx:311
 AliTrackResidualsFast.cxx:312