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> AliTPCBoundaryVoltError class   </h2>       
//   This class calculates the space point distortions due to residual voltage errors on 
//   the main boundaries of the TPC. For example, the inner vessel of the TPC is shifted 
//   by a certain amount, whereas the ROCs on the A side and the C side follow this mechanical 
//   shift (at the inner vessel) in z direction. This example can be named "conical deformation"
//   of the TPC field cage (see example below).                    
//   <p>
//   The boundary conditions can be set via two arrays (A and C side) which contain the 
//   residual voltage setting modeling a possible shift or an inhomogeneity on the TPC field 
//   cage. In order to avoid that the user splits the Central Electrode (CE), the settings for 
//   the C side is taken from the array on the A side (points: A.6 and A.7). The region betweem
//    the points is interpolated linearly onto the boundaries.
//   <p>
//   The class uses the PoissonRelaxation2D (see AliTPCCorrection) to calculate the resulting 
//   electrical field inhomogeneities in the (r,z)-plane. Then, the Langevin-integral formalism 
//   is used to calculate the space point distortions. <br>
//   Note: This class assumes a homogeneous magnetic field. 
//   <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. This emulates for example a Gating-Grid Voltage offset 
//   without moving the ROCs. When the flag is set to TRUE, the ROCs are assumed to be misaligned 
//   and the corresponding offset in z is added.
// End_Html
//
// Begin_Macro(source)
//   {
//   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
//   TCanvas *c2 = new TCanvas("cAliTPCBoundaryVoltError","cAliTPCBoundaryVoltError",500,300); 
//   AliTPCBoundaryVoltError bve;     
//   Float_t val = 40;// [Volt]; 40V corresponds to 1mm
//   /* IFC shift, CE follows, ROC follows by factor half */
//   Float_t boundA[8] = { val, val, val,0,0,0,0,val}; // voltages A-side
//   Float_t boundC[6] = {-val,-val,-val,0,0,0};       // voltages C-side  
//   bve.SetBoundariesA(boundA);                                        
//   bve.SetBoundariesC(boundC);                                        
//   bve.SetOmegaTauT1T2(-0.32,1,1);
//   bve.SetROCDisplacement(kTRUE); // include the chamber offset in z when calculating the dz distortions
//   bve.CreateHistoDRinZR(0)->Draw("surf2"); 
//   return c2;
//   } 
// End_Macro
//
// Begin_Html
//   <p>
// Date: 01/06/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 "TMath.h"
#include "AliTPCROC.h"
#include "AliTPCBoundaryVoltError.h"

ClassImp(AliTPCBoundaryVoltError)

AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
  : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
    fC0(0.),fC1(0.),
    fROCdisplacement(kTRUE),
    fInitLookUp(kFALSE)
{
  //
  // default constructor
  //
  for (Int_t i=0; i<8; i++){
    fBoundariesA[i]= 0;  
    if (i<6) fBoundariesC[i]= 0;
  }
}

AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
  //
  // default destructor
  //
}




Bool_t AliTPCBoundaryVoltError::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;
  }  
  AliTPCBoundaryVoltError* corrC = dynamic_cast<AliTPCBoundaryVoltError *>(corr);
  if (corrC == NULL) {
    AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
    return kFALSE;
  }
  if (fROCdisplacement!=corrC->fROCdisplacement){
    AliError(TString::Format("Inconsistent fROCdisplacement : %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
    return kFALSE;    
  }
  for (Int_t i=0;i <8; i++){
    fBoundariesA[i]+= corrC->fBoundariesA[i]*weight;
    fBoundariesC[i]+= corrC->fBoundariesC[i]*weight;
  }
  //
  return kTRUE;
}




void AliTPCBoundaryVoltError::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);

  InitBoundaryVoltErrorDistortion();
}

void AliTPCBoundaryVoltError::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 AliTPCBoundaryVoltError::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
  //   

  if (!fInitLookUp) {
    AliInfo("Lookup table was not initialized!  Perform the inizialisation now ...");
    InitBoundaryVoltErrorDistortion();
  }

  Int_t   order     = 1 ;               // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
                                        // note that the poisson solution was linearly mirroed on this grid!
  Double_t intEr, intEphi, intdEz ;
  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
  

  intEphi = 0.0;  // Efield is symmetric in phi - 2D calculation

  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 E field integral
  Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
  // Get DeltaEz field integral
  Interpolate2DEdistortion( order, r, z, fLookUpDeltaEz, intdEz );
  
  // 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] = intdEz;  // z distortion - (internally scaled with driftvelocity dependency 
                   // on the Ez field plus the actual ROC misalignment (if set TRUE)


}

void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
  //
  // Initialization of the Lookup table which contains the solutions of the 
  // Dirichlet boundary problem
  //

  const Float_t  gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
  const Float_t  gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;

  TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
  TMatrixD chargeDensity(kRows,kColumns);                              // dummy charge
  TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
  TMatrixD arrayDeltaEzA(kRows,kColumns),  arrayDeltaEzC(kRows,kColumns); // solution

  Double_t  rList[kRows], zedList[kColumns] ;
  
  // Fill arrays with initial conditions.  V on the boundary and ChargeDensity in the volume.      
  for ( Int_t j = 0 ; j < kColumns ; j++ ) {
    Double_t zed = j*gridSizeZ ;
    zedList[j] = zed ;
    for ( Int_t i = 0 ; i < kRows ; i++ )  {
      Double_t radius = fgkIFCRadius + i*gridSizeR ;
      rList[i]           = radius ;
      voltArrayA(i,j)        = 0;  // Initialize voltArrayA to zero
      voltArrayC(i,j)        = 0;  // Initialize voltArrayC to zero
      chargeDensity(i,j)     = 0;  // Initialize ChargeDensity to zero - not used in this class
    }
  }      


  // check if boundary values are the same for the C side (for later, saving some calculation time)

  Int_t symmetry = -1; // assume that  A and C side are identical (but anti-symmetric!) // e.g conical deformation
  Int_t sVec[8];

  // check if boundaries are different (regardless the sign)
  for (Int_t i=0; i<8; i++) { 
    if (TMath::Abs(TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5) 
      symmetry = 0;  
    sVec[i] = (Int_t)( TMath::Sign((Float_t)1.,fBoundariesA[i]) * TMath::Sign((Float_t)1.,fBoundariesC[i])); 
  }
  if (symmetry==-1) { // still the same values?
    // check the kind of symmetry , if even ...
    if (sVec[0]==1 && sVec[1]==1 && sVec[2]==1 && sVec[3]==1 && sVec[4]==1 && sVec[5]==1 && sVec[6]==1 && sVec[7]==1 ) 
      symmetry =  1;
    else if (sVec[0]==-1 && sVec[1]==-1 && sVec[2]==-1 && sVec[3]==-1 && sVec[4]==-1 && sVec[5]==-1 && sVec[6]==-1 && sVec[7]==-1 ) 
      symmetry = -1;
    else
      symmetry =  0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
  }



  // Solve the electrosatic problem in 2D 

  // Fill the complete Boundary vectors
  // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
  for ( Int_t j = 0 ; j < kColumns ; j++ ) {
    Double_t zed = j*gridSizeZ ;
    for ( Int_t i = 0 ; i < kRows ; i++ ) { 
      Double_t radius = fgkIFCRadius + i*gridSizeR ;

      // A side boundary vectors
      if ( i == 0 ) voltArrayA(i,j) += zed   *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
	+ fBoundariesA[0] ; // IFC
      if ( j == kColumns-1 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR))
	+ fBoundariesA[2] ; // ROC
      if ( i == kRows-1 ) voltArrayA(i,j) += zed   *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
	+ fBoundariesA[5] ; // OFC
      if ( j == 0 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR))
	+ fBoundariesA[7] ; // CE
      
      if (symmetry==0) {
	// C side boundary vectors
	if ( i == 0 ) voltArrayC(i,j) += zed   *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
	  + fBoundariesC[0] ; // IFC
	if ( j == kColumns-1 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR))
	  + fBoundariesC[2] ; // ROC
	if ( i == kRows-1 ) voltArrayC(i,j) += zed   *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
	  + fBoundariesC[5] ; // OFC
	if ( j == 0 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR))
	  + fBoundariesC[7] ; // CE
      }

    }
  }

  voltArrayA(0,0)               *= 0.5 ; // Use average boundary condition at corner
  voltArrayA(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
  voltArrayA(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
  voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner

  if (symmetry==0) {
    voltArrayC(0,0)               *= 0.5 ; // Use average boundary condition at corner
    voltArrayC(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
    voltArrayC(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
    voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
  }


  // always solve the problem on the A side
  PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, arrayDeltaEzA, 
		       kRows, kColumns, kIterations, fROCdisplacement ) ;

  if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
    for ( Int_t j = 0 ; j < kColumns ; j++ ) {
      for ( Int_t i = 0 ; i < kRows ; i++ ) { 
	arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
	arrayDeltaEzC(i,j) = -symmetry*arrayDeltaEzA(i,j);
      }
    }
  } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
    PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, arrayDeltaEzC,
			 kRows, kColumns, kIterations, fROCdisplacement ) ;
    for ( Int_t j = 0 ; j < kColumns ; j++ ) {
      for ( Int_t i = 0 ; i < kRows ; i++ ) { 
	arrayDeltaEzC(i,j) = -arrayDeltaEzC(i,j); // negative z coordinate!
      }
    }
  }

  // Interpolate results onto standard grid for Electric Fields
  Int_t ilow=0, jlow=0 ;
  Double_t z,r;
  Float_t saveEr[2] ;	      
  Float_t saveEz[2] ;	      
  for ( Int_t i = 0 ; i < kNZ ; ++i )  {
    z = TMath::Abs( fgkZList[i] ) ;
    for ( Int_t j = 0 ; j < kNR ; ++j ) {
      // Linear interpolation !!
      r = fgkRList[j] ;
      Search( kRows,   rList, r, ilow ) ;          // Note switch - R in rows and Z in columns
      Search( kColumns, zedList, z, jlow ) ;
      if ( ilow < 0 ) ilow = 0 ;                   // check if out of range
      if ( jlow < 0 ) jlow = 0 ;   
      if ( ilow + 1  >=  kRows - 1 ) ilow =  kRows - 2 ;	      
      if ( jlow + 1  >=  kColumns - 1 ) jlow =  kColumns - 2 ; 
      
      if (fgkZList[i]>0) {         // A side solution
	saveEr[0] = arrayErOverEzA(ilow,jlow) + 
	  (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEr[1] = arrayErOverEzA(ilow+1,jlow) + 
	  (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEz[0] = arrayDeltaEzA(ilow,jlow) + 
	  (arrayDeltaEzA(ilow,jlow+1)-arrayDeltaEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEz[1] = arrayDeltaEzA(ilow+1,jlow) + 
	  (arrayDeltaEzA(ilow+1,jlow+1)-arrayDeltaEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;

      } else if (fgkZList[i]<0) {  // C side solution
	saveEr[0] = arrayErOverEzC(ilow,jlow) + 
	  (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEr[1] = arrayErOverEzC(ilow+1,jlow) + 
	  (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEz[0] = arrayDeltaEzC(ilow,jlow) + 
	  (arrayDeltaEzC(ilow,jlow+1)-arrayDeltaEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
	saveEz[1] = arrayDeltaEzC(ilow+1,jlow) + 
	  (arrayDeltaEzC(ilow+1,jlow+1)-arrayDeltaEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;

      } else {
	AliWarning("Field calculation at z=0 (CE) is not allowed!");
	saveEr[0]=0; saveEr[1]=0;
	saveEz[0]=0; saveEz[1]=0;
      }
      fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
      fLookUpDeltaEz[i][j]  = saveEz[0] + (saveEz[1]-saveEz[0])*(r-rList[ilow])/gridSizeR ;
    }
  }
  
  voltArrayA.Clear();
  voltArrayC.Clear();
  chargeDensity.Clear();
  arrayErOverEzA.Clear();
  arrayErOverEzC.Clear();
  arrayDeltaEzA.Clear();
  arrayDeltaEzC.Clear();
  
  fInitLookUp = kTRUE;

}

void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
  //
  // Print function to check the settings of the boundary vectors
  // option=="a" prints the C0 and C1 coefficents for calibration purposes
  //

  TString opt = option; opt.ToLower();
  printf("%s\n",GetTitle());
  printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
  printf("  : A-side (anti-clockwise)\n"); 
  printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
  printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
  printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
  printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
  printf("  : C-side (clockwise)\n"); 
  printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
  printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
  printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
  printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);

  // Check wether the settings of the Central Electrode agree (on the A and C side)
  // Note: they have to be antisymmetric!
  if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
    AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
    AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
    AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
  }

  if (opt.Contains("a")) { // Print all details
    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 InitBoundaryVoltErrorDistortion() ...");

}


void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
  //
  // set voltage errors on the TPC boundaries - A side 
  //
  // Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through 
  // IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear 
  // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.  
  // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
  // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
  // correspond to ~ 40 V
  // 
  // Note: The setting for the CE will be passed to the C side!
  
  for (Int_t i=0; i<8; i++) {
    fBoundariesA[i]= boundariesA[i];  
    if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
  }
  fInitLookUp=kFALSE;
}
void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
  //
  // set voltage errors on the TPC boundaries - C side 
  //
  // Start at IFC at the Central electrode and work clockwise (for C side) through 
  // IFC, ROC and OFC. The boundary conditions are currently defined to be a linear 
  // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.  
  // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
  // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
  // correspond to ~ 40 V
  // 
  // Note: The setting for the CE will be taken from the A side (pos 6 and 7)!

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