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

/*  
  
  Unlinearities fitting:

  Model for Outer field cage:
  Unlinearities at the edge aproximated using two exponential decays.
 
  dz = dz0(r,z) +dr(r,z)*tan(theta) 
  dy =          +dr(r,z)*tan(phi)

   
  
    
*/

#include "AliTPCcalibDB.h"
#include "TLinearFitter.h"
#include "Riostream.h"
#include "TList.h"
#include "TMath.h"
#include "TCanvas.h"
#include "TFile.h"
#include "TF1.h"
#include "TVectorD.h"
#include "AliLog.h"
#include "AliTPCROC.h"
#include "AliTPCClusterParam.h"
#include "AliTPCPointCorrection.h"

AliTPCPointCorrection* AliTPCPointCorrection::fgInstance = 0;

ClassImp(AliTPCPointCorrection)

AliTPCPointCorrection::AliTPCPointCorrection():
  TNamed(),
  fParamsOutR(38),
  fParamsOutZ(38),
  fParamOutRVersion(0),
  fErrorsOutR(38),
  fErrorsOutZ(38),
  fParamOutZVersion(0),
  //
  fXIO(0),
  fXmiddle(0),
  fXquadrant(0),
  fArraySectorIntParam(36), // array of sector alignment parameters
  fArraySectorIntCovar(36), // array of sector alignment covariances 
  //
  // Kalman filter for global alignment
  //
  fSectorParam(0),     // Kalman parameter   
  fSectorCovar(0)     // Kalman covariance  

{
  //
  // Default constructor
  //
  AliTPCROC * roc = AliTPCROC::Instance();
  fXquadrant = roc->GetPadRowRadii(36,53);
  fXmiddle   = ( roc->GetPadRowRadii(0,0)+roc->GetPadRowRadii(36,roc->GetNRows(36)-1))*0.5;
  fXIO       = ( roc->GetPadRowRadii(0,roc->GetNRows(0)-1)+roc->GetPadRowRadii(36,0))*0.5;
}

AliTPCPointCorrection::AliTPCPointCorrection(const Text_t *name, const Text_t *title):
  TNamed(name,title),
  fParamsOutR(38),
  fParamsOutZ(38),
  fParamOutRVersion(0),
  fErrorsOutR(38),
  fErrorsOutZ(38),
  fParamOutZVersion(0),
  //
  // 
  //
  fXIO(0),
  fXmiddle(0),
  fXquadrant(0),
  fArraySectorIntParam(36), // array of sector alignment parameters
  fArraySectorIntCovar(36), // array of sector alignment covariances 
  //
  // Kalman filter for global alignment
  //
  fSectorParam(0),     // Kalman parameter   for A side
  fSectorCovar(0)     // Kalman covariance  for A side 
  
{
  //
  // Non default constructor
  //
  AliTPCROC * roc = AliTPCROC::Instance();
  fXquadrant = roc->GetPadRowRadii(36,53);
  fXmiddle   = ( roc->GetPadRowRadii(0,0)+roc->GetPadRowRadii(36,roc->GetNRows(36)-1))*0.5;
  fXIO       = ( roc->GetPadRowRadii(0,roc->GetNRows(0)-1)+roc->GetPadRowRadii(36,0))*0.5;

}

AliTPCPointCorrection::~AliTPCPointCorrection(){
  //
  //
  //
}


AliTPCPointCorrection* AliTPCPointCorrection::Instance()
{
  //
  // Singleton implementation
  // Returns an instance of this class, it is created if neccessary
  //
  if (fgInstance == 0){
    fgInstance = new AliTPCPointCorrection();
  }
  return fgInstance;
}



Double_t AliTPCPointCorrection::GetDrOut(Bool_t isGlobal, Bool_t type, Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  //  return radial correction
  //
  if (fParamOutRVersion==0) return CorrectionOutR0(isGlobal, type,cx,cy,cz,sector);
  return 0;
}

Double_t      AliTPCPointCorrection::SGetDrOut(Bool_t isGlobal, Bool_t type, Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  // return radial correction - static function
  // 
  return fgInstance->GetDrOut(isGlobal, type,cx,cy,cz,sector);
}




Double_t AliTPCPointCorrection::GetDzOut(Bool_t isGlobal, Bool_t type, Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  //
  //
  if (fParamOutZVersion==0) return CorrectionOutZ0(isGlobal, type,cx,cy,cz,sector);
  return 0;
}


