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

// _________________________________________________________________
//
// Begin_Html
//   <h2> AliTPCROCVoltError3D class   </h2>    
//   The class calculates the space point distortions due to z offsets of the 
//   ROCs via the residual voltage technique (Poisson relaxation) in 3D. 
//   Since the GG (part of the ROCs) represents the closure of the FC in z direction,
//   every misalignment in z produces not only dz distortions but also electrical 
//   field inhomogeneities throughout the volume, which produces additional dr and rd$\phi$ distortions.
//   <p>
//   Each ROC can be misaligned (in z direction) in three ways. A general z0 offset, 
//   an inclination along the x and an inclination along the y axis. The z-misalignment's
//   can be set via the function SetROCData(TMatrixD *mat) for each single chamber (ROC). 
//   The array size has to be (72,3) representing the 72 chambers according to the 
//   offline numbering scheme (IROC: roc$<$36; OROC: roc$\geq$36) and the three misalignment's
//   (see the source code for further details).
//   <p>
//   Internally, these z offsets (unit is cm)  are recalculated into residual voltage 
//   equivalents in order to make use of the relaxation technique. 
//   <p>
//   One has two possibilities when calculating the $dz$ distortions. The resulting 
//   distortions are purely due to the change of the drift velocity (along with the 
//   change of the drift field) when the SetROCDisplacement is FALSE. <br>
//   For this class, this is a rather unphysical setting and should be avoided. <br>
//   When the flag is set to TRUE, the corresponding offset in z is added to the dz 
//   calculation of the outer ROCs. <br>
//   For the Alice TPC gas, both effects are of similar magnitude. This means, if the 
//   drift length is sufficiently large, a z-offset of a chamber appears to have (approx.) 
//   twice the magnitude when one looks at the actual dz distortions.
//   <p>
//   In addition, this class allows a correction regarding the different arrival times 
//   of the electrons due to the geometrical difference of the inner and outer chambers.
//   The impact was simulated via Garfield. If the flag is set via the 
//   function SetElectronArrivalCorrection, the electron-arrival correction is added to the dz calculation.
// End_Html
//
// Begin_Macro(source)
//   {
//   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
//   TCanvas *c2 = new TCanvas("cAliTPCROCVoltError3D","cAliTPCROCVoltError3D",500,400); 
//   AliTPCROCVoltError3D roc; 
//   roc.SetROCDataFileName("$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root");
//   roc.SetElectronArrivalCorrection(kFALSE);  // Correction for electron arrival offset, IROC vs OROC
//   roc.SetROCDisplacement(kTRUE);   // include the chamber offset in z when calculating the dz 
//   roc.SetOmegaTauT1T2(0,1,1); // B=0
//   roc.CreateHistoDZinXY(1.,300,300)->Draw("colz"); 
//   return c2;
//   } 
// End_Macro
//
// Begin_Html
//   <p>
//   Date: 08/08/2010    <br>                                                 
//   Authors: Jim Thomas, Stefan Rossegger                                
// End_Html 
// _________________________________________________________________


#include "AliMagF.h"
#include "TGeoGlobalMagField.h"
#include "AliTPCcalibDB.h"
#include "AliTPCParam.h"
#include "AliLog.h"
#include "TMatrixD.h"
#include "TFile.h"

#include "TMath.h"
#include "AliTPCROC.h"
#include "AliTPCROCVoltError3D.h"

ClassImp(AliTPCROCVoltError3D)

AliTPCROCVoltError3D::AliTPCROCVoltError3D()
  : AliTPCCorrection("ROCVoltErrors","ROC z alignment Errors"),
    fC0(0.),fC1(0.),
    fROCdisplacement(kTRUE),
    fElectronArrivalCorrection(kTRUE),
    fInitLookUp(kFALSE),
    fROCDataFileName(""),  
    fdzDataLinFit(0)
{
  //
  // default constructor
  //

  // Array which will contain the solution according to the setted boundary conditions
  // main input: z alignment of the Read Out chambers
  // see InitROCVoltError3D() function
  for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
    fLookUpErOverEz[k]   =  new TMatrixF(kNR,kNZ);  
    fLookUpEphiOverEz[k] =  new TMatrixF(kNR,kNZ);
    fLookUpDeltaEz[k]    =  new TMatrixF(kNR,kNZ);   
  }
  //  fROCDataFileName="$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root";
  //  SetROCDataFileName(fROCDataFileName.Data()); // initialization of fdzDataLinFit is included

}

