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

#include <TString.h>
#include <TSystem.h>
#include <TROOT.h>
#include <TRandom.h>
#include <stdio.h>
#include <TMethodCall.h>
#include <TMath.h>
#include <TH1.h>
#include "AliCheb3D.h"
#include "AliCheb3DCalc.h"
#include "AliLog.h"

ClassImp(AliCheb3D)

const Float_t AliCheb3D::fgkMinPrec = 1.e-12f;

//__________________________________________________________________________________________
AliCheb3D::AliCheb3D() : 
  fDimOut(0), 
  fPrec(0), 
  fChebCalc(1), 
  fMaxCoefs(0), 
  fResTmp(0), 
  fGrid(0), 
  fUsrFunName(""), 
  fUsrMacro(0) 
{
// Default constructor
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = fArgsTmp[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0;
  }
}

//__________________________________________________________________________________________
AliCheb3D::AliCheb3D(const AliCheb3D& src) : 
  TNamed(src),
  fDimOut(src.fDimOut), 
  fPrec(src.fPrec), 
  fChebCalc(1), 
  fMaxCoefs(src.fMaxCoefs), 
  fResTmp(0),
  fGrid(0), 
  fUsrFunName(src.fUsrFunName), 
  fUsrMacro(0)
{
  // read coefs from text file
  for (int i=3;i--;) {
    fBMin[i]    = src.fBMin[i];
    fBMax[i]    = src.fBMax[i];
    fBScale[i]  = src.fBScale[i];
    fBOffset[i] = src.fBOffset[i];
    fNPoints[i] = src.fNPoints[i];
    fGridOffs[i] = src.fGridOffs[i];
    fArgsTmp[i]  = 0;
  }
  for (int i=0;i<fDimOut;i++) {
    AliCheb3DCalc* cbc = src.GetChebCalc(i);
    if (cbc) fChebCalc.AddAtAndExpand(new AliCheb3DCalc(*cbc),i);
  }
}

//__________________________________________________________________________________________
AliCheb3D::AliCheb3D(const char* inpFile) : 
  fDimOut(0), 
  fPrec(0),  
  fChebCalc(1),
  fMaxCoefs(0),  
  fResTmp(0),
  fGrid(0), 
  fUsrFunName(""), 
  fUsrMacro(0)
{
  // read coefs from text file
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0;
    fArgsTmp[i]  = 0;
  }
  LoadData(inpFile);
}

//__________________________________________________________________________________________
AliCheb3D::AliCheb3D(FILE* stream) : 
  fDimOut(0), 
  fPrec(0), 
  fChebCalc(1), 
  fMaxCoefs(0),
  fResTmp(0),
  fGrid(0),
  fUsrFunName(""),
  fUsrMacro(0)
{
  // read coefs from stream
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0;
    fArgsTmp[i]  = 0;
  }
  LoadData(stream);
}

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
AliCheb3D::AliCheb3D(const char* funName, int DimOut, const Float_t  *bmin, const Float_t  *bmax, Int_t *npoints, Float_t prec, const Float_t* precD) : 
  TNamed(funName,funName), 
  fDimOut(0), 
  fPrec(TMath::Max(fgkMinPrec,prec)), 
  fChebCalc(1), 
  fMaxCoefs(0), 
  fResTmp(0), 
  fGrid(0), 
  fUsrFunName("") ,
  fUsrMacro(0)
{
  // Construct the parameterization for the function
  // funName : name of the file containing the function: void funName(Float_t * inp,Float_t * out)
  // DimOut  : dimension of the vector computed by the user function
  // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
  // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
  // npoints : array of 3 elements with the number of points to compute in each of 3 dimension
  // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
  // precD   : optional array with precisions per output dimension (if >fgkMinPrec will override common prec)
  if (DimOut<1) AliFatalF("Requested output dimension is %d",fDimOut);
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0.;
    fArgsTmp[i]  = 0;
  }
  SetDimOut(DimOut,precD);
  PrepareBoundaries(bmin,bmax);
  DefineGrid(npoints);
  SetUsrFunction(funName);
  ChebFit();
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Int_t *npoints, Float_t prec, const Float_t* precD) : 
  TNamed("Cheb3D","Cheb3D"),
  fDimOut(0), 
  fPrec(TMath::Max(fgkMinPrec,prec)), 
  fChebCalc(1), 
  fMaxCoefs(0), 
  fResTmp(0), 
  fGrid(0), 
  fUsrFunName(""),
  fUsrMacro(0)
{
  // Construct the parameterization for the function
  // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
  // DimOut  : dimension of the vector computed by the user function
  // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
  // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
  // npoints : array of 3 elements with the number of points to compute in each of 3 dimension
  // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
  // precD   : optional array with precisions per output dimension (if >fgkMinPrec will override common prec)
  //
  if (DimOut<1) AliFatalF("Requested output dimension is %d",fDimOut);
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0.;
    fArgsTmp[i]  = 0;
  }
  SetDimOut(DimOut,precD);
  PrepareBoundaries(bmin,bmax);
  DefineGrid(npoints);
  SetUsrFunction(ptr);
  ChebFit();
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Int_t *npX,Int_t *npY,Int_t *npZ, Float_t prec, const Float_t* precD) : 
  TNamed("Cheb3D","Cheb3D"),
  fDimOut(0), 
  fPrec(TMath::Max(fgkMinPrec,prec)), 
  fChebCalc(1), 
  fMaxCoefs(0), 
  fResTmp(0), 
  fGrid(0), 
  fUsrFunName(""),
  fUsrMacro(0)
{
  // Construct very economic  parameterization for the function
  // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
  // DimOut  : dimension of the vector computed by the user function
  // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
  // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
  // npX     : array of 3 elements with the number of points to compute in each dimension for 1st component 
  // npY     : array of 3 elements with the number of points to compute in each dimension for 2nd component 
  // npZ     : array of 3 elements with the number of points to compute in each dimension for 3d  component 
  // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
  // precD   : optional array with precisions per output dimension (if >fgkMinPrec will override common prec)
  //
  if (DimOut<1) AliFatalF("Requested output dimension is %d",fDimOut);
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0.;
    fArgsTmp[i]  = 0;
  }
  SetDimOut(DimOut,precD);
  PrepareBoundaries(bmin,bmax);
  SetUsrFunction(ptr);
  //
  DefineGrid(npX);
  ChebFit(0);
  DefineGrid(npY);
  ChebFit(1);
  DefineGrid(npZ);
  ChebFit(2);
  //
}
#endif


//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Float_t prec, Bool_t run, const Float_t* precD) : 
  TNamed("Cheb3D","Cheb3D"),
  fDimOut(0), 
  fPrec(TMath::Max(fgkMinPrec,prec)), 
  fChebCalc(1), 
  fMaxCoefs(0), 
  fResTmp(0), 
  fGrid(0), 
  fUsrFunName(""),
  fUsrMacro(0)
{
  // Construct very economic  parameterization for the function with automatic calculation of the root's grid
  // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
  // DimOut  : dimension of the vector computed by the user function
  // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
  // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
  // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
  // precD   : optional array with precisions per output dimension (if >fgkMinPrec will override common prec)
  //
  if (DimOut!=3) AliFatalF("This constructor works only for 3D fits, %dD fit was requested",fDimOut);
  if (DimOut<1)  AliFatalF("Requested output dimension is %d",fDimOut);
  for (int i=3;i--;) {
    fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
    fNPoints[i] = 0;
    fGridOffs[i] = 0.;
    fArgsTmp[i]  = 0;
  }
  SetDimOut(DimOut,precD);
  PrepareBoundaries(bmin,bmax);
  SetUsrFunction(ptr);
  //
  if (run) {
    int gridNC[3][3];
    EstimateNPoints(prec,gridNC);
    DefineGrid(gridNC[0]);
    ChebFit(0);
    DefineGrid(gridNC[1]);
    ChebFit(1);
    DefineGrid(gridNC[2]);
    ChebFit(2);
  }
  //
}
#endif


//__________________________________________________________________________________________
AliCheb3D& AliCheb3D::operator=(const AliCheb3D& rhs)
{
  // assignment operator
  //
  if (this != &rhs) {
    Clear();
    fDimOut   = rhs.fDimOut;
    fPrec     = rhs.fPrec;
    fMaxCoefs = rhs.fMaxCoefs;
    fUsrFunName = rhs.fUsrFunName;
    fUsrMacro   = 0;
    for (int i=3;i--;) {
      fBMin[i]    = rhs.fBMin[i];
      fBMax[i]    = rhs.fBMax[i];
      fBScale[i]  = rhs.fBScale[i];
      fBOffset[i] = rhs.fBOffset[i];
      fNPoints[i] = rhs.fNPoints[i];
    } 
    for (int i=0;i<fDimOut;i++) {
      AliCheb3DCalc* cbc = rhs.GetChebCalc(i);
      if (cbc) fChebCalc.AddAtAndExpand(new AliCheb3DCalc(*cbc),i);
    }    
  }
  return *this;
  //
}