Double_t      AliTPCPointCorrection::SGetDzOut(Bool_t isGlobal, Bool_t type, Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  //
  //
  return fgInstance->GetDzOut(isGlobal, type,cx,cy,cz,sector);
}




Double_t AliTPCPointCorrection::CorrectionOutR0(Bool_t isGlobal, Bool_t type,  Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  // return dR correction - for correction version 0 
  // Parameters:
  // isGlobal   - is the x in global frame
  // type       - kTRUE - use sectors - kFALSE - only Side param
  // cx, cy,cz  - cluster position
  // sector     - sector number
  if (isGlobal){
    // recalculate sector if in global frame
    Double_t alpha    = TMath::ATan2(cy,cx);
    if (alpha<0)  alpha+=TMath::Pi()*2;
    sector = Int_t(18*alpha/(TMath::Pi()*2));
  }

  if (type==kFALSE) sector=36+(sector%36>=18);  // side level parameters
  // dout
  Double_t radius = TMath::Sqrt(cx*cx+cy*cy);  
  AliTPCROC * roc = AliTPCROC::Instance();
  Double_t xout = roc->GetPadRowRadiiUp(roc->GetNRows(37)-1);
  Double_t dout = xout-radius;  
  if (dout<0) return 0;
  //drift
  Double_t dr   = 0.5 - TMath::Abs(cz/250.);
  //
  //
  TVectorD * vec = GetParamOutR(sector);
  if (!vec) return 0;
  Double_t eout10 = TMath::Exp(-dout/10.);
  Double_t eout5 = TMath::Exp(-dout/5.);		    
  Double_t eoutp  = (eout10+eout5)*0.5;
  Double_t eoutm  = (eout10-eout5)*0.5;
  //
  Double_t result=0;
  result+=(*vec)[1]*eoutp;
  result+=(*vec)[2]*eoutm;
  result+=(*vec)[3]*eoutp*dr;
  result+=(*vec)[4]*eoutm*dr;
  result+=(*vec)[5]*eoutp*dr*dr;
  result+=(*vec)[6]*eoutm*dr*dr;
  return result;
}

Double_t AliTPCPointCorrection::RPhiCOGCorrection(Int_t sector, Int_t padrow, Float_t pad, Float_t cy, Float_t y, Float_t z, Float_t ky,Float_t qMax, Float_t threshold){
  //
  // Calculates COG corection in RPHI direction
  // cluster and track position  y is supposed to be corrected before for other effects
  // (e.g ExB and alignemnt)
  // Rphi distortion dependeds on the distance to the edge-pad, distance to the wire edge and
  // relative distance to the center of the pad. Therefore the y position is trnsfromed to the 
  // pad coordiante frame (correction offset (ExB alignemnt) substracted). 
  //   
  // Input parameters:
  //
  // sector - sector number - 0-71  - cl.GetDetector()
  // padrow - padrow number -0-63 -IROC 0-95 OROC - cl->GetRow()
  // pad    - mean pad number  - cl->GetPad()
  // cy     - cluster y        - cl->GetY()
  //  y     - track -or cluster y
  //  qMax  - cluster max charge - cl-.GetMax()
  //  threshold - clusterer threshold
  //
  AliTPCClusterParam * clparam = AliTPCcalibDB::Instance()->GetClusterParam();
  if (!clparam) {
    AliFatal("TPC OCDB not initialized"); 
    return 0;
  }
  Int_t padtype=0;
  if (sector>=36) padtype = (padrow>64)?2:1;
  Double_t padwidth=(padtype==0)? 0.4:0.6;

  Double_t sigma = clparam->GetRMS0(0,padtype,250-TMath::Abs(z),ky)/padwidth;
  //
  Int_t   npads   =  AliTPCROC::Instance()->GetNPads(sector,padrow);
  Float_t yshift  =  TMath::Abs(cy)-TMath::Abs((-npads*0.5+pad)*padwidth);   // the clusters can be corrected before
                                            // shift to undo correction
  
  y -= yshift*((y>0)?1.:-1.);                             // y in the sector coordinate
  Double_t y0     = npads*0.5*padwidth;
  Double_t dy     = (TMath::Abs(y0)-TMath::Abs(y))/padwidth-0.5;  // distance to  the center of pad0   
  //
  Double_t padcenter = TMath::Nint(dy);
  Double_t sumw=0;
  Double_t sumwi=0;
  for (Double_t ip=padcenter-2.5; ip<=padcenter+2.5;ip++){
    Double_t weightGaus = qMax*TMath::Exp(-(dy-ip)*(dy-ip)/(2*sigma*sigma));
    Double_t ypad       = (ip-npads*0.5)*padwidth;
    Double_t weightGain = GetEdgeQ0(sector,padrow,ypad);
    //
    Double_t weight = TMath::Nint(weightGaus*weightGain);
    if (ip<0 &&weight> threshold) weight = threshold;  // this is done in cl finder
    if (weight<0) continue;
    sumw+=weight;
    sumwi+=weight*(ip);    
  }
  Double_t result =0;
  if (sumw>0) result = sumwi/sumw;
  result = (result-dy)*padwidth;
  return result;
}