AliTPCROCVoltError3D::~AliTPCROCVoltError3D() {
  //
  // destructor
  //
  
  for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
    delete fLookUpErOverEz[k];
    delete fLookUpEphiOverEz[k];
    delete fLookUpDeltaEz[k];
  }

  delete fdzDataLinFit;
}



Bool_t AliTPCROCVoltError3D::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
  //
  // Add correction  and make them compact
  // Assumptions:
  //  - origin of distortion/correction are additive
  //  - only correction ot the same type supported ()
  if (corr==NULL) {
    AliError("Zerro pointer - correction");
    return kFALSE;
  }  
  AliTPCROCVoltError3D * corrC = dynamic_cast<AliTPCROCVoltError3D *>(corr);
  if (corrC == NULL) {
    AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
    return kFALSE;
  }
  //
  TMatrixD matrixDz=*(corrC->fdzDataLinFit);
  matrixDz*=weight;
  if (fdzDataLinFit==NULL) {
    fdzDataLinFit=new TMatrixD(matrixDz);
  }
  else{
    (*fdzDataLinFit)+=matrixDz;
  }
  return kTRUE;
}





void AliTPCROCVoltError3D::SetROCData(TMatrixD * matrix){
  //
  // Set a z alignment map of the chambers not via a file, but
  // directly via a TMatrix(72,3), where dz = p0 + p1*(lx-133.4) + p2*ly (all in cm)
  //
  if (!fdzDataLinFit) fdzDataLinFit=new TMatrixD(*matrix);
  else *fdzDataLinFit = *matrix;
}


void AliTPCROCVoltError3D::Init() {
  //
  // Initialization funtion
  //
  
  AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
  if (!magF) AliError("Magneticd field - not initialized");
  Double_t bzField = magF->SolenoidField()/10.; //field in T
  AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
  if (!param) AliError("Parameters - not initialized");
  Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
  Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
  Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
  // Correction Terms for effective omegaTau; obtained by a laser calibration run
  SetOmegaTauT1T2(wt,fT1,fT2);

  if (!fInitLookUp) InitROCVoltError3D();
}

void AliTPCROCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
  //
  // Update function 
  //
  AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
  if (!magF) AliError("Magneticd field - not initialized");
  Double_t bzField = magF->SolenoidField()/10.; //field in T
  AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
  if (!param) AliError("Parameters - not initialized");
  Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]  // From dataBase: to be updated: per second (ideally)
  Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
  Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
  // Correction Terms for effective omegaTau; obtained by a laser calibration run
  SetOmegaTauT1T2(wt,fT1,fT2);

}

void  AliTPCROCVoltError3D::SetROCDataFileName(const char * fname) {
  //
  // Set / load the ROC data (linear fit of ROC misalignments)
  //

  fROCDataFileName = fname;
  
  TFile f(fROCDataFileName.Data(),"READ");
  TMatrixD *m = (TMatrixD*) f.Get("dzSurveyLinFitData");
  TMatrixD &mf = *m;

  // prepare some space

  if (fdzDataLinFit) delete fdzDataLinFit;
  fdzDataLinFit = new TMatrixD(72,3);
  TMatrixD &dataIntern = *fdzDataLinFit;
  
  for (Int_t iroc=0;iroc<72;iroc++) {
    dataIntern(iroc,0) = mf(iroc,0);  // z0 offset
    dataIntern(iroc,1) = mf(iroc,1);  // slope in x
    dataIntern(iroc,2) = mf(iroc,2);  // slope in y
  }

  f.Close();

  fInitLookUp = kFALSE;

}

void AliTPCROCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
  //
  // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
  //   
  const Double_t kEpsilon=Double_t(FLT_MIN);
  if (!fInitLookUp) {
    AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
    InitROCVoltError3D();
  }
  static Bool_t forceInit=kTRUE; //temporary needed for back compatibility with old OCDB entries
  if (forceInit&&fLookUpErOverEz[0]){
    if (TMath::Abs(fLookUpErOverEz[0]->Sum())<kEpsilon){//temporary needed for back compatibility with old OCDB entries
      ForceInitROCVoltError3D();
    }
    forceInit=kFALSE;
  }

  
  Int_t   order     = 1 ;    // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         

  Float_t intEr, intEphi, intDeltaEz;
  Double_t r, phi, z ;
  Int_t    sign;

  r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
  phi    =  TMath::ATan2(x[1],x[0]) ;
  if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
  z      =  x[2] ;                                         // Create temporary copy of x[2]

  if ( (roc%36) < 18 ) {
    sign =  1;       // (TPC A side)
  } else {
    sign = -1;       // (TPC C side)
  }
  
  if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
  if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
  

  if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
    AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");

  // Get the Er and Ephi field integrals plus the integral over DeltaEz 
  intEr      = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
				  fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz  );
  intEphi    = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
				  fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz);
  intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
				  fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz   );

  //  printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);

  // Calculate distorted position
  if ( r > 0.0 ) {
    phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
    r   =  r   + ( fC0*intEr   + fC1*intEphi );  
  }
  
  // Calculate correction in cartesian coordinates
  dx[0] = r * TMath::Cos(phi) - x[0];
  dx[1] = r * TMath::Sin(phi) - x[1]; 
  dx[2] = intDeltaEz;  // z distortion - (internally scaled with driftvelocity dependency 
                       // on the Ez field plus the actual ROC misalignment (if set TRUE)


  if (fElectronArrivalCorrection) {

    // correction for the OROC (in average, a 0.014usec longer drift time
    // due to different position of the anode wires) -> vd*dt -> 2.64*0.014 = 0.0369 cm
    // FIXME: correction are token from Magboltz simulations
    //        should be loaded from a database
 
    AliTPCROC * rocInfo = AliTPCROC::Instance();
    Double_t rCrossingROC  =  (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
  
    if (r>rCrossingROC) {
      if (sign==1)
	dx[2] += 0.0369; // A side - negative correction
      else
	dx[2] -= 0.0369; // C side - positive correction
    }
    
  }

}