//__________________________________________________________________________________________
void AliCheb3D::Clear(const Option_t*)
{
  // clear all dynamic structures
  //
  if (fResTmp)        { delete[] fResTmp; fResTmp = 0; }
  if (fGrid)          { delete[] fGrid;   fGrid   = 0; }
  if (fUsrMacro)      { delete fUsrMacro; fUsrMacro = 0;}
  fChebCalc.SetOwner(kTRUE);
  fChebCalc.Delete();
  //
}

//__________________________________________________________________________________________
void AliCheb3D::Print(const Option_t* opt) const
{
  // print info
  //
  printf("%s: Chebyshev parameterization for 3D->%dD function. Precision: %e\n",GetName(),fDimOut,fPrec);
  printf("Region of validity: [%+.5e:%+.5e] [%+.5e:%+.5e] [%+.5e:%+.5e]\n",fBMin[0],fBMax[0],fBMin[1],fBMax[1],fBMin[2],fBMax[2]);
  TString opts = opt; opts.ToLower();
  if (opts.Contains("l")) for (int i=0;i<fDimOut;i++) {printf("Output dimension %d:\n",i+1); GetChebCalc(i)->Print();}
  //
}

//__________________________________________________________________________________________
void AliCheb3D::PrepareBoundaries(const Float_t  *bmin, const Float_t  *bmax)
{
  // Set and check boundaries defined by user, prepare coefficients for their conversion to [-1:1] interval
  //
  for (int i=3;i--;) {
    fBMin[i]   = bmin[i];
    fBMax[i]   = bmax[i];
    fBScale[i] = bmax[i]-bmin[i];
    if (fBScale[i]<=0) { 
      AliFatalF("Boundaries for %d-th dimension are not increasing: %+.4e %+.4e\nStop\n",i,fBMin[i],fBMax[i]);
    }
    fBOffset[i] = bmin[i] + fBScale[i]/2.0;
    fBScale[i] = 2./fBScale[i];
  }
  //
}


//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_

// Pointer on user function (faster altrnative to TMethodCall)
void (*gUsrFunAliCheb3D) (float* ,float* );

void AliCheb3D::EvalUsrFunction() 
{
  // call user supplied function
  if   (gUsrFunAliCheb3D) gUsrFunAliCheb3D(fArgsTmp,fResTmp);
  else fUsrMacro->Execute(); 
}