Double_t AliTPCPointCorrection::CorrectionOutZ0(Bool_t isGlobal, Bool_t type,  Double_t cx, Double_t cy, Double_t cz, Int_t sector){
  //
  // return dR correction - for correction version 0 
  // Parameters:
  // isGlobal   - is the x in global frame
  // type       - kTRUE - use sectors - kFALSE - only Side param
  // cx, cy,cz  - cluster position
  // sector     - sector number
  if (isGlobal){
    // recalculate sector if in global frame
    Double_t alpha    = TMath::ATan2(cy,cx);
    if (alpha<0)  alpha+=TMath::Pi()*2;
    sector = Int_t(18*alpha/(TMath::Pi()*2));
  }

  if (type==kFALSE) sector=36+(sector%36>=18);  // side level parameters
  // dout
  Double_t radius = TMath::Sqrt(cx*cx+cy*cy);  
  AliTPCROC * roc = AliTPCROC::Instance();
  Double_t xout = roc->GetPadRowRadiiUp(roc->GetNRows(37)-1);
  Double_t dout = xout-radius;
  if (dout<0) return 0;
  //drift
  Double_t dr   = 0.5 - TMath::Abs(cz/250.);
  //
  //
  TVectorD * vec = GetParamOutR(sector);
  if (!vec) return 0;
  Double_t eout10 = TMath::Exp(-dout/10.);
  Double_t eout5 = TMath::Exp(-dout/5.);		    
  Double_t eoutp  = (eout10+eout5)*0.5;
  Double_t eoutm  = (eout10-eout5)*0.5;
  //
  Double_t result=0;
  result+=(*vec)[1]*eoutp;
  result+=(*vec)[2]*eoutm;
  result+=(*vec)[3]*eoutp*dr;
  result+=(*vec)[4]*eoutm*dr;
  result+=(*vec)[5]*eoutp*dr*dr;
  result+=(*vec)[6]*eoutm*dr*dr;
  return result;

}

Double_t  AliTPCPointCorrection::GetEdgeQ0(Int_t sector, Int_t padrow, Float_t y){
  //
  /* TF1 fexp("fexp","1-exp(-[0]*(x-[1]))",0,20)
          | param [0] | param [1]
     IROC | 4.71413   | 1.39558
     OROC1| 2.11437   | 1.52643
     OROC2| 2.15082   | 1.53537 
  */

  Double_t params[2]={100,0};
  if (sector<36){
    params[0]=4.71413; params[1]=1.39558;
  }
  if (sector>=36){
    if (padrow<64) {  params[0]=2.114; params[1]=1.526;}
    if (padrow>=64){  params[0]=2.15; params[1]=1.535;}
  }
  Double_t result= 1;
  Double_t xrow  = AliTPCROC::Instance()->GetPadRowRadii(sector,padrow);
  Double_t ymax  = TMath::Tan(TMath::Pi()/18.)*xrow;
  Double_t dedge = ymax-TMath::Abs(y);
  if (dedge-params[1]<0)             return 0;
  if (dedge>10.*params[1]) return 1;
  result = 1.-TMath::Exp(-params[0]*(dedge-params[1]));
  return result;
}

Double_t AliTPCPointCorrection::SRPhiCOGCorrection(Int_t sector, Int_t padrow, Float_t pad, Float_t cy, Float_t y, Float_t z, Float_t ky,Float_t qMax, Float_t threshold){
  return fgInstance->RPhiCOGCorrection(sector, padrow, pad, cy, y,z, ky, qMax, threshold);
}

Double_t AliTPCPointCorrection::SGetEdgeQ0(Int_t sector, Int_t padrow, Float_t y){
  //
  //
  return fgInstance->GetEdgeQ0(sector, padrow, y);
}