void AliTPCROCVoltError3D::InitROCVoltError3D() {
  //
  // Initialization of the Lookup table which contains the solutions of the 
  // Dirichlet boundary problem
  // Calculation of the single 3D-Poisson solver is done just if needed
  // (see basic lookup tables in header file)
  //

  const Int_t   order       =    1  ;  // Linear interpolation = 1, Quadratic = 2  
  const Float_t gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
  const Float_t gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
  const Float_t gridSizePhi =  TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);

  // temporary arrays to create the boundary conditions
  TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ; 
  TMatrixD *arrayofEroverEz[kPhiSlices], *arrayofEphioverEz[kPhiSlices], *arrayofDeltaEz[kPhiSlices] ; 

  for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
    arrayofArrayV[k]     =   new TMatrixD(kRows,kColumns) ;
    arrayofCharge[k]     =   new TMatrixD(kRows,kColumns) ;
    arrayofEroverEz[k]   =   new TMatrixD(kRows,kColumns) ;
    arrayofEphioverEz[k] =   new TMatrixD(kRows,kColumns) ;
    arrayofDeltaEz[k]    =   new TMatrixD(kRows,kColumns) ;
  }
  
  // list of point as used in the poisson relation and the interpolation (during sum up)
  Double_t  rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
  for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
    philist[k] =  gridSizePhi * k;
    for ( Int_t i = 0 ; i < kRows ; i++ )    {
      rlist[i] = fgkIFCRadius + i*gridSizeR ;
      for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
	zedlist[j]  = j * gridSizeZ ;
      }
    }
  }

  // ==========================================================================
  // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
  // Allow for different size grid spacing in R and Z directions
  
  const Int_t   symmetry = 0;
 
  // Set bondaries and solve Poisson's equation --------------------------
  
  if ( !fInitLookUp ) {
    
    AliInfo(Form("Solving the poisson equation (~ %d sec)",2*10*(int)(kPhiSlices/10)));
    
    for ( Int_t side = 0 ; side < 2 ; side++ ) {  // Solve Poisson3D twice; once for +Z and once for -Z
      
      for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
	TMatrixD &arrayV    =  *arrayofArrayV[k] ;
	TMatrixD &charge    =  *arrayofCharge[k] ;
	
	//Fill arrays with initial conditions.  V on the boundary and Charge in the volume.
	for ( Int_t i = 0 ; i < kRows ; i++ ) {
	  for ( Int_t j = 0 ; j < kColumns ; j++ ) {  // Fill Vmatrix with Boundary Conditions
	    arrayV(i,j) = 0.0 ; 
	    charge(i,j) = 0.0 ;

	    Float_t radius0 = rlist[i] ;
	    Float_t phi0    = gridSizePhi * k ;
	    
	    // To avoid problems at sector boundaries, use an average of +- 1 degree from actual phi location
	    if ( j == (kColumns-1) ) {
	      arrayV(i,j) = 0.5*  ( GetROCVoltOffset( side, radius0, phi0+0.02 ) + GetROCVoltOffset( side, radius0, phi0-0.02 ) ) ;

	      if (side==1) // C side
		arrayV(i,j) = -arrayV(i,j); // minus sign on the C side to allow a consistent usage of global z when setting the boundaries
	    }
	  }
	}      
	
	for ( Int_t i = 1 ; i < kRows-1 ; i++ ) { 
	  for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
	    charge(i,j)  =  0.0 ;
	  }
	}
      }      
      
      // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
      // Allow for different size grid spacing in R and Z directions
      
      PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
			   arrayofEroverEz, arrayofEphioverEz, arrayofDeltaEz,
			   kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, 
			   symmetry, fROCdisplacement) ;
      
      
      //Interpolate results onto a custom grid which is used just for these calculations.
      Double_t  r, phi, z ;
      for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
	phi = fgkPhiList[k] ;
	
	TMatrixF &erOverEz   =  *fLookUpErOverEz[k]  ;
	TMatrixF &ephiOverEz =  *fLookUpEphiOverEz[k];
	TMatrixF &deltaEz    =  *fLookUpDeltaEz[k]   ;
	
	for ( Int_t j = 0 ; j < kNZ ; j++ ) {

	  z = TMath::Abs(fgkZList[j]) ;  // Symmetric solution in Z that depends only on ABS(Z)
  
	  if ( side == 0 &&  fgkZList[j] < 0 ) continue; // Skip rest of this loop if on the wrong side
	  if ( side == 1 &&  fgkZList[j] > 0 ) continue; // Skip rest of this loop if on the wrong side
	  
	  for ( Int_t i = 0 ; i < kNR ; i++ ) { 
	    r = fgkRList[i] ;

	    // Interpolate basicLookup tables; once for each rod, then sum the results
	    erOverEz(i,j)   = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices, 
						 rlist, zedlist, philist, arrayofEroverEz  );
	    ephiOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
						 rlist, zedlist, philist, arrayofEphioverEz);
	    deltaEz(i,j)    = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
						 rlist, zedlist, philist, arrayofDeltaEz  );

	    if (side == 1)  deltaEz(i,j) = -  deltaEz(i,j); // negative coordinate system on C side

	  } // end r loop
	}// end z loop
      }// end phi loop

      if ( side == 0 ) AliInfo(" A side done");
      if ( side == 1 ) AliInfo(" C side done");
    } // end side loop
  }
  
  // clear the temporary arrays lists
  for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
    delete arrayofArrayV[k];
    delete arrayofCharge[k];
    delete arrayofEroverEz[k];  
    delete arrayofEphioverEz[k];
    delete arrayofDeltaEz[k];
  }
 

  fInitLookUp = kTRUE;

}