void AliCheb3D::SetUsrFunction(const char* name)
{
  // load user macro with function definition and compile it
  gUsrFunAliCheb3D = 0; 
  fUsrFunName = name;
  gSystem->ExpandPathName(fUsrFunName);
  if (fUsrMacro) delete fUsrMacro;
  TString tmpst = fUsrFunName;
  tmpst += "+"; // prepare filename to compile
  if (gROOT->LoadMacro(tmpst.Data())) AliFatalF("Failed to load user function from %s",name);
  fUsrMacro = new TMethodCall();        
  tmpst = tmpst.Data() + tmpst.Last('/')+1; //Strip away any path preceding the macro file name
  int dot = tmpst.Last('.');
  if (dot>0) tmpst.Resize(dot);
  fUsrMacro->InitWithPrototype(tmpst.Data(),"Float_t *,Float_t *");
  long args[2];
  args[0] = (long)fArgsTmp;
  args[1] = (long)fResTmp;
  fUsrMacro->SetParamPtrs(args); 
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::SetUsrFunction(void (*ptr)(float*,float*))
{
  // assign user training function
  //
  if (fUsrMacro) delete fUsrMacro;
  fUsrMacro = 0;
  fUsrFunName = "";
  gUsrFunAliCheb3D = ptr;
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::EvalUsrFunction(const Float_t  *x, Float_t  *res) 
{
  // evaluate user function value
  //
  for (int i=3;i--;) fArgsTmp[i] = x[i];
  if   (gUsrFunAliCheb3D) gUsrFunAliCheb3D(fArgsTmp,fResTmp);
  else fUsrMacro->Execute(); 
  for (int i=fDimOut;i--;) res[i] = fResTmp[i];
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
Int_t AliCheb3D::CalcChebCoefs(const Float_t  *funval,int np, Float_t  *outCoefs, Float_t  prec)
{
  // Calculate Chebyshev coeffs using precomputed function values at np roots.
  // If prec>0, estimate the highest coeff number providing the needed precision
  //
  double sm;                 // do summations in double to minimize the roundoff error
  for (int ic=0;ic<np;ic++) { // compute coeffs
    sm = 0;          
    for (int ir=0;ir<np;ir++) {
      float  rt = TMath::Cos( ic*(ir+0.5)*TMath::Pi()/np);
      sm += funval[ir]*rt;
    }
    outCoefs[ic] = Float_t( sm * ((ic==0) ? 1./np : 2./np) );
  }
  //
  if (prec<=0) return np;
  //
  sm = 0;
  int cfMax = 0;
  for (cfMax=np;cfMax--;) {
    sm += TMath::Abs(outCoefs[cfMax]);
    if (sm>=prec) break;
  }
  if (++cfMax==0) cfMax=1;
  return cfMax;
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::DefineGrid(Int_t* npoints)
{
  // prepare the grid of Chebyshev roots in each dimension
  const int kMinPoints = 1;
  int ntot = 0;
  fMaxCoefs = 1;
  for (int id=3;id--;) { 
    fNPoints[id] = npoints[id];
    if (fNPoints[id]<kMinPoints) AliFatalF("at %d-th dimension %d point is requested, at least %d is needed",id,fNPoints[id],kMinPoints);
    ntot += fNPoints[id];
    fMaxCoefs *= fNPoints[id];
  }
  printf("Computing Chebyshev nodes on [%2d/%2d/%2d] grid\n",npoints[0],npoints[1],npoints[2]);
  if (fGrid) delete[] fGrid;
  fGrid = new Float_t [ntot];
  //
  int curp = 0;
  for (int id=3;id--;) { 
    int np = fNPoints[id];
    fGridOffs[id] = curp;
    for (int ip=0;ip<np;ip++) {
      Float_t x = TMath::Cos( TMath::Pi()*(ip+0.5)/np );
      fGrid[curp++] = MapToExternal(x,id);
    }
  }
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
Int_t AliCheb3D::ChebFit()
{
  // prepare parameterization for all output dimensions
  int ir=0; 
  for (int i=fDimOut;i--;) ir+=ChebFit(i); 
  return ir;
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
Int_t AliCheb3D::ChebFit(int dmOut)
{
  // prepare paramaterization of 3D function for dmOut-th dimension 
  int maxDim = 0;
  for (int i=0;i<3;i++) if (maxDim<fNPoints[i]) maxDim = fNPoints[i];
  Float_t  *fvals      = new Float_t [ fNPoints[0] ];
  Float_t  *tmpCoef3D  = new Float_t [ fNPoints[0]*fNPoints[1]*fNPoints[2] ]; 
  Float_t  *tmpCoef2D  = new Float_t [ fNPoints[0]*fNPoints[1] ]; 
  Float_t  *tmpCoef1D  = new Float_t [ maxDim ];
  //
  // 1D Cheb.fit for 0-th dimension at current steps of remaining dimensions
  int ncmax = 0;
  //
  printf("Dim%d : 00.00%% Done",dmOut);fflush(stdout);
  AliCheb3DCalc* cheb =  GetChebCalc(dmOut);
  //
  Float_t prec = cheb->GetPrecision(); 
  if (prec<fgkMinPrec) prec = fPrec;         // no specific precision for this dim.
  //
  Float_t rTiny = 0.1*prec/Float_t(maxDim); // neglect coefficient below this threshold
  //
  float ncals2count = fNPoints[2]*fNPoints[1]*fNPoints[0];
  float ncals = 0;
  float frac = 0;
  float fracStep = 0.001;
  //
  for (int id2=fNPoints[2];id2--;) {
    fArgsTmp[2] = fGrid[ fGridOffs[2]+id2 ];
    //
    for (int id1=fNPoints[1];id1--;) {
      fArgsTmp[1] = fGrid[ fGridOffs[1]+id1 ];
      //
      for (int id0=fNPoints[0];id0--;) {
	fArgsTmp[0] = fGrid[ fGridOffs[0]+id0 ];
	EvalUsrFunction();     // compute function values at Chebyshev roots of 0-th dimension
	fvals[id0] =  fResTmp[dmOut];
	float fr = (++ncals)/ncals2count;
	if (fr-frac>=fracStep) {
	  frac = fr;
	  printf("\b\b\b\b\b\b\b\b\b\b\b");
	  printf("%05.2f%% Done",fr*100);
	  fflush(stdout);
	}
	//
      }
      int nc = CalcChebCoefs(fvals,fNPoints[0], tmpCoef1D, prec);
      for (int id0=fNPoints[0];id0--;) tmpCoef2D[id1 + id0*fNPoints[1]] = tmpCoef1D[id0];
      if (ncmax<nc) ncmax = nc;              // max coefs to be kept in dim0 to guarantee needed precision
    }
    //
    // once each 1d slice of given 2d slice is parametrized, parametrize the Cheb.coeffs
    for (int id0=fNPoints[0];id0--;) {
      CalcChebCoefs( tmpCoef2D+id0*fNPoints[1], fNPoints[1], tmpCoef1D, -1);
      for (int id1=fNPoints[1];id1--;) tmpCoef3D[id2 + fNPoints[2]*(id1+id0*fNPoints[1])] = tmpCoef1D[id1];
    }
  }
  //
  // now fit the last dimensions Cheb.coefs
  for (int id0=fNPoints[0];id0--;) {
    for (int id1=fNPoints[1];id1--;) {
      CalcChebCoefs( tmpCoef3D+ fNPoints[2]*(id1+id0*fNPoints[1]), fNPoints[2], tmpCoef1D, -1);
      for (int id2=fNPoints[2];id2--;) tmpCoef3D[id2+ fNPoints[2]*(id1+id0*fNPoints[1])] = tmpCoef1D[id2]; // store on place
    }
  }
  //
  // now find 2D surface which separates significant coefficients of 3D matrix from nonsignificant ones (up to prec)
  UShort_t *tmpCoefSurf = new UShort_t[ fNPoints[0]*fNPoints[1] ];
  for (int id0=fNPoints[0];id0--;) for (int id1=fNPoints[1];id1--;) tmpCoefSurf[id1+id0*fNPoints[1]]=0;  
  Double_t resid = 0;
  for (int id0=fNPoints[0];id0--;) {
    for (int id1=fNPoints[1];id1--;) {
      for (int id2=fNPoints[2];id2--;) {
	int id = id2 + fNPoints[2]*(id1+id0*fNPoints[1]);
	Float_t  cfa = TMath::Abs(tmpCoef3D[id]);
	if (cfa < rTiny) {tmpCoef3D[id] = 0; continue;} // neglect coefs below the threshold
	resid += cfa;
	if (resid<prec) continue; // this coeff is negligible
	// otherwise go back 1 step
	resid -= cfa;
	tmpCoefSurf[id1+id0*fNPoints[1]] = id2+1; // how many coefs to keep
	break;
      }
    }
  }
  /*
  printf("\n\nCoeffs\n");  
  int cnt = 0;
  for (int id0=0;id0<fNPoints[0];id0++) {
    for (int id1=0;id1<fNPoints[1];id1++) {
      for (int id2=0;id2<fNPoints[2];id2++) {
	printf("%2d%2d%2d %+.4e |",id0,id1,id2,tmpCoef3D[cnt++]);
      }
      printf("\n");
    }
    printf("\n");
  }
  */
  // see if there are rows to reject, find max.significant column at each row
  int nRows = fNPoints[0];
  UShort_t *tmpCols = new UShort_t[nRows]; 
  for (int id0=fNPoints[0];id0--;) {
    int id1 = fNPoints[1];
    while (id1>0 && tmpCoefSurf[(id1-1)+id0*fNPoints[1]]==0) id1--;
    tmpCols[id0] = id1;
  }
  // find max significant row
  for (int id0=nRows;id0--;) {if (tmpCols[id0]>0) break; nRows--;}
  // find max significant column and fill the permanent storage for the max sigificant column of each row
  cheb->InitRows(nRows);                  // create needed arrays;
  UShort_t *nColsAtRow = cheb->GetNColsAtRow();
  UShort_t *colAtRowBg = cheb->GetColAtRowBg();
  int nCols = 0;
  int nElemBound2D = 0;
  for (int id0=0;id0<nRows;id0++) {
    nColsAtRow[id0] = tmpCols[id0];     // number of columns to store for this row
    colAtRowBg[id0] = nElemBound2D;     // begining of this row in 2D boundary surface
    nElemBound2D += tmpCols[id0];
    if (nCols<nColsAtRow[id0]) nCols = nColsAtRow[id0];
  }
  cheb->InitCols(nCols);
  delete[] tmpCols;
  //  
  // create the 2D matrix defining the boundary of significance for 3D coeffs.matrix 
  // and count the number of siginifacnt coefficients
  //
  cheb->InitElemBound2D(nElemBound2D);
  UShort_t *coefBound2D0 = cheb->GetCoefBound2D0();
  UShort_t *coefBound2D1 = cheb->GetCoefBound2D1();
  fMaxCoefs = 0; // redefine number of coeffs
  for (int id0=0;id0<nRows;id0++) {
    int nCLoc = nColsAtRow[id0];
    int col0  = colAtRowBg[id0];
    for (int id1=0;id1<nCLoc;id1++) {
      coefBound2D0[col0 + id1] = tmpCoefSurf[id1+id0*fNPoints[1]];  // number of coefs to store for 3-d dimension
      coefBound2D1[col0 + id1] = fMaxCoefs;
      fMaxCoefs += coefBound2D0[col0 + id1];
    }
  }
  //
  // create final compressed 3D matrix for significant coeffs
  cheb->InitCoefs(fMaxCoefs);
  Float_t  *coefs = cheb->GetCoefs();
  int count = 0;
  for (int id0=0;id0<nRows;id0++) {
    int ncLoc = nColsAtRow[id0];
    int col0  = colAtRowBg[id0];
    for (int id1=0;id1<ncLoc;id1++) {
      int ncf2 = coefBound2D0[col0 + id1];
      for (int id2=0;id2<ncf2;id2++) {
	coefs[count++] = tmpCoef3D[id2 + fNPoints[2]*(id1+id0*fNPoints[1])];
      }
    }
  }
  /*
  printf("\n\nNewSurf\n");
  for (int id0=0;id0<fNPoints[0];id0++) {
    for (int id1=0;id1<fNPoints[1];id1++) {
      printf("(%2d %2d) %2d |",id0,id1,tmpCoefSurf[id1+id0*fNPoints[1]]);  
    }
    printf("\n");
  }
  */
  //
  delete[] tmpCoefSurf;
  delete[] tmpCoef1D;
  delete[] tmpCoef2D;
  delete[] tmpCoef3D;
  delete[] fvals;
  //
  printf("\b\b\b\b\b\b\b\b\b\b\b\b");
  printf("100.00%% Done\n");
  return 1;
}
#endif

//_______________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::SaveData(const char* outfile,Bool_t append) const
{
  // writes coefficients data to output text file, optionallt appending on the end of existing file
  TString strf = outfile;
  gSystem->ExpandPathName(strf);
  FILE* stream = fopen(strf,append ? "a":"w");
  SaveData(stream);
  fclose(stream);
  //
}
#endif

//_______________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::SaveData(FILE* stream) const
{
  // writes coefficients data to existing output stream
  //
  fprintf(stream,"\n# These are automatically generated data for the Chebyshev interpolation of 3D->%dD function\n",fDimOut); 
  fprintf(stream,"#\nSTART %s\n",GetName());
  fprintf(stream,"# Dimensionality of the output\n%d\n",fDimOut);
  fprintf(stream,"# Interpolation abs. precision\n%+.8e\n",fPrec);
  //
  fprintf(stream,"# Lower boundaries of interpolation region\n");
  for (int i=0;i<3;i++) fprintf(stream,"%+.8e\n",fBMin[i]);
  fprintf(stream,"# Upper boundaries of interpolation region\n");
  for (int i=0;i<3;i++) fprintf(stream,"%+.8e\n",fBMax[i]);
  fprintf(stream,"# Parameterization for each output dimension follows:\n");
  //
  for (int i=0;i<fDimOut;i++) GetChebCalc(i)->SaveData(stream);
  fprintf(stream,"#\nEND %s\n#\n",GetName());
  //
}
#endif

//__________________________________________________________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
void AliCheb3D::InvertSign()
{
  // invert the sign of all parameterizations
  for (int i=fDimOut;i--;) {
    AliCheb3DCalc* par =  GetChebCalc(i);
    int ncf = par->GetNCoefs();
    float *coefs = par->GetCoefs();
    for (int j=ncf;j--;) coefs[j] = -coefs[j];
  }
}
#endif


//_______________________________________________
void AliCheb3D::LoadData(const char* inpFile)
{
  // load coefficients data from txt file
  //
  TString strf = inpFile;
  gSystem->ExpandPathName(strf);
  FILE* stream = fopen(strf.Data(),"r");
  LoadData(stream);
  fclose(stream);
  //
}

//_______________________________________________
void AliCheb3D::LoadData(FILE* stream)
{
  // load coefficients data from stream
  //
  if (!stream) AliFatal("No stream provided");
  TString buffs;
  Clear();
  AliCheb3DCalc::ReadLine(buffs,stream);
  if (!buffs.BeginsWith("START")) AliFatalF("Expected: \"START <fit_name>\", found \"%s\"",buffs.Data());
  SetName(buffs.Data()+buffs.First(' ')+1);
  //
  AliCheb3DCalc::ReadLine(buffs,stream); // N output dimensions
  fDimOut = buffs.Atoi(); 
  if (fDimOut<1) AliFatalF("Expected: '<number_of_output_dimensions>', found \"%s\"",buffs.Data());
  //
  SetDimOut(fDimOut);
  //
  AliCheb3DCalc::ReadLine(buffs,stream); // Interpolation abs. precision
  fPrec = buffs.Atof();
  if (fPrec<=0) AliFatalF("Expected: '<abs.precision>', found \"%s\"",buffs.Data());
  //
  for (int i=0;i<3;i++) { // Lower boundaries of interpolation region
    AliCheb3DCalc::ReadLine(buffs,stream);
    fBMin[i] = buffs.Atof(); 
  }
  for (int i=0;i<3;i++) { // Upper boundaries of interpolation region
    AliCheb3DCalc::ReadLine(buffs,stream);
    fBMax[i] = buffs.Atof(); 
  }
  PrepareBoundaries(fBMin,fBMax);
  //
  // data for each output dimension
  for (int i=0;i<fDimOut;i++) GetChebCalc(i)->LoadData(stream);
  //
  // check end_of_data record
  AliCheb3DCalc::ReadLine(buffs,stream);
  if (!buffs.BeginsWith("END") || !buffs.Contains(GetName())) {
    AliFatalF("Expected \"END %s\", found \"%s\"",GetName(),buffs.Data());
  }
  //
}

//_______________________________________________
void AliCheb3D::SetDimOut(const int d, const float* prec)
{
  // init output dimensions
  fDimOut = d;
  if (fResTmp) delete fResTmp;
  fResTmp = new Float_t[fDimOut];
  fChebCalc.Delete();
  for (int i=0;i<d;i++) {
    AliCheb3DCalc* clc = new AliCheb3DCalc();
    clc->SetPrecision(prec && prec[i]>fgkMinPrec ? prec[i] : fPrec);
    fChebCalc.AddAtAndExpand(clc,i);
  }
}

//_______________________________________________
void AliCheb3D::ShiftBound(int id,float dif)
{
  // modify the bounds of the grid
  //
  if (id<0||id>2) {printf("Maximum 3 dimensions are supported\n"); return;}
  fBMin[id] += dif;
  fBMax[id] += dif;
  fBOffset[id] += dif;
}

//_______________________________________________
#ifdef _INC_CREATION_ALICHEB3D_
TH1* AliCheb3D::TestRMS(int idim,int npoints,TH1* histo)
{
  // fills the difference between the original function and parameterization (for idim-th component of the output)
  // to supplied histogram. Calculations are done in npoints random points. 
  // If the hostgram was not supplied, it will be created. It is up to the user to delete it! 
  if (!fUsrMacro) {
    printf("No user function is set\n");
    return 0;
  }
  float prc = GetChebCalc(idim)->GetPrecision();
  if (prc<fgkMinPrec) prc = fPrec;   // no dimension specific precision
  if (!histo) histo = new TH1D(GetName(),"Control: Function - Parametrization",100,-2*prc,2*prc);
  for (int ip=npoints;ip--;) {
    gRandom->RndmArray(3,(Float_t *)fArgsTmp);
    for (int i=3;i--;) fArgsTmp[i] = fBMin[i] + fArgsTmp[i]*(fBMax[i]-fBMin[i]);
    EvalUsrFunction();
    Float_t valFun = fResTmp[idim];
    Eval(fArgsTmp,fResTmp);
    Float_t valPar = fResTmp[idim];
    histo->Fill(valFun - valPar);
  }
  return histo;
  //
}
#endif

//_______________________________________________
#ifdef _INC_CREATION_ALICHEB3D_

void AliCheb3D::EstimateNPoints(float prec, int gridBC[3][3],Int_t npd1,Int_t npd2,Int_t npd3)
{
  // Estimate number of points to generate a training data
  //
  const int kScp = 9;
  const float kScl[9] = {0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9};
  //
  const float sclDim[2] = {0.001,0.999};
  const int   compDim[3][2] = { {1,2}, {2,0}, {0,1} };
  static float xyz[3];
  Int_t npdTst[3] = {npd1,npd2,npd3};
  //

  for (int i=3;i--;)for (int j=3;j--;) gridBC[i][j] = -1;
  //
  for (int idim=0;idim<3;idim++) {
    float dimMN = fBMin[idim] + sclDim[0]*(fBMax[idim]-fBMin[idim]);
    float dimMX = fBMin[idim] + sclDim[1]*(fBMax[idim]-fBMin[idim]);
    //
    int id1 = compDim[idim][0]; // 1st fixed dim
    int id2 = compDim[idim][1]; // 2nd fixed dim
    for (int i1=0;i1<kScp;i1++) {
      xyz[ id1 ] = fBMin[id1] + kScl[i1]*( fBMax[id1]-fBMin[id1] );
      for (int i2=0;i2<kScp;i2++) {
	xyz[ id2 ] = fBMin[id2] + kScl[i2]*( fBMax[id2]-fBMin[id2] );
	int* npt = GetNCNeeded(xyz,idim, dimMN,dimMX, prec, npdTst[idim]); // npoints for Bx,By,Bz
	for (int ib=0;ib<3;ib++) if (npt[ib]>gridBC[ib][idim]) gridBC[ib][idim] = npt[ib];
      }
    }
  }
}

/*
void AliCheb3D::EstimateNPoints(float prec, int gridBC[3][3])
{
  // Estimate number of points to generate a training data
  //
  const float sclA[9] = {0.1, 0.5, 0.9, 0.1, 0.5, 0.9, 0.1, 0.5, 0.9} ;
  const float sclB[9] = {0.1, 0.1, 0.1, 0.5, 0.5, 0.5, 0.9, 0.9, 0.9} ;
  const float sclDim[2] = {0.01,0.99};
  const int   compDim[3][2] = { {1,2}, {2,0}, {0,1} };
  static float xyz[3];
  //
  for (int i=3;i--;)for (int j=3;j--;) gridBC[i][j] = -1;
  //
  for (int idim=0;idim<3;idim++) {
    float dimMN = fBMin[idim] + sclDim[0]*(fBMax[idim]-fBMin[idim]);
    float dimMX = fBMin[idim] + sclDim[1]*(fBMax[idim]-fBMin[idim]);
    //
    for (int it=0;it<9;it++) { // test in 9 points
      int id1 = compDim[idim][0]; // 1st fixed dim
      int id2 = compDim[idim][1]; // 2nd fixed dim
      xyz[ id1 ] = fBMin[id1] + sclA[it]*( fBMax[id1]-fBMin[id1] );
      xyz[ id2 ] = fBMin[id2] + sclB[it]*( fBMax[id2]-fBMin[id2] );
      //
      int* npt = GetNCNeeded(xyz,idim, dimMN,dimMX, prec); // npoints for Bx,By,Bz
      for (int ib=0;ib<3;ib++) if (npt[ib]>gridBC[ib][idim]) gridBC[ib][idim] = npt[ib];//+2;
      //
    }
  }
}


int* AliCheb3D::GetNCNeeded(float xyz[3],int DimVar, float mn,float mx, float prec)
{
  // estimate needed number of chebyshev coefs for given function description in DimVar dimension
  // The values for two other dimensions must be set beforehand
  //
  static int curNC[3];
  static int retNC[3];
  const int kMaxPoint = 400;
  float* gridVal = new float[3*kMaxPoint];
  float* coefs   = new float[3*kMaxPoint];
  //
  float scale = mx-mn;
  float offs  = mn + scale/2.0;
  scale = 2./scale;
  // 
  int curNP;
  int maxNC=-1;
  int maxNCPrev=-1;
  for (int i=0;i<3;i++) retNC[i] = -1;
  for (int i=0;i<3;i++) fArgsTmp[i] = xyz[i];
  //
  for (curNP=3; curNP<kMaxPoint; curNP+=3) { 
    maxNCPrev = maxNC;
    //
    for (int i=0;i<curNP;i++) { // get function values on Cheb. nodes
      float x = TMath::Cos( TMath::Pi()*(i+0.5)/curNP );
      fArgsTmp[DimVar] =  x/scale+offs; // map to requested interval
      EvalUsrFunction();
      for (int ib=3;ib--;) gridVal[ib*kMaxPoint + i] = fResTmp[ib];
    }
    //
    for (int ib=0;ib<3;ib++) {
      curNC[ib] = AliCheb3D::CalcChebCoefs(&gridVal[ib*kMaxPoint], curNP, &coefs[ib*kMaxPoint],prec);
      if (maxNC < curNC[ib]) maxNC = curNC[ib];
      if (retNC[ib] < curNC[ib]) retNC[ib] = curNC[ib];
    }
    if ( (curNP-maxNC)>3 &&  (maxNC-maxNCPrev)<1 ) break;
    maxNCPrev = maxNC;
    //
  }
  delete[] gridVal;
  delete[] coefs;
  return retNC;
  //
}
*/


int* AliCheb3D::GetNCNeeded(float xyz[3],int DimVar, float mn,float mx, float prec, Int_t npCheck)
{
  // estimate needed number of chebyshev coefs for given function description in DimVar dimension
  // The values for two other dimensions must be set beforehand
  //
  static int retNC[3];
  static int npChLast = 0;
  static float *gridVal=0,*coefs=0;
  if (npCheck<3) npCheck = 3;
  if (npChLast<npCheck) {
    if (gridVal) delete[] gridVal;
    if (coefs)   delete[] coefs;
    gridVal = new float[3*npCheck];
    coefs   = new float[3*npCheck];
    npChLast = npCheck;
  }
  //
  float scale = mx-mn;
  float offs  = mn + scale/2.0;
  scale = 2./scale;
  //
  for (int i=0;i<3;i++) fArgsTmp[i] = xyz[i];
  for (int i=0;i<npCheck;i++) {
    fArgsTmp[DimVar] =  TMath::Cos( TMath::Pi()*(i+0.5)/npCheck)/scale+offs; // map to requested interval
    EvalUsrFunction();
    for (int ib=3;ib--;) gridVal[ib*npCheck + i] = fResTmp[ib];
  } 
  //
  for (int ib=0;ib<3;ib++) retNC[ib] = AliCheb3D::CalcChebCoefs(&gridVal[ib*npCheck], npCheck, &coefs[ib*npCheck],prec);
  return retNC;
  //
}


#endif
 AliCheb3D.cxx:1
 AliCheb3D.cxx:2
 AliCheb3D.cxx:3
 AliCheb3D.cxx:4
 AliCheb3D.cxx:5
 AliCheb3D.cxx:6
 AliCheb3D.cxx:7
 AliCheb3D.cxx:8
 AliCheb3D.cxx:9
 AliCheb3D.cxx:10
 AliCheb3D.cxx:11
 AliCheb3D.cxx:12
 AliCheb3D.cxx:13
 AliCheb3D.cxx:14
 AliCheb3D.cxx:15
 AliCheb3D.cxx:16
 AliCheb3D.cxx:17
 AliCheb3D.cxx:18
 AliCheb3D.cxx:19
 AliCheb3D.cxx:20
 AliCheb3D.cxx:21
 AliCheb3D.cxx:22
 AliCheb3D.cxx:23
 AliCheb3D.cxx:24
 AliCheb3D.cxx:25
 AliCheb3D.cxx:26
 AliCheb3D.cxx:27
 AliCheb3D.cxx:28
 AliCheb3D.cxx:29
 AliCheb3D.cxx:30
 AliCheb3D.cxx:31
 AliCheb3D.cxx:32
 AliCheb3D.cxx:33
 AliCheb3D.cxx:34
 AliCheb3D.cxx:35
 AliCheb3D.cxx:36
 AliCheb3D.cxx:37
 AliCheb3D.cxx:38
 AliCheb3D.cxx:39
 AliCheb3D.cxx:40
 AliCheb3D.cxx:41
 AliCheb3D.cxx:42
 AliCheb3D.cxx:43
 AliCheb3D.cxx:44
 AliCheb3D.cxx:45
 AliCheb3D.cxx:46
 AliCheb3D.cxx:47
 AliCheb3D.cxx:48
 AliCheb3D.cxx:49
 AliCheb3D.cxx:50
 AliCheb3D.cxx:51
 AliCheb3D.cxx:52
 AliCheb3D.cxx:53
 AliCheb3D.cxx:54
 AliCheb3D.cxx:55
 AliCheb3D.cxx:56
 AliCheb3D.cxx:57
 AliCheb3D.cxx:58
 AliCheb3D.cxx:59
 AliCheb3D.cxx:60
 AliCheb3D.cxx:61
 AliCheb3D.cxx:62
 AliCheb3D.cxx:63
 AliCheb3D.cxx:64
 AliCheb3D.cxx:65
 AliCheb3D.cxx:66
 AliCheb3D.cxx:67
 AliCheb3D.cxx:68
 AliCheb3D.cxx:69
 AliCheb3D.cxx:70
 AliCheb3D.cxx:71
 AliCheb3D.cxx:72
 AliCheb3D.cxx:73
 AliCheb3D.cxx:74
 AliCheb3D.cxx:75
 AliCheb3D.cxx:76
 AliCheb3D.cxx:77
 AliCheb3D.cxx:78
 AliCheb3D.cxx:79
 AliCheb3D.cxx:80
 AliCheb3D.cxx:81
 AliCheb3D.cxx:82
 AliCheb3D.cxx:83
 AliCheb3D.cxx:84
 AliCheb3D.cxx:85
 AliCheb3D.cxx:86
 AliCheb3D.cxx:87
 AliCheb3D.cxx:88
 AliCheb3D.cxx:89
 AliCheb3D.cxx:90
 AliCheb3D.cxx:91
 AliCheb3D.cxx:92
 AliCheb3D.cxx:93
 AliCheb3D.cxx:94
 AliCheb3D.cxx:95
 AliCheb3D.cxx:96
 AliCheb3D.cxx:97
 AliCheb3D.cxx:98
 AliCheb3D.cxx:99
 AliCheb3D.cxx:100
 AliCheb3D.cxx:101
 AliCheb3D.cxx:102
 AliCheb3D.cxx:103
 AliCheb3D.cxx:104
 AliCheb3D.cxx:105
 AliCheb3D.cxx:106
 AliCheb3D.cxx:107
 AliCheb3D.cxx:108
 AliCheb3D.cxx:109
 AliCheb3D.cxx:110
 AliCheb3D.cxx:111
 AliCheb3D.cxx:112
 AliCheb3D.cxx:113
 AliCheb3D.cxx:114
 AliCheb3D.cxx:115
 AliCheb3D.cxx:116
 AliCheb3D.cxx:117
 AliCheb3D.cxx:118
 AliCheb3D.cxx:119
 AliCheb3D.cxx:120
 AliCheb3D.cxx:121
 AliCheb3D.cxx:122
 AliCheb3D.cxx:123
 AliCheb3D.cxx:124
 AliCheb3D.cxx:125
 AliCheb3D.cxx:126
 AliCheb3D.cxx:127
 AliCheb3D.cxx:128
 AliCheb3D.cxx:129
 AliCheb3D.cxx:130
 AliCheb3D.cxx:131
 AliCheb3D.cxx:132
 AliCheb3D.cxx:133
 AliCheb3D.cxx:134
 AliCheb3D.cxx:135
 AliCheb3D.cxx:136
 AliCheb3D.cxx:137
 AliCheb3D.cxx:138
 AliCheb3D.cxx:139
 AliCheb3D.cxx:140
 AliCheb3D.cxx:141
 AliCheb3D.cxx:142
 AliCheb3D.cxx:143
 AliCheb3D.cxx:144
 AliCheb3D.cxx:145
 AliCheb3D.cxx:146
 AliCheb3D.cxx:147
 AliCheb3D.cxx:148
 AliCheb3D.cxx:149
 AliCheb3D.cxx:150
 AliCheb3D.cxx:151
 AliCheb3D.cxx:152
 AliCheb3D.cxx:153
 AliCheb3D.cxx:154
 AliCheb3D.cxx:155
 AliCheb3D.cxx:156
 AliCheb3D.cxx:157
 AliCheb3D.cxx:158
 AliCheb3D.cxx:159
 AliCheb3D.cxx:160
 AliCheb3D.cxx:161
 AliCheb3D.cxx:162
 AliCheb3D.cxx:163
 AliCheb3D.cxx:164
 AliCheb3D.cxx:165
 AliCheb3D.cxx:166
 AliCheb3D.cxx:167
 AliCheb3D.cxx:168
 AliCheb3D.cxx:169
 AliCheb3D.cxx:170
 AliCheb3D.cxx:171
 AliCheb3D.cxx:172
 AliCheb3D.cxx:173
 AliCheb3D.cxx:174
 AliCheb3D.cxx:175
 AliCheb3D.cxx:176
 AliCheb3D.cxx:177
 AliCheb3D.cxx:178
 AliCheb3D.cxx:179
 AliCheb3D.cxx:180
 AliCheb3D.cxx:181
 AliCheb3D.cxx:182
 AliCheb3D.cxx:183
 AliCheb3D.cxx:184
 AliCheb3D.cxx:185
 AliCheb3D.cxx:186
 AliCheb3D.cxx:187
 AliCheb3D.cxx:188
 AliCheb3D.cxx:189
 AliCheb3D.cxx:190
 AliCheb3D.cxx:191
 AliCheb3D.cxx:192
 AliCheb3D.cxx:193
 AliCheb3D.cxx:194
 AliCheb3D.cxx:195
 AliCheb3D.cxx:196
 AliCheb3D.cxx:197
 AliCheb3D.cxx:198
 AliCheb3D.cxx:199
 AliCheb3D.cxx:200
 AliCheb3D.cxx:201
 AliCheb3D.cxx:202
 AliCheb3D.cxx:203
 AliCheb3D.cxx:204
 AliCheb3D.cxx:205
 AliCheb3D.cxx:206
 AliCheb3D.cxx:207
 AliCheb3D.cxx:208
 AliCheb3D.cxx:209
 AliCheb3D.cxx:210
 AliCheb3D.cxx:211
 AliCheb3D.cxx:212
 AliCheb3D.cxx:213
 AliCheb3D.cxx:214
 AliCheb3D.cxx:215
 AliCheb3D.cxx:216
 AliCheb3D.cxx:217
 AliCheb3D.cxx:218
 AliCheb3D.cxx:219
 AliCheb3D.cxx:220
 AliCheb3D.cxx:221
 AliCheb3D.cxx:222
 AliCheb3D.cxx:223
 AliCheb3D.cxx:224
 AliCheb3D.cxx:225
 AliCheb3D.cxx:226
 AliCheb3D.cxx:227
 AliCheb3D.cxx:228
 AliCheb3D.cxx:229
 AliCheb3D.cxx:230
 AliCheb3D.cxx:231
 AliCheb3D.cxx:232
 AliCheb3D.cxx:233
 AliCheb3D.cxx:234
 AliCheb3D.cxx:235
 AliCheb3D.cxx:236
 AliCheb3D.cxx:237
 AliCheb3D.cxx:238
 AliCheb3D.cxx:239
 AliCheb3D.cxx:240
 AliCheb3D.cxx:241
 AliCheb3D.cxx:242
 AliCheb3D.cxx:243
 AliCheb3D.cxx:244
 AliCheb3D.cxx:245
 AliCheb3D.cxx:246
 AliCheb3D.cxx:247
 AliCheb3D.cxx:248
 AliCheb3D.cxx:249
 AliCheb3D.cxx:250
 AliCheb3D.cxx:251
 AliCheb3D.cxx:252
 AliCheb3D.cxx:253
 AliCheb3D.cxx:254
 AliCheb3D.cxx:255
 AliCheb3D.cxx:256
 AliCheb3D.cxx:257
 AliCheb3D.cxx:258
 AliCheb3D.cxx:259
 AliCheb3D.cxx:260
 AliCheb3D.cxx:261
 AliCheb3D.cxx:262
 AliCheb3D.cxx:263
 AliCheb3D.cxx:264
 AliCheb3D.cxx:265
 AliCheb3D.cxx:266
 AliCheb3D.cxx:267
 AliCheb3D.cxx:268
 AliCheb3D.cxx:269
 AliCheb3D.cxx:270
 AliCheb3D.cxx:271
 AliCheb3D.cxx:272
 AliCheb3D.cxx:273
 AliCheb3D.cxx:274
 AliCheb3D.cxx:275
 AliCheb3D.cxx:276
 AliCheb3D.cxx:277
 AliCheb3D.cxx:278
 AliCheb3D.cxx:279
 AliCheb3D.cxx:280
 AliCheb3D.cxx:281
 AliCheb3D.cxx:282
 AliCheb3D.cxx:283
 AliCheb3D.cxx:284
 AliCheb3D.cxx:285
 AliCheb3D.cxx:286
 AliCheb3D.cxx:287
 AliCheb3D.cxx:288
 AliCheb3D.cxx:289
 AliCheb3D.cxx:290
 AliCheb3D.cxx:291
 AliCheb3D.cxx:292
 AliCheb3D.cxx:293
 AliCheb3D.cxx:294
 AliCheb3D.cxx:295
 AliCheb3D.cxx:296
 AliCheb3D.cxx:297
 AliCheb3D.cxx:298
 AliCheb3D.cxx:299
 AliCheb3D.cxx:300
 AliCheb3D.cxx:301
 AliCheb3D.cxx:302
 AliCheb3D.cxx:303
 AliCheb3D.cxx:304
 AliCheb3D.cxx:305
 AliCheb3D.cxx:306
 AliCheb3D.cxx:307
 AliCheb3D.cxx:308
 AliCheb3D.cxx:309
 AliCheb3D.cxx:310
 AliCheb3D.cxx:311
 AliCheb3D.cxx:312
 AliCheb3D.cxx:313
 AliCheb3D.cxx:314
 AliCheb3D.cxx:315
 AliCheb3D.cxx:316
 AliCheb3D.cxx:317
 AliCheb3D.cxx:318
 AliCheb3D.cxx:319
 AliCheb3D.cxx:320
 AliCheb3D.cxx:321
 AliCheb3D.cxx:322
 AliCheb3D.cxx:323
 AliCheb3D.cxx:324
 AliCheb3D.cxx:325
 AliCheb3D.cxx:326
 AliCheb3D.cxx:327
 AliCheb3D.cxx:328
 AliCheb3D.cxx:329
 AliCheb3D.cxx:330
 AliCheb3D.cxx:331
 AliCheb3D.cxx:332
 AliCheb3D.cxx:333
 AliCheb3D.cxx:334
 AliCheb3D.cxx:335
 AliCheb3D.cxx:336
 AliCheb3D.cxx:337
 AliCheb3D.cxx:338
 AliCheb3D.cxx:339
 AliCheb3D.cxx:340
 AliCheb3D.cxx:341
 AliCheb3D.cxx:342
 AliCheb3D.cxx:343
 AliCheb3D.cxx:344
 AliCheb3D.cxx:345
 AliCheb3D.cxx:346
 AliCheb3D.cxx:347
 AliCheb3D.cxx:348
 AliCheb3D.cxx:349
 AliCheb3D.cxx:350
 AliCheb3D.cxx:351
 AliCheb3D.cxx:352
 AliCheb3D.cxx:353
 AliCheb3D.cxx:354
 AliCheb3D.cxx:355
 AliCheb3D.cxx:356
 AliCheb3D.cxx:357
 AliCheb3D.cxx:358
 AliCheb3D.cxx:359
 AliCheb3D.cxx:360
 AliCheb3D.cxx:361
 AliCheb3D.cxx:362
 AliCheb3D.cxx:363
 AliCheb3D.cxx:364
 AliCheb3D.cxx:365
 AliCheb3D.cxx:366
 AliCheb3D.cxx:367
 AliCheb3D.cxx:368
 AliCheb3D.cxx:369
 AliCheb3D.cxx:370
 AliCheb3D.cxx:371
 AliCheb3D.cxx:372
 AliCheb3D.cxx:373
 AliCheb3D.cxx:374
 AliCheb3D.cxx:375
 AliCheb3D.cxx:376
 AliCheb3D.cxx:377
 AliCheb3D.cxx:378
 AliCheb3D.cxx:379
 AliCheb3D.cxx:380
 AliCheb3D.cxx:381
 AliCheb3D.cxx:382
 AliCheb3D.cxx:383
 AliCheb3D.cxx:384
 AliCheb3D.cxx:385
 AliCheb3D.cxx:386
 AliCheb3D.cxx:387
 AliCheb3D.cxx:388
 AliCheb3D.cxx:389
 AliCheb3D.cxx:390
 AliCheb3D.cxx:391
 AliCheb3D.cxx:392
 AliCheb3D.cxx:393
 AliCheb3D.cxx:394
 AliCheb3D.cxx:395
 AliCheb3D.cxx:396
 AliCheb3D.cxx:397
 AliCheb3D.cxx:398
 AliCheb3D.cxx:399
 AliCheb3D.cxx:400
 AliCheb3D.cxx:401
 AliCheb3D.cxx:402
 AliCheb3D.cxx:403
 AliCheb3D.cxx:404
 AliCheb3D.cxx:405
 AliCheb3D.cxx:406
 AliCheb3D.cxx:407
 AliCheb3D.cxx:408
 AliCheb3D.cxx:409
 AliCheb3D.cxx:410
 AliCheb3D.cxx:411
 AliCheb3D.cxx:412
 AliCheb3D.cxx:413
 AliCheb3D.cxx:414
 AliCheb3D.cxx:415
 AliCheb3D.cxx:416
 AliCheb3D.cxx:417
 AliCheb3D.cxx:418
 AliCheb3D.cxx:419
 AliCheb3D.cxx:420
 AliCheb3D.cxx:421
 AliCheb3D.cxx:422
 AliCheb3D.cxx:423
 AliCheb3D.cxx:424
 AliCheb3D.cxx:425
 AliCheb3D.cxx:426
 AliCheb3D.cxx:427
 AliCheb3D.cxx:428
 AliCheb3D.cxx:429
 AliCheb3D.cxx:430
 AliCheb3D.cxx:431
 AliCheb3D.cxx:432
 AliCheb3D.cxx:433
 AliCheb3D.cxx:434
 AliCheb3D.cxx:435
 AliCheb3D.cxx:436
 AliCheb3D.cxx:437
 AliCheb3D.cxx:438
 AliCheb3D.cxx:439
 AliCheb3D.cxx:440
 AliCheb3D.cxx:441
 AliCheb3D.cxx:442
 AliCheb3D.cxx:443
 AliCheb3D.cxx:444
 AliCheb3D.cxx:445
 AliCheb3D.cxx:446
 AliCheb3D.cxx:447
 AliCheb3D.cxx:448
 AliCheb3D.cxx:449
 AliCheb3D.cxx:450
 AliCheb3D.cxx:451
 AliCheb3D.cxx:452
 AliCheb3D.cxx:453
 AliCheb3D.cxx:454
 AliCheb3D.cxx:455
 AliCheb3D.cxx:456
 AliCheb3D.cxx:457
 AliCheb3D.cxx:458
 AliCheb3D.cxx:459
 AliCheb3D.cxx:460
 AliCheb3D.cxx:461
 AliCheb3D.cxx:462
 AliCheb3D.cxx:463
 AliCheb3D.cxx:464
 AliCheb3D.cxx:465
 AliCheb3D.cxx:466
 AliCheb3D.cxx:467
 AliCheb3D.cxx:468
 AliCheb3D.cxx:469
 AliCheb3D.cxx:470
 AliCheb3D.cxx:471
 AliCheb3D.cxx:472
 AliCheb3D.cxx:473
 AliCheb3D.cxx:474
 AliCheb3D.cxx:475
 AliCheb3D.cxx:476
 AliCheb3D.cxx:477
 AliCheb3D.cxx:478
 AliCheb3D.cxx:479
 AliCheb3D.cxx:480
 AliCheb3D.cxx:481
 AliCheb3D.cxx:482
 AliCheb3D.cxx:483
 AliCheb3D.cxx:484
 AliCheb3D.cxx:485
 AliCheb3D.cxx:486
 AliCheb3D.cxx:487
 AliCheb3D.cxx:488
 AliCheb3D.cxx:489
 AliCheb3D.cxx:490
 AliCheb3D.cxx:491
 AliCheb3D.cxx:492
 AliCheb3D.cxx:493
 AliCheb3D.cxx:494
 AliCheb3D.cxx:495
 AliCheb3D.cxx:496
 AliCheb3D.cxx:497
 AliCheb3D.cxx:498
 AliCheb3D.cxx:499
 AliCheb3D.cxx:500
 AliCheb3D.cxx:501
 AliCheb3D.cxx:502
 AliCheb3D.cxx:503
 AliCheb3D.cxx:504
 AliCheb3D.cxx:505
 AliCheb3D.cxx:506
 AliCheb3D.cxx:507
 AliCheb3D.cxx:508
 AliCheb3D.cxx:509
 AliCheb3D.cxx:510
 AliCheb3D.cxx:511
 AliCheb3D.cxx:512
 AliCheb3D.cxx:513
 AliCheb3D.cxx:514
 AliCheb3D.cxx:515
 AliCheb3D.cxx:516
 AliCheb3D.cxx:517
 AliCheb3D.cxx:518
 AliCheb3D.cxx:519
 AliCheb3D.cxx:520
 AliCheb3D.cxx:521
 AliCheb3D.cxx:522
 AliCheb3D.cxx:523
 AliCheb3D.cxx:524
 AliCheb3D.cxx:525
 AliCheb3D.cxx:526
 AliCheb3D.cxx:527
 AliCheb3D.cxx:528
 AliCheb3D.cxx:529
 AliCheb3D.cxx:530
 AliCheb3D.cxx:531
 AliCheb3D.cxx:532
 AliCheb3D.cxx:533
 AliCheb3D.cxx:534
 AliCheb3D.cxx:535
 AliCheb3D.cxx:536
 AliCheb3D.cxx:537
 AliCheb3D.cxx:538
 AliCheb3D.cxx:539
 AliCheb3D.cxx:540
 AliCheb3D.cxx:541
 AliCheb3D.cxx:542
 AliCheb3D.cxx:543
 AliCheb3D.cxx:544
 AliCheb3D.cxx:545
 AliCheb3D.cxx:546
 AliCheb3D.cxx:547
 AliCheb3D.cxx:548
 AliCheb3D.cxx:549
 AliCheb3D.cxx:550
 AliCheb3D.cxx:551
 AliCheb3D.cxx:552
 AliCheb3D.cxx:553
 AliCheb3D.cxx:554
 AliCheb3D.cxx:555
 AliCheb3D.cxx:556
 AliCheb3D.cxx:557
 AliCheb3D.cxx:558
 AliCheb3D.cxx:559
 AliCheb3D.cxx:560
 AliCheb3D.cxx:561
 AliCheb3D.cxx:562
 AliCheb3D.cxx:563
 AliCheb3D.cxx:564
 AliCheb3D.cxx:565
 AliCheb3D.cxx:566
 AliCheb3D.cxx:567
 AliCheb3D.cxx:568
 AliCheb3D.cxx:569
 AliCheb3D.cxx:570
 AliCheb3D.cxx:571
 AliCheb3D.cxx:572
 AliCheb3D.cxx:573
 AliCheb3D.cxx:574
 AliCheb3D.cxx:575
 AliCheb3D.cxx:576
 AliCheb3D.cxx:577
 AliCheb3D.cxx:578
 AliCheb3D.cxx:579
 AliCheb3D.cxx:580
 AliCheb3D.cxx:581
 AliCheb3D.cxx:582
 AliCheb3D.cxx:583
 AliCheb3D.cxx:584
 AliCheb3D.cxx:585
 AliCheb3D.cxx:586
 AliCheb3D.cxx:587
 AliCheb3D.cxx:588
 AliCheb3D.cxx:589
 AliCheb3D.cxx:590
 AliCheb3D.cxx:591
 AliCheb3D.cxx:592
 AliCheb3D.cxx:593
 AliCheb3D.cxx:594
 AliCheb3D.cxx:595
 AliCheb3D.cxx:596
 AliCheb3D.cxx:597
 AliCheb3D.cxx:598
 AliCheb3D.cxx:599
 AliCheb3D.cxx:600
 AliCheb3D.cxx:601
 AliCheb3D.cxx:602
 AliCheb3D.cxx:603
 AliCheb3D.cxx:604
 AliCheb3D.cxx:605
 AliCheb3D.cxx:606
 AliCheb3D.cxx:607
 AliCheb3D.cxx:608
 AliCheb3D.cxx:609
 AliCheb3D.cxx:610
 AliCheb3D.cxx:611
 AliCheb3D.cxx:612
 AliCheb3D.cxx:613
 AliCheb3D.cxx:614
 AliCheb3D.cxx:615
 AliCheb3D.cxx:616
 AliCheb3D.cxx:617
 AliCheb3D.cxx:618
 AliCheb3D.cxx:619
 AliCheb3D.cxx:620
 AliCheb3D.cxx:621
 AliCheb3D.cxx:622
 AliCheb3D.cxx:623
 AliCheb3D.cxx:624
 AliCheb3D.cxx:625
 AliCheb3D.cxx:626
 AliCheb3D.cxx:627
 AliCheb3D.cxx:628
 AliCheb3D.cxx:629
 AliCheb3D.cxx:630
 AliCheb3D.cxx:631
 AliCheb3D.cxx:632
 AliCheb3D.cxx:633
 AliCheb3D.cxx:634
 AliCheb3D.cxx:635
 AliCheb3D.cxx:636
 AliCheb3D.cxx:637
 AliCheb3D.cxx:638
 AliCheb3D.cxx:639
 AliCheb3D.cxx:640
 AliCheb3D.cxx:641
 AliCheb3D.cxx:642
 AliCheb3D.cxx:643
 AliCheb3D.cxx:644
 AliCheb3D.cxx:645
 AliCheb3D.cxx:646
 AliCheb3D.cxx:647
 AliCheb3D.cxx:648
 AliCheb3D.cxx:649
 AliCheb3D.cxx:650
 AliCheb3D.cxx:651
 AliCheb3D.cxx:652
 AliCheb3D.cxx:653
 AliCheb3D.cxx:654
 AliCheb3D.cxx:655
 AliCheb3D.cxx:656
 AliCheb3D.cxx:657
 AliCheb3D.cxx:658
 AliCheb3D.cxx:659
 AliCheb3D.cxx:660
 AliCheb3D.cxx:661
 AliCheb3D.cxx:662
 AliCheb3D.cxx:663
 AliCheb3D.cxx:664
 AliCheb3D.cxx:665
 AliCheb3D.cxx:666
 AliCheb3D.cxx:667
 AliCheb3D.cxx:668
 AliCheb3D.cxx:669
 AliCheb3D.cxx:670
 AliCheb3D.cxx:671
 AliCheb3D.cxx:672
 AliCheb3D.cxx:673
 AliCheb3D.cxx:674
 AliCheb3D.cxx:675
 AliCheb3D.cxx:676
 AliCheb3D.cxx:677
 AliCheb3D.cxx:678
 AliCheb3D.cxx:679
 AliCheb3D.cxx:680
 AliCheb3D.cxx:681
 AliCheb3D.cxx:682
 AliCheb3D.cxx:683
 AliCheb3D.cxx:684
 AliCheb3D.cxx:685
 AliCheb3D.cxx:686
 AliCheb3D.cxx:687
 AliCheb3D.cxx:688
 AliCheb3D.cxx:689
 AliCheb3D.cxx:690
 AliCheb3D.cxx:691
 AliCheb3D.cxx:692
 AliCheb3D.cxx:693
 AliCheb3D.cxx:694
 AliCheb3D.cxx:695
 AliCheb3D.cxx:696
 AliCheb3D.cxx:697
 AliCheb3D.cxx:698
 AliCheb3D.cxx:699
 AliCheb3D.cxx:700
 AliCheb3D.cxx:701
 AliCheb3D.cxx:702
 AliCheb3D.cxx:703
 AliCheb3D.cxx:704
 AliCheb3D.cxx:705
 AliCheb3D.cxx:706
 AliCheb3D.cxx:707
 AliCheb3D.cxx:708
 AliCheb3D.cxx:709
 AliCheb3D.cxx:710
 AliCheb3D.cxx:711
 AliCheb3D.cxx:712
 AliCheb3D.cxx:713
 AliCheb3D.cxx:714
 AliCheb3D.cxx:715
 AliCheb3D.cxx:716
 AliCheb3D.cxx:717
 AliCheb3D.cxx:718
 AliCheb3D.cxx:719
 AliCheb3D.cxx:720
 AliCheb3D.cxx:721
 AliCheb3D.cxx:722
 AliCheb3D.cxx:723
 AliCheb3D.cxx:724
 AliCheb3D.cxx:725
 AliCheb3D.cxx:726
 AliCheb3D.cxx:727
 AliCheb3D.cxx:728
 AliCheb3D.cxx:729
 AliCheb3D.cxx:730
 AliCheb3D.cxx:731
 AliCheb3D.cxx:732
 AliCheb3D.cxx:733
 AliCheb3D.cxx:734
 AliCheb3D.cxx:735
 AliCheb3D.cxx:736
 AliCheb3D.cxx:737
 AliCheb3D.cxx:738
 AliCheb3D.cxx:739
 AliCheb3D.cxx:740
 AliCheb3D.cxx:741
 AliCheb3D.cxx:742
 AliCheb3D.cxx:743
 AliCheb3D.cxx:744
 AliCheb3D.cxx:745
 AliCheb3D.cxx:746
 AliCheb3D.cxx:747
 AliCheb3D.cxx:748
 AliCheb3D.cxx:749
 AliCheb3D.cxx:750
 AliCheb3D.cxx:751
 AliCheb3D.cxx:752
 AliCheb3D.cxx:753
 AliCheb3D.cxx:754
 AliCheb3D.cxx:755
 AliCheb3D.cxx:756
 AliCheb3D.cxx:757
 AliCheb3D.cxx:758
 AliCheb3D.cxx:759
 AliCheb3D.cxx:760
 AliCheb3D.cxx:761
 AliCheb3D.cxx:762
 AliCheb3D.cxx:763
 AliCheb3D.cxx:764
 AliCheb3D.cxx:765
 AliCheb3D.cxx:766
 AliCheb3D.cxx:767
 AliCheb3D.cxx:768
 AliCheb3D.cxx:769
 AliCheb3D.cxx:770
 AliCheb3D.cxx:771
 AliCheb3D.cxx:772
 AliCheb3D.cxx:773
 AliCheb3D.cxx:774
 AliCheb3D.cxx:775
 AliCheb3D.cxx:776
 AliCheb3D.cxx:777
 AliCheb3D.cxx:778
 AliCheb3D.cxx:779
 AliCheb3D.cxx:780
 AliCheb3D.cxx:781
 AliCheb3D.cxx:782
 AliCheb3D.cxx:783
 AliCheb3D.cxx:784
 AliCheb3D.cxx:785
 AliCheb3D.cxx:786
 AliCheb3D.cxx:787
 AliCheb3D.cxx:788
 AliCheb3D.cxx:789
 AliCheb3D.cxx:790
 AliCheb3D.cxx:791
 AliCheb3D.cxx:792
 AliCheb3D.cxx:793
 AliCheb3D.cxx:794
 AliCheb3D.cxx:795
 AliCheb3D.cxx:796
 AliCheb3D.cxx:797
 AliCheb3D.cxx:798
 AliCheb3D.cxx:799
 AliCheb3D.cxx:800
 AliCheb3D.cxx:801
 AliCheb3D.cxx:802
 AliCheb3D.cxx:803
 AliCheb3D.cxx:804
 AliCheb3D.cxx:805
 AliCheb3D.cxx:806
 AliCheb3D.cxx:807
 AliCheb3D.cxx:808
 AliCheb3D.cxx:809
 AliCheb3D.cxx:810
 AliCheb3D.cxx:811
 AliCheb3D.cxx:812
 AliCheb3D.cxx:813
 AliCheb3D.cxx:814
 AliCheb3D.cxx:815
 AliCheb3D.cxx:816
 AliCheb3D.cxx:817
 AliCheb3D.cxx:818
 AliCheb3D.cxx:819
 AliCheb3D.cxx:820
 AliCheb3D.cxx:821
 AliCheb3D.cxx:822
 AliCheb3D.cxx:823
 AliCheb3D.cxx:824
 AliCheb3D.cxx:825
 AliCheb3D.cxx:826
 AliCheb3D.cxx:827
 AliCheb3D.cxx:828
 AliCheb3D.cxx:829
 AliCheb3D.cxx:830
 AliCheb3D.cxx:831
 AliCheb3D.cxx:832
 AliCheb3D.cxx:833
 AliCheb3D.cxx:834
 AliCheb3D.cxx:835
 AliCheb3D.cxx:836
 AliCheb3D.cxx:837
 AliCheb3D.cxx:838
 AliCheb3D.cxx:839
 AliCheb3D.cxx:840
 AliCheb3D.cxx:841
 AliCheb3D.cxx:842
 AliCheb3D.cxx:843
 AliCheb3D.cxx:844
 AliCheb3D.cxx:845
 AliCheb3D.cxx:846
 AliCheb3D.cxx:847
 AliCheb3D.cxx:848
 AliCheb3D.cxx:849
 AliCheb3D.cxx:850
 AliCheb3D.cxx:851
 AliCheb3D.cxx:852
 AliCheb3D.cxx:853
 AliCheb3D.cxx:854
 AliCheb3D.cxx:855
 AliCheb3D.cxx:856
 AliCheb3D.cxx:857
 AliCheb3D.cxx:858
 AliCheb3D.cxx:859
 AliCheb3D.cxx:860
 AliCheb3D.cxx:861
 AliCheb3D.cxx:862
 AliCheb3D.cxx:863
 AliCheb3D.cxx:864
 AliCheb3D.cxx:865
 AliCheb3D.cxx:866
 AliCheb3D.cxx:867
 AliCheb3D.cxx:868
 AliCheb3D.cxx:869
 AliCheb3D.cxx:870
 AliCheb3D.cxx:871
 AliCheb3D.cxx:872
 AliCheb3D.cxx:873
 AliCheb3D.cxx:874
 AliCheb3D.cxx:875
 AliCheb3D.cxx:876
 AliCheb3D.cxx:877
 AliCheb3D.cxx:878
 AliCheb3D.cxx:879
 AliCheb3D.cxx:880
 AliCheb3D.cxx:881
 AliCheb3D.cxx:882
 AliCheb3D.cxx:883
 AliCheb3D.cxx:884
 AliCheb3D.cxx:885
 AliCheb3D.cxx:886
 AliCheb3D.cxx:887
 AliCheb3D.cxx:888
 AliCheb3D.cxx:889
 AliCheb3D.cxx:890
 AliCheb3D.cxx:891
 AliCheb3D.cxx:892
 AliCheb3D.cxx:893
 AliCheb3D.cxx:894
 AliCheb3D.cxx:895
 AliCheb3D.cxx:896
 AliCheb3D.cxx:897
 AliCheb3D.cxx:898
 AliCheb3D.cxx:899
 AliCheb3D.cxx:900
 AliCheb3D.cxx:901
 AliCheb3D.cxx:902
 AliCheb3D.cxx:903
 AliCheb3D.cxx:904
 AliCheb3D.cxx:905
 AliCheb3D.cxx:906
 AliCheb3D.cxx:907
 AliCheb3D.cxx:908
 AliCheb3D.cxx:909
 AliCheb3D.cxx:910
 AliCheb3D.cxx:911
 AliCheb3D.cxx:912
 AliCheb3D.cxx:913
 AliCheb3D.cxx:914
 AliCheb3D.cxx:915
 AliCheb3D.cxx:916
 AliCheb3D.cxx:917
 AliCheb3D.cxx:918
 AliCheb3D.cxx:919
 AliCheb3D.cxx:920
 AliCheb3D.cxx:921
 AliCheb3D.cxx:922
 AliCheb3D.cxx:923
 AliCheb3D.cxx:924
 AliCheb3D.cxx:925
 AliCheb3D.cxx:926
 AliCheb3D.cxx:927
 AliCheb3D.cxx:928
 AliCheb3D.cxx:929
 AliCheb3D.cxx:930
 AliCheb3D.cxx:931
 AliCheb3D.cxx:932
 AliCheb3D.cxx:933
 AliCheb3D.cxx:934
 AliCheb3D.cxx:935
 AliCheb3D.cxx:936
 AliCheb3D.cxx:937
 AliCheb3D.cxx:938
 AliCheb3D.cxx:939
 AliCheb3D.cxx:940
 AliCheb3D.cxx:941
 AliCheb3D.cxx:942
 AliCheb3D.cxx:943
 AliCheb3D.cxx:944
 AliCheb3D.cxx:945
 AliCheb3D.cxx:946
 AliCheb3D.cxx:947
 AliCheb3D.cxx:948
 AliCheb3D.cxx:949
 AliCheb3D.cxx:950
 AliCheb3D.cxx:951
 AliCheb3D.cxx:952
 AliCheb3D.cxx:953
 AliCheb3D.cxx:954
 AliCheb3D.cxx:955
 AliCheb3D.cxx:956
 AliCheb3D.cxx:957
 AliCheb3D.cxx:958
 AliCheb3D.cxx:959
 AliCheb3D.cxx:960
 AliCheb3D.cxx:961
 AliCheb3D.cxx:962
 AliCheb3D.cxx:963
 AliCheb3D.cxx:964
 AliCheb3D.cxx:965
 AliCheb3D.cxx:966
 AliCheb3D.cxx:967
 AliCheb3D.cxx:968
 AliCheb3D.cxx:969
 AliCheb3D.cxx:970
 AliCheb3D.cxx:971
 AliCheb3D.cxx:972
 AliCheb3D.cxx:973
 AliCheb3D.cxx:974
 AliCheb3D.cxx:975
 AliCheb3D.cxx:976
 AliCheb3D.cxx:977
 AliCheb3D.cxx:978
 AliCheb3D.cxx:979
 AliCheb3D.cxx:980
 AliCheb3D.cxx:981
 AliCheb3D.cxx:982
 AliCheb3D.cxx:983
 AliCheb3D.cxx:984