void     AliTPCPointCorrection::AddCorrectionSector(TObjArray & sideAPar, TObjArray &sideCPar, TObjArray & sideACov, TObjArray &sideCCov, Bool_t reset){
  //
  //
  //
  for (Int_t isec=0;isec<36;isec++){
    TMatrixD * mat1 = (TMatrixD*)(sideAPar.At(isec));
    TMatrixD * mat1Covar = (TMatrixD*)(sideACov.At(isec));
    if (!mat1) continue;
    TMatrixD * mat0 = (TMatrixD*)(fArraySectorIntParam.At(isec));
    TMatrixD * mat0Covar = (TMatrixD*)(fArraySectorIntCovar.At(isec));
    if (!mat0) {
      fArraySectorIntParam.AddAt(mat1->Clone(),isec); 
      fArraySectorIntCovar.AddAt(mat1Covar->Clone(),isec); 
      continue;
    }
    (*mat0Covar)=(*mat1Covar);      
    if (reset)   (*mat0)=(*mat1);
    if (!reset) (*mat0)+=(*mat1);
  }

  for (Int_t isec=0;isec<36;isec++){
    TMatrixD * mat1 = (TMatrixD*)(sideCPar.At(isec));
    TMatrixD * mat1Covar = (TMatrixD*)(sideCCov.At(isec));
    if (!mat1) continue;
    TMatrixD * mat0 = (TMatrixD*)(fArraySectorIntParam.At(isec));
    TMatrixD * mat0Covar = (TMatrixD*)(fArraySectorIntCovar.At(isec));
    if (!mat0) {
      fArraySectorIntParam.AddAt(mat1->Clone(),isec); 
      fArraySectorIntCovar.AddAt(mat1Covar->Clone(),isec); 
      continue;
    }
    (*mat0Covar)=(*mat1Covar);      
    if (reset)   (*mat0)=(*mat1);
    if (!reset) (*mat0)+=(*mat1);
  }
}


Double_t AliTPCPointCorrection::GetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t /*lz*/, Int_t quadrant){
  //
  // Get position correction for given sector
  //

  TMatrixD * param = (TMatrixD*)fArraySectorIntParam.At(sector%36);
  if (!param) return 0;
  if (quadrant<0){   //recalc quadrant if not specified
    if (lx<fXIO) quadrant=0;
    if(lx>fXIO){
      if (lx<fXquadrant) {
	if (ly<0) quadrant=1;
	if (ly>0) quadrant=2;
      }
      if (lx>=fXquadrant) {
	if (ly<0) quadrant=3;
	if (ly>0) quadrant=4;
      }
    }    
  }
  Double_t a10 = (*param)(quadrant*6+0,0);
  Double_t a20 = (*param)(quadrant*6+1,0);
  Double_t a21 = (*param)(quadrant*6+2,0);
  Double_t dx  = (*param)(quadrant*6+3,0);
  Double_t dy  = (*param)(quadrant*6+4,0);
  Double_t dz  = (*param)(quadrant*6+5,0);
  Double_t deltaX = lx-fXmiddle;
  if (coord==0) return dx;
  if (coord==1) return dy+deltaX*a10;
  if (coord==2) return dz+deltaX*a20+ly*a21;
  if (coord==3) return a10;
  if (coord==4) return a20;
  if (coord==5) return a21;
  //
  return 0;
}

Double_t AliTPCPointCorrection::SGetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz, Int_t quadrant){
  //
  //
  //
  if (!Instance()) return 0;
  return Instance()->GetCorrectionSector(coord,sector,lx,ly,lz, quadrant);
}



Double_t AliTPCPointCorrection::GetCorrection(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t /*lz*/){
  //
  // Get position correction for given sector
  //
  if (!fSectorParam) return 0;
  
  Double_t a10 = (*fSectorParam)(sector*6+0,0);
  Double_t a20 = (*fSectorParam)(sector*6+1,0);
  Double_t a21 = (*fSectorParam)(sector*6+2,0);
  Double_t dx  = (*fSectorParam)(sector*6+3,0);
  Double_t dy  = (*fSectorParam)(sector*6+4,0);
  Double_t dz  = (*fSectorParam)(sector*6+5,0);
  Double_t deltaX = lx-fXIO;
  //
  if (coord==0) return dx;
  if (coord==1) return dy+deltaX*a10;
  if (coord==2) return dz+deltaX*a20+ly*a21;
  if (coord==3) return a10;
  if (coord==4) return a20;
  if (coord==5) return a21;
  return 0;
}

Double_t AliTPCPointCorrection::SGetCorrection(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz){
  //
  //
  //
  if (!Instance()) return 0;
  return Instance()->GetCorrection(coord,sector,lx,ly,lz);
}







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