Float_t AliTPCROCVoltError3D::GetROCVoltOffset(Int_t side, Float_t r0, Float_t phi0) {
  // 
  // Returns the dz alignment data (in voltage equivalents) at 
  // the given position
  //

  Float_t xp = r0*TMath::Cos(phi0);
  Float_t yp = r0*TMath::Sin(phi0);
  
  // phi0 should be between 0 and 2pi 
  if (phi0<0) phi0+=TMath::TwoPi();
  Int_t roc = (Int_t)TMath::Floor((TMath::RadToDeg()*phi0)/20);
  if (side==1) roc+=18; // C side
  if (r0>132) roc+=36;  // OROC 
  
  // linear-plane data:  z = z0 + kx*lx + ky*ly (rotation in local coordinates)
  TMatrixD &fitData = *fdzDataLinFit;

  // local coordinates                                                          
  Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)roc)%18));
  Float_t lx = xp*TMath::Cos(secAlpha)+yp*TMath::Sin(secAlpha);
  Float_t ly = yp*TMath::Cos(secAlpha)-xp*TMath::Sin(secAlpha);

  // reference of rotation in lx is at the intersection between OROC and IROC
  // necessary, since the Fitprozedure is otherwise useless
  
  AliTPCROC * rocInfo = AliTPCROC::Instance();
  Double_t lxRef  = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
  
  Float_t dz = fitData(roc,0)+fitData(roc,1)*(lx-lxRef) + fitData(roc,2)*ly; // value in cm

  // aproximated Voltage-offset-aquivalent to the z misalignment
  // (linearly scaled with the z position)
  Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
  Float_t voltOff = dz*ezField;            // values in "Volt equivalents"

  return voltOff;
}

TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment(Int_t side, Int_t nx, Int_t ny) {
  //
  // return a simple histogramm containing the input to the poisson solver
  // (z positions of the Read-out chambers, linearly interpolated)

  char hname[100];
  if (side==0) snprintf(hname,100,"survey_dz_Aside");
  if (side==1) snprintf(hname,100,"survey_dz_Cside");

  TH2F *h = new TH2F(hname,hname,nx,-250.,250.,ny,-250.,250.);

  for (Int_t iy=1;iy<=ny;++iy) {
    Double_t yp = h->GetYaxis()->GetBinCenter(iy);
    for (Int_t ix=1;ix<=nx;++ix) {
      Double_t xp = h->GetXaxis()->GetBinCenter(ix);
    
      Float_t r0 = TMath::Sqrt(xp*xp+yp*yp);
      Float_t phi0 = TMath::ATan2(yp,xp); 
   
      Float_t dz = GetROCVoltOffset(side,r0,phi0); // in [volt]

      Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
      dz = dz/ezField;    // in [cm]

      if (85.<=r0 && r0<=245.) {
	h->SetBinContent(ix,iy,dz); 
      } else {
	h->SetBinContent(ix,iy,0.);
      }
    }
  }
  
  h->GetXaxis()->SetTitle("x [cm]");
  h->GetYaxis()->SetTitle("y [cm]");
  h->GetZaxis()->SetTitle("dz [cm]");
  h->SetStats(0);
  //  h->DrawCopy("colz");

  return h;
} 

void AliTPCROCVoltError3D::Print(const Option_t* option) const {
  //
  // Print function to check the settings of the Rod shifts and the rotated clips
  // option=="a" prints the C0 and C1 coefficents for calibration purposes
  //

  TString opt = option; opt.ToLower();
  printf("%s\n",GetTitle());
  printf(" - z aligmnet of the TPC Read-Out chambers \n");
  printf("   (linear interpolation within the chamber:  dz = z0 + kx*(lx-133) + ky*ly [cm] ) \n");
  printf("   Info: Check the following data-file for more details: %s \n",fROCDataFileName.Data());

  if (opt.Contains("a")) { // Print all details
    TMatrixD &fitData = *fdzDataLinFit;
    printf(" A side:  IROC   ROCX=(z0,kx,ky): \n");
    for (Int_t roc = 0; roc<18; roc++) 
      printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
    printf("\n A side:  OROC   ROCX=(z0,kx,ky): \n");
    for (Int_t roc = 36; roc<54; roc++) 
      printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
    printf("\n C side:  IROC   ROCX=(z0,kx,ky): \n");
    for (Int_t roc = 18; roc<36; roc++) 
      printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
    printf("\n C side:  OROC   ROCX=(z0,kx,ky): \n");
    for (Int_t roc = 54; roc<72; roc++) 
      printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
    printf("\n\n");
    printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
    printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
  }

  if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitROCVoltError3D() ...");

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