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


///////////////////////////////////////////////////////////////////////////
// Class AliMathBase
// 
// Subset of  matheamtical functions  not included in the TMath
//

///////////////////////////////////////////////////////////////////////////
#include "TMath.h"
#include "AliMathBase.h"
#include "Riostream.h"
#include "TH1F.h"
#include "TH3.h"
#include "TF1.h"
#include "TLinearFitter.h"
#include "AliLog.h"

#include "AliExternalTrackParam.h"

//
// includes neccessary for test functions
//

#include "TSystem.h"
#include "TRandom.h"
#include "TStopwatch.h"
#include "TTreeStream.h"
 
ClassImp(AliMathBase) // Class implementation to enable ROOT I/O
 
AliMathBase::AliMathBase() : TObject()
{
  //
  // Default constructor
  //
}
///////////////////////////////////////////////////////////////////////////
AliMathBase::~AliMathBase()
{
  //
  // Destructor
  //
}


//_____________________________________________________________________________
void AliMathBase::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
                           , Double_t &sigma, Int_t hh)
{
  //
  // Robust estimator in 1D case MI version - (faster than ROOT version)
  //
  // For the univariate case
  // estimates of location and scatter are returned in mean and sigma parameters
  // the algorithm works on the same principle as in multivariate case -
  // it finds a subset of size hh with smallest sigma, and then returns mean and
  // sigma of this subset
  //

  if (nvectors<2) {
    AliErrorClass(Form("nvectors = %d, should be > 1",nvectors));
    return;
  }
  if (hh<2)
    hh=(nvectors+2)/2;
  Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
  Int_t *index=new Int_t[nvectors];
  TMath::Sort(nvectors, data, index, kFALSE);
  
  Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
  Double_t factor = faclts[TMath::Max(0,nquant-1)];
  
  Double_t sumx  =0;
  Double_t sumx2 =0;
  Int_t    bestindex = -1;
  Double_t bestmean  = 0; 
  Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.);   // maximal possible sigma
  bestsigma *=bestsigma;

  for (Int_t i=0; i<hh; i++){
    sumx  += data[index[i]];
    sumx2 += data[index[i]]*data[index[i]];
  }
  
  Double_t norm = 1./Double_t(hh);
  Double_t norm2 = 1./Double_t(hh-1);
  for (Int_t i=hh; i<nvectors; i++){
    Double_t cmean  = sumx*norm;
    Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
    if (csigma<bestsigma){
      bestmean  = cmean;
      bestsigma = csigma;
      bestindex = i-hh;
    }
    
    sumx  += data[index[i]]-data[index[i-hh]];
    sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
  }
  
  Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
  mean  = bestmean;
  sigma = bstd;
  delete [] index;

}



void AliMathBase::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh,  Float_t externalfactor)
{
  // Modified version of ROOT robust EvaluateUni
  // robust estimator in 1D case MI version
  // added external factor to include precision of external measurement
  // 

  if (hh==0)
    hh=(nvectors+2)/2;
  Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
  Int_t *index=new Int_t[nvectors];
  TMath::Sort(nvectors, data, index, kFALSE);
  //
  Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
  Double_t factor = faclts[0];
  if (nquant>0){
    // fix proper normalization - Anja
    factor = faclts[nquant-1];
  }

  //
  //
  Double_t sumx  =0;
  Double_t sumx2 =0;
  Int_t    bestindex = -1;
  Double_t bestmean  = 0; 
  Double_t bestsigma = -1;
  for (Int_t i=0; i<hh; i++){
    sumx  += data[index[i]];
    sumx2 += data[index[i]]*data[index[i]];
  }
  //   
  Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
  Double_t norm = 1./Double_t(hh);
  for (Int_t i=hh; i<nvectors; i++){
    Double_t cmean  = sumx*norm;
    Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
    if (csigma<bestsigma ||  bestsigma<0){
      bestmean  = cmean;
      bestsigma = csigma;
      bestindex = i-hh;
    }
    //
    //
    sumx  += data[index[i]]-data[index[i-hh]];
    sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
  }
  
  Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
  mean  = bestmean;
  sigma = bstd;
  delete [] index;
}


//_____________________________________________________________________________
Int_t AliMathBase::Freq(Int_t n, const Int_t *inlist
                        , Int_t *outlist, Bool_t down)
{    
  //
  //  Sort eleements according occurancy 
  //  The size of output array has is 2*n 
  //

  Int_t * sindexS = new Int_t[n];     // temp array for sorting
  Int_t * sindexF = new Int_t[2*n];   
  for (Int_t i=0;i<n;i++) sindexF[i]=0;
  //
  TMath::Sort(n,inlist, sindexS, down);  
  Int_t last      = inlist[sindexS[0]];
  Int_t val       = last;
  sindexF[0]      = 1;
  sindexF[0+n]    = last;
  Int_t countPos  = 0;
  //
  //  find frequency
  for(Int_t i=1;i<n; i++){
    val = inlist[sindexS[i]];
    if (last == val)   sindexF[countPos]++;
    else{      
      countPos++;
      sindexF[countPos+n] = val;
      sindexF[countPos]++;
      last =val;
    }
  }
  if (last==val) countPos++;
  // sort according frequency
  TMath::Sort(countPos, sindexF, sindexS, kTRUE);
  for (Int_t i=0;i<countPos;i++){
    outlist[2*i  ] = sindexF[sindexS[i]+n];
    outlist[2*i+1] = sindexF[sindexS[i]];
  }
  delete [] sindexS;
  delete [] sindexF;
  
  return countPos;

}

//___AliMathBase__________________________________________________________________________
void AliMathBase::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
  //
  //
  //
  Int_t nbins    = his->GetNbinsX();
  Float_t nentries = his->GetEntries();
  Float_t sum      =0;
  Float_t mean   = 0;
  Float_t sigma2 = 0;
  Float_t ncumul=0;  
  for (Int_t ibin=1;ibin<nbins; ibin++){
    ncumul+= his->GetBinContent(ibin);
    Float_t fraction = Float_t(ncumul)/Float_t(nentries);
    if (fraction>down && fraction<up){
      sum+=his->GetBinContent(ibin);
      mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
      sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);      
    }
  }
  mean/=sum;
  sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
  if (param){
    (*param)[0] = his->GetMaximum();
    (*param)[1] = mean;
    (*param)[2] = sigma2;
    
  }
  if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
}

void AliMathBase::LTM(TH1F * his, TVectorD *param , Float_t fraction,  Bool_t verbose){
  //
  // LTM
  //
  Int_t nbins    = his->GetNbinsX();
  Int_t nentries = (Int_t)his->GetEntries();
  Double_t *data  = new Double_t[nentries];
  Int_t npoints=0;
  for (Int_t ibin=1;ibin<nbins; ibin++){
    Float_t entriesI = his->GetBinContent(ibin);
    Float_t xcenter= his->GetBinCenter(ibin);
    for (Int_t ic=0; ic<entriesI; ic++){
      if (npoints<nentries){
	data[npoints]= xcenter;
	npoints++;
      }
    }
  }
  Double_t mean, sigma;
  Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
  npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
  AliMathBase::EvaluateUni(npoints, data, mean,sigma,npoints2);
  delete [] data;
  if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
    (*param)[0] = his->GetMaximum();
    (*param)[1] = mean;
    (*param)[2] = sigma;    
  }
}

Double_t  AliMathBase::FitGaus(TH1F* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
  //
  //  Fit histogram with gaussian function
  //  
  //  Prameters:
  //       return value- chi2 - if negative ( not enough points)
  //       his        -  input histogram
  //       param      -  vector with parameters 
  //       xmin, xmax -  range to fit - if xmin=xmax=0 - the full histogram range used
  //  Fitting:
  //  1. Step - make logarithm
  //  2. Linear  fit (parabola) - more robust - always converge
  //  3. In case of small statistic bins are averaged
  //  
  static TLinearFitter fitter(3,"pol2");
  TVectorD  par(3);
  TVectorD  sigma(3);
  TMatrixD mat(3,3);
  if (his->GetMaximum()<4) return -1;  
  if (his->GetEntries()<12) return -1;  
  if (his->GetRMS()<mat.GetTol()) return -1;
  Float_t maxEstimate   = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
  Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));

  if (maxEstimate<1) return -1;
  Int_t nbins    = his->GetNbinsX();
  Int_t npoints=0;
  //


  if (xmin>=xmax){
    xmin = his->GetXaxis()->GetXmin();
    xmax = his->GetXaxis()->GetXmax();
  }
  for (Int_t iter=0; iter<2; iter++){
    fitter.ClearPoints();
    npoints=0;
    for (Int_t ibin=1;ibin<nbins+1; ibin++){
      Int_t countB=1;
      Float_t entriesI =  his->GetBinContent(ibin);
      for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
	if (ibin+delta>1 &&ibin+delta<nbins-1){
	  entriesI +=  his->GetBinContent(ibin+delta);
	  countB++;
	}
      }
      entriesI/=countB;
      Double_t xcenter= his->GetBinCenter(ibin);
      if (xcenter<xmin || xcenter>xmax) continue;
      Double_t error=1./TMath::Sqrt(countB);
      Float_t   cont=2;
      if (iter>0){
	if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
	cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
	if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
      }
      if (entriesI>1&&cont>1){
	fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
	npoints++;
      }
    }  
    if (npoints>3){
      fitter.Eval();
      fitter.GetParameters(par);
    }else{
      break;
    }
  }
  if (npoints<=3){
    return -1;
  }
  fitter.GetParameters(par);
  fitter.GetCovarianceMatrix(mat);
  if (TMath::Abs(par[1])<mat.GetTol()) return -1;
  if (TMath::Abs(par[2])<mat.GetTol()) return -1;
  Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
  //fitter.GetParameters();
  if (!param)  param  = new TVectorD(3);
  //if (!matrix) matrix = new TMatrixD(3,3);
  (*param)[1] = par[1]/(-2.*par[2]);
  (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
  (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1]);
  if (verbose){
    par.Print();
    mat.Print();
    param->Print();
    printf("Chi2=%f\n",chi2);
    TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
    f1->SetParameter(0, (*param)[0]);
    f1->SetParameter(1, (*param)[1]);
    f1->SetParameter(2, (*param)[2]);    
    f1->Draw("same");
  }
  return chi2;
}

Double_t  AliMathBase::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
  //
  //  Fit histogram with gaussian function
  //  
  //  Prameters:
  //     nbins: size of the array and number of histogram bins
  //     xMin, xMax: histogram range
  //     param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma, 3-Sum)
  //     matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
  //
  //  Return values:
  //    >0: the chi2 returned by TLinearFitter
  //    -3: only three points have been used for the calculation - no fitter was used
  //    -2: only two points have been used for the calculation - center of gravity was uesed for calculation
  //    -1: only one point has been used for the calculation - center of gravity was uesed for calculation
  //    -4: invalid result!!
  //
  //  Fitting:
  //  1. Step - make logarithm
  //  2. Linear  fit (parabola) - more robust - always converge
  //  
  static TLinearFitter fitter(3,"pol2");
  static TMatrixD mat(3,3);
  static Double_t kTol = mat.GetTol();
  fitter.StoreData(kFALSE);
  fitter.ClearPoints();
  TVectorD  par(3);
  TVectorD  sigma(3);
  TMatrixD A(3,3);
  TMatrixD b(3,1);
  Float_t rms = TMath::RMS(nBins,arr);
  Float_t max = TMath::MaxElement(nBins,arr);
  Float_t binWidth = (xMax-xMin)/(Float_t)nBins;

  Float_t meanCOG = 0;
  Float_t rms2COG = 0;
  Float_t sumCOG  = 0;

  Float_t entries = 0;
  Int_t nfilled=0;

  if (!param)  param  = new TVectorD(4);

  for (Int_t i=0; i<nBins; i++){
      entries+=arr[i];
      if (arr[i]>0) nfilled++;
  }
  (*param)[0] = 0;
  (*param)[1] = 0;
  (*param)[2] = 0;
  (*param)[3] = 0;

  if (max<4) return -4;
  if (entries<12) return -4;
  if (rms<kTol) return -4;

  (*param)[3] = entries;

  Int_t npoints=0;
  for (Int_t ibin=0;ibin<nBins; ibin++){
      Float_t entriesI = arr[ibin];
    if (entriesI>1){
      Double_t xcenter = xMin+(ibin+0.5)*binWidth;
      Float_t error    = 1./TMath::Sqrt(entriesI);
      Float_t val = TMath::Log(Float_t(entriesI));
      fitter.AddPoint(&xcenter,val,error);
      if (npoints<3){
	  A(npoints,0)=1;
	  A(npoints,1)=xcenter;
	  A(npoints,2)=xcenter*xcenter;
	  b(npoints,0)=val;
	  meanCOG+=xcenter*entriesI;
	  rms2COG +=xcenter*entriesI*xcenter;
	  sumCOG +=entriesI;
      }
      npoints++;
    }
  }
  
  Double_t chi2 = 0;
  if (npoints>=3){
      if ( npoints == 3 ){
	  //analytic calculation of the parameters for three points
	  A.Invert();
	  TMatrixD res(1,3);
	  res.Mult(A,b);
	  par[0]=res(0,0);
	  par[1]=res(0,1);
	  par[2]=res(0,2);
          chi2 = -3.;
      } else {
          // use fitter for more than three points
	  fitter.Eval();
	  fitter.GetParameters(par);
	  fitter.GetCovarianceMatrix(mat);
	  chi2 = fitter.GetChisquare()/Float_t(npoints);
      }
      if (TMath::Abs(par[1])<kTol) return -4;
      if (TMath::Abs(par[2])<kTol) return -4;

      //if (!param)  param  = new TVectorD(4);
      if ( param->GetNrows()<4 ) param->ResizeTo(4);
      //if (!matrix) matrix = new TMatrixD(3,3);  // !!!!might be a memory leek. use dummy matrix pointer to call this function!

      (*param)[1] = par[1]/(-2.*par[2]);
      (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
      Double_t lnparam0 = par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1];
      if ( lnparam0>307 ) return -4;
      (*param)[0] = TMath::Exp(lnparam0);
      if (verbose){
	  par.Print();
	  mat.Print();
	  param->Print();
	  printf("Chi2=%f\n",chi2);
	  TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
	  f1->SetParameter(0, (*param)[0]);
	  f1->SetParameter(1, (*param)[1]);
	  f1->SetParameter(2, (*param)[2]);
	  f1->Draw("same");
      }
      return chi2;
  }

  if (npoints == 2){
      //use center of gravity for 2 points
      meanCOG/=sumCOG;
      rms2COG /=sumCOG;
      (*param)[0] = max;
      (*param)[1] = meanCOG;
      (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
      chi2=-2.;
  }
  if ( npoints == 1 ){
      meanCOG/=sumCOG;
      (*param)[0] = max;
      (*param)[1] = meanCOG;
      (*param)[2] = binWidth/TMath::Sqrt(12);
      chi2=-1.;
  }
  return chi2;

}


Float_t AliMathBase::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
{
    //
    //  calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
    //  return COG; in case of failure return xMin
    //
    Float_t meanCOG = 0;
    Float_t rms2COG = 0;
    Float_t sumCOG  = 0;
    Int_t npoints   = 0;

    Float_t binWidth = (xMax-xMin)/(Float_t)nBins;

    for (Int_t ibin=0; ibin<nBins; ibin++){
	Float_t entriesI = (Float_t)arr[ibin];
	Double_t xcenter = xMin+(ibin+0.5)*binWidth;
	if ( entriesI>0 ){
	    meanCOG += xcenter*entriesI;
	    rms2COG += xcenter*entriesI*xcenter;
	    sumCOG  += entriesI;
	    npoints++;
	}
    }
    if ( sumCOG == 0 ) return xMin;
    meanCOG/=sumCOG;

    if ( rms ){
	rms2COG /=sumCOG;
	(*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
	if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
    }

    if ( sum )
        (*sum) = sumCOG;

    return meanCOG;
}


Double_t AliMathBase::ErfcFast(Double_t x){
  // Fast implementation of the complementary error function
  // The error of the approximation is |eps(x)| < 5E-4
  // See Abramowitz and Stegun, p.299, 7.1.27

  Double_t z = TMath::Abs(x);
  Double_t ans = 1+z*(0.278393+z*(0.230389+z*(0.000972+z*0.078108)));
  ans = 1.0/ans;
  ans *= ans;
  ans *= ans;

  return (x>=0.0 ? ans : 2.0 - ans);
}

///////////////////////////////////////////////////////////////
//////////////         TEST functions /////////////////////////
///////////////////////////////////////////////////////////////





void AliMathBase::TestGausFit(Int_t nhistos){
  //
  // Test performance of the parabolic - gaussian fit - compare it with 
  // ROOT gauss fit
  //  nhistos - number of histograms to be used for test
  //
  TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
  
  Float_t  *xTrue = new Float_t[nhistos];
  Float_t  *sTrue = new Float_t[nhistos];
  TVectorD **par1  = new TVectorD*[nhistos];
  TVectorD **par2  = new TVectorD*[nhistos];
  TMatrixD dummy(3,3);
  
  
  TH1F **h1f = new TH1F*[nhistos];
  TF1  *myg = new TF1("myg","gaus");
  TF1  *fit = new TF1("fit","gaus");
  gRandom->SetSeed(0);
  
  //init
  for (Int_t i=0;i<nhistos; i++){
    par1[i] = new TVectorD(3);
    par2[i] = new TVectorD(3);
    h1f[i]  = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
    xTrue[i]= gRandom->Rndm();
    gSystem->Sleep(2);
    sTrue[i]= .75+gRandom->Rndm()*.5;
    myg->SetParameters(1,xTrue[i],sTrue[i]);
    h1f[i]->FillRandom("myg");
  }
  
  TStopwatch s;
  s.Start();
  //standard gaus fit
  for (Int_t i=0; i<nhistos; i++){
    h1f[i]->Fit(fit,"0q");
    (*par1[i])(0) = fit->GetParameter(0);
    (*par1[i])(1) = fit->GetParameter(1);
    (*par1[i])(2) = fit->GetParameter(2);
  }
  s.Stop();
  printf("Gaussian fit\t");
  s.Print();
  
  s.Start();
  //AliMathBase gaus fit
  for (Int_t i=0; i<nhistos; i++){
    AliMathBase::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
  }
  
  s.Stop();
  printf("Parabolic fit\t");
  s.Print();
  //write stream
  for (Int_t i=0;i<nhistos; i++){
    Float_t xt  = xTrue[i];
    Float_t st  = sTrue[i];
    (*pcstream)<<"data"
	       <<"xTrue="<<xt
	       <<"sTrue="<<st
	       <<"pg.="<<(par1[i])
	       <<"pa.="<<(par2[i])
	       <<"\n";
  }    
  //delete pointers
  for (Int_t i=0;i<nhistos; i++){
    delete par1[i];
    delete par2[i];
    delete h1f[i];
  }
  delete pcstream;
  delete []h1f;
  delete []xTrue;
  delete []sTrue;
  //
  delete []par1;
  delete []par2;

}



TGraph2D * AliMathBase::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
  //
  //
  //
  // delta - number of bins to integrate
  // type - 0 - mean value

  TAxis * xaxis  = his->GetXaxis();
  TAxis * yaxis  = his->GetYaxis();
  //  TAxis * zaxis  = his->GetZaxis();
  Int_t   nbinx  = xaxis->GetNbins();
  Int_t   nbiny  = yaxis->GetNbins();
  const Int_t nc=1000;
  char name[nc];
  Int_t icount=0;
  TGraph2D  *graph = new TGraph2D(nbinx*nbiny);
  TF1 f1("f1","gaus");
  for (Int_t ix=0; ix<nbinx;ix++)
    for (Int_t iy=0; iy<nbiny;iy++){
      Float_t xcenter = xaxis->GetBinCenter(ix); 
      Float_t ycenter = yaxis->GetBinCenter(iy); 
      snprintf(name,nc,"%s_%d_%d",his->GetName(), ix,iy);
      TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
      Float_t stat= 0;
      if (type==0) stat = projection->GetMean();
      if (type==1) stat = projection->GetRMS();
      if (type==2 || type==3){
	TVectorD vec(3);
	AliMathBase::LTM((TH1F*)projection,&vec,0.7);
	if (type==2) stat= vec[1];
	if (type==3) stat= vec[0];	
      }
      if (type==4|| type==5){
	projection->Fit(&f1);
	if (type==4) stat= f1.GetParameter(1);
	if (type==5) stat= f1.GetParameter(2);
      }
      //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
      graph->SetPoint(icount,xcenter, ycenter, stat);
      icount++;
    }
  return graph;
}

TGraph * AliMathBase::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
  //
  //
  //
  // delta - number of bins to integrate
  // type - 0 - mean value

  TAxis * xaxis  = his->GetXaxis();
  TAxis * yaxis  = his->GetYaxis();
  //  TAxis * zaxis  = his->GetZaxis();
  Int_t   nbinx  = xaxis->GetNbins();
  Int_t   nbiny  = yaxis->GetNbins();
  const Int_t nc=1000;
  char name[nc];
  Int_t icount=0;
  TGraph  *graph = new TGraph(nbinx);
  TF1 f1("f1","gaus");
  for (Int_t ix=0; ix<nbinx;ix++){
    Float_t xcenter = xaxis->GetBinCenter(ix); 
    //    Float_t ycenter = yaxis->GetBinCenter(iy); 
    snprintf(name,nc,"%s_%d",his->GetName(), ix);
    TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
    Float_t stat= 0;
    if (type==0) stat = projection->GetMean();
    if (type==1) stat = projection->GetRMS();
    if (type==2 || type==3){
      TVectorD vec(3);
	AliMathBase::LTM((TH1F*)projection,&vec,0.7);
	if (type==2) stat= vec[1];
	if (type==3) stat= vec[0];	
    }
    if (type==4|| type==5){
      projection->Fit(&f1);
      if (type==4) stat= f1.GetParameter(1);
      if (type==5) stat= f1.GetParameter(2);
    }
      //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
    graph->SetPoint(icount,xcenter, stat);
    icount++;
  }
  return graph;
}

Double_t AliMathBase::TruncatedGaus(Double_t mean, Double_t sigma, Double_t cutat)
{
  // return number generated according to a gaussian distribution N(mean,sigma) truncated at cutat
  //
  Double_t value;
  do{
    value=gRandom->Gaus(mean,sigma);
  }while(TMath::Abs(value-mean)>cutat);
  return value;
}

Double_t AliMathBase::TruncatedGaus(Double_t mean, Double_t sigma, Double_t leftCut, Double_t rightCut)
{
  // return number generated according to a gaussian distribution N(mean,sigma)
  // truncated at leftCut and rightCut
  //
  Double_t value;
  do{
    value=gRandom->Gaus(mean,sigma);
  }while((value-mean)<-leftCut || (value-mean)>rightCut);
  return value;
}

Double_t AliMathBase::BetheBlochAleph(Double_t bg,
         Double_t kp1,
         Double_t kp2,
         Double_t kp3,
         Double_t kp4,
         Double_t kp5) {
  //
  // This is the empirical ALEPH parameterization of the Bethe-Bloch formula.
  // It is normalized to 1 at the minimum.
  //
  // bg - beta*gamma
  // 
  // The default values for the kp* parameters are for ALICE TPC.
  // The returned value is in MIP units
  //

  return AliExternalTrackParam::BetheBlochAleph(bg,kp1,kp2,kp3,kp4,kp5);
}

Double_t AliMathBase::Gamma(Double_t k)
{
  // from
  // Hisashi Tanizaki
  // http://citeseerx.ist.psu.edu/viewdoc/download?doi=10.1.1.158.3866&rep=rep1&type=pdf
  // A. Morsch 14/01/2014
  static Double_t n=0;
  static Double_t c1=0;
  static Double_t c2=0;
  static Double_t b1=0;
  static Double_t b2=0;
  if (k > 0) {
    if (k < 0.4) 
      n = 1./k;
    else if (k >= 0.4 && k < 4) 
      n = 1./k + (k - 0.4)/k/3.6;
    else if (k >= 4.) 
      n = 1./TMath::Sqrt(k);
    b1 = k - 1./n;
    b2 = k + 1./n;
    c1 = (k < 0.4)? 0 : b1 * (TMath::Log(b1) - 1.)/2.;
    c2 = b2 * (TMath::Log(b2) - 1.)/2.;
  }
  Double_t x;
  Double_t y = -1.;
  while (1) {
    Double_t nu1 = gRandom->Rndm();
    Double_t nu2 = gRandom->Rndm();
    Double_t w1 = c1 + TMath::Log(nu1);
    Double_t w2 = c2 + TMath::Log(nu2);
    y = n * (b1 * w2 - b2 * w1);
    if (y < 0) continue;
    x = n * (w2 - w1);
    if (TMath::Log(y) >= x) break;
  }
  return TMath::Exp(x);
}
 AliMathBase.cxx:1
 AliMathBase.cxx:2
 AliMathBase.cxx:3
 AliMathBase.cxx:4
 AliMathBase.cxx:5
 AliMathBase.cxx:6
 AliMathBase.cxx:7
 AliMathBase.cxx:8
 AliMathBase.cxx:9
 AliMathBase.cxx:10
 AliMathBase.cxx:11
 AliMathBase.cxx:12
 AliMathBase.cxx:13
 AliMathBase.cxx:14
 AliMathBase.cxx:15
 AliMathBase.cxx:16
 AliMathBase.cxx:17
 AliMathBase.cxx:18
 AliMathBase.cxx:19
 AliMathBase.cxx:20
 AliMathBase.cxx:21
 AliMathBase.cxx:22
 AliMathBase.cxx:23
 AliMathBase.cxx:24
 AliMathBase.cxx:25
 AliMathBase.cxx:26
 AliMathBase.cxx:27
 AliMathBase.cxx:28
 AliMathBase.cxx:29
 AliMathBase.cxx:30
 AliMathBase.cxx:31
 AliMathBase.cxx:32
 AliMathBase.cxx:33
 AliMathBase.cxx:34
 AliMathBase.cxx:35
 AliMathBase.cxx:36
 AliMathBase.cxx:37
 AliMathBase.cxx:38
 AliMathBase.cxx:39
 AliMathBase.cxx:40
 AliMathBase.cxx:41
 AliMathBase.cxx:42
 AliMathBase.cxx:43
 AliMathBase.cxx:44
 AliMathBase.cxx:45
 AliMathBase.cxx:46
 AliMathBase.cxx:47
 AliMathBase.cxx:48
 AliMathBase.cxx:49
 AliMathBase.cxx:50
 AliMathBase.cxx:51
 AliMathBase.cxx:52
 AliMathBase.cxx:53
 AliMathBase.cxx:54
 AliMathBase.cxx:55
 AliMathBase.cxx:56
 AliMathBase.cxx:57
 AliMathBase.cxx:58
 AliMathBase.cxx:59
 AliMathBase.cxx:60
 AliMathBase.cxx:61
 AliMathBase.cxx:62
 AliMathBase.cxx:63
 AliMathBase.cxx:64
 AliMathBase.cxx:65
 AliMathBase.cxx:66
 AliMathBase.cxx:67
 AliMathBase.cxx:68
 AliMathBase.cxx:69
 AliMathBase.cxx:70
 AliMathBase.cxx:71
 AliMathBase.cxx:72
 AliMathBase.cxx:73
 AliMathBase.cxx:74
 AliMathBase.cxx:75
 AliMathBase.cxx:76
 AliMathBase.cxx:77
 AliMathBase.cxx:78
 AliMathBase.cxx:79
 AliMathBase.cxx:80
 AliMathBase.cxx:81
 AliMathBase.cxx:82
 AliMathBase.cxx:83
 AliMathBase.cxx:84
 AliMathBase.cxx:85
 AliMathBase.cxx:86
 AliMathBase.cxx:87
 AliMathBase.cxx:88
 AliMathBase.cxx:89
 AliMathBase.cxx:90
 AliMathBase.cxx:91
 AliMathBase.cxx:92
 AliMathBase.cxx:93
 AliMathBase.cxx:94
 AliMathBase.cxx:95
 AliMathBase.cxx:96
 AliMathBase.cxx:97
 AliMathBase.cxx:98
 AliMathBase.cxx:99
 AliMathBase.cxx:100
 AliMathBase.cxx:101
 AliMathBase.cxx:102
 AliMathBase.cxx:103
 AliMathBase.cxx:104
 AliMathBase.cxx:105
 AliMathBase.cxx:106
 AliMathBase.cxx:107
 AliMathBase.cxx:108
 AliMathBase.cxx:109
 AliMathBase.cxx:110
 AliMathBase.cxx:111
 AliMathBase.cxx:112
 AliMathBase.cxx:113
 AliMathBase.cxx:114
 AliMathBase.cxx:115
 AliMathBase.cxx:116
 AliMathBase.cxx:117
 AliMathBase.cxx:118
 AliMathBase.cxx:119
 AliMathBase.cxx:120
 AliMathBase.cxx:121
 AliMathBase.cxx:122
 AliMathBase.cxx:123
 AliMathBase.cxx:124
 AliMathBase.cxx:125
 AliMathBase.cxx:126
 AliMathBase.cxx:127
 AliMathBase.cxx:128
 AliMathBase.cxx:129
 AliMathBase.cxx:130
 AliMathBase.cxx:131
 AliMathBase.cxx:132
 AliMathBase.cxx:133
 AliMathBase.cxx:134
 AliMathBase.cxx:135
 AliMathBase.cxx:136
 AliMathBase.cxx:137
 AliMathBase.cxx:138
 AliMathBase.cxx:139
 AliMathBase.cxx:140
 AliMathBase.cxx:141
 AliMathBase.cxx:142
 AliMathBase.cxx:143
 AliMathBase.cxx:144
 AliMathBase.cxx:145
 AliMathBase.cxx:146
 AliMathBase.cxx:147
 AliMathBase.cxx:148
 AliMathBase.cxx:149
 AliMathBase.cxx:150
 AliMathBase.cxx:151
 AliMathBase.cxx:152
 AliMathBase.cxx:153
 AliMathBase.cxx:154
 AliMathBase.cxx:155
 AliMathBase.cxx:156
 AliMathBase.cxx:157
 AliMathBase.cxx:158
 AliMathBase.cxx:159
 AliMathBase.cxx:160
 AliMathBase.cxx:161
 AliMathBase.cxx:162
 AliMathBase.cxx:163
 AliMathBase.cxx:164
 AliMathBase.cxx:165
 AliMathBase.cxx:166
 AliMathBase.cxx:167
 AliMathBase.cxx:168
 AliMathBase.cxx:169
 AliMathBase.cxx:170
 AliMathBase.cxx:171
 AliMathBase.cxx:172
 AliMathBase.cxx:173
 AliMathBase.cxx:174
 AliMathBase.cxx:175
 AliMathBase.cxx:176
 AliMathBase.cxx:177
 AliMathBase.cxx:178
 AliMathBase.cxx:179
 AliMathBase.cxx:180
 AliMathBase.cxx:181
 AliMathBase.cxx:182
 AliMathBase.cxx:183
 AliMathBase.cxx:184
 AliMathBase.cxx:185
 AliMathBase.cxx:186
 AliMathBase.cxx:187
 AliMathBase.cxx:188
 AliMathBase.cxx:189
 AliMathBase.cxx:190
 AliMathBase.cxx:191
 AliMathBase.cxx:192
 AliMathBase.cxx:193
 AliMathBase.cxx:194
 AliMathBase.cxx:195
 AliMathBase.cxx:196
 AliMathBase.cxx:197
 AliMathBase.cxx:198
 AliMathBase.cxx:199
 AliMathBase.cxx:200
 AliMathBase.cxx:201
 AliMathBase.cxx:202
 AliMathBase.cxx:203
 AliMathBase.cxx:204
 AliMathBase.cxx:205
 AliMathBase.cxx:206
 AliMathBase.cxx:207
 AliMathBase.cxx:208
 AliMathBase.cxx:209
 AliMathBase.cxx:210
 AliMathBase.cxx:211
 AliMathBase.cxx:212
 AliMathBase.cxx:213
 AliMathBase.cxx:214
 AliMathBase.cxx:215
 AliMathBase.cxx:216
 AliMathBase.cxx:217
 AliMathBase.cxx:218
 AliMathBase.cxx:219
 AliMathBase.cxx:220
 AliMathBase.cxx:221
 AliMathBase.cxx:222
 AliMathBase.cxx:223
 AliMathBase.cxx:224
 AliMathBase.cxx:225
 AliMathBase.cxx:226
 AliMathBase.cxx:227
 AliMathBase.cxx:228
 AliMathBase.cxx:229
 AliMathBase.cxx:230
 AliMathBase.cxx:231
 AliMathBase.cxx:232
 AliMathBase.cxx:233
 AliMathBase.cxx:234
 AliMathBase.cxx:235
 AliMathBase.cxx:236
 AliMathBase.cxx:237
 AliMathBase.cxx:238
 AliMathBase.cxx:239
 AliMathBase.cxx:240
 AliMathBase.cxx:241
 AliMathBase.cxx:242
 AliMathBase.cxx:243
 AliMathBase.cxx:244
 AliMathBase.cxx:245
 AliMathBase.cxx:246
 AliMathBase.cxx:247
 AliMathBase.cxx:248
 AliMathBase.cxx:249
 AliMathBase.cxx:250
 AliMathBase.cxx:251
 AliMathBase.cxx:252
 AliMathBase.cxx:253
 AliMathBase.cxx:254
 AliMathBase.cxx:255
 AliMathBase.cxx:256
 AliMathBase.cxx:257
 AliMathBase.cxx:258
 AliMathBase.cxx:259
 AliMathBase.cxx:260
 AliMathBase.cxx:261
 AliMathBase.cxx:262
 AliMathBase.cxx:263
 AliMathBase.cxx:264
 AliMathBase.cxx:265
 AliMathBase.cxx:266
 AliMathBase.cxx:267
 AliMathBase.cxx:268
 AliMathBase.cxx:269
 AliMathBase.cxx:270
 AliMathBase.cxx:271
 AliMathBase.cxx:272
 AliMathBase.cxx:273
 AliMathBase.cxx:274
 AliMathBase.cxx:275
 AliMathBase.cxx:276
 AliMathBase.cxx:277
 AliMathBase.cxx:278
 AliMathBase.cxx:279
 AliMathBase.cxx:280
 AliMathBase.cxx:281
 AliMathBase.cxx:282
 AliMathBase.cxx:283
 AliMathBase.cxx:284
 AliMathBase.cxx:285
 AliMathBase.cxx:286
 AliMathBase.cxx:287
 AliMathBase.cxx:288
 AliMathBase.cxx:289
 AliMathBase.cxx:290
 AliMathBase.cxx:291
 AliMathBase.cxx:292
 AliMathBase.cxx:293
 AliMathBase.cxx:294
 AliMathBase.cxx:295
 AliMathBase.cxx:296
 AliMathBase.cxx:297
 AliMathBase.cxx:298
 AliMathBase.cxx:299
 AliMathBase.cxx:300
 AliMathBase.cxx:301
 AliMathBase.cxx:302
 AliMathBase.cxx:303
 AliMathBase.cxx:304
 AliMathBase.cxx:305
 AliMathBase.cxx:306
 AliMathBase.cxx:307
 AliMathBase.cxx:308
 AliMathBase.cxx:309
 AliMathBase.cxx:310
 AliMathBase.cxx:311
 AliMathBase.cxx:312
 AliMathBase.cxx:313
 AliMathBase.cxx:314
 AliMathBase.cxx:315
 AliMathBase.cxx:316
 AliMathBase.cxx:317
 AliMathBase.cxx:318
 AliMathBase.cxx:319
 AliMathBase.cxx:320
 AliMathBase.cxx:321
 AliMathBase.cxx:322
 AliMathBase.cxx:323
 AliMathBase.cxx:324
 AliMathBase.cxx:325
 AliMathBase.cxx:326
 AliMathBase.cxx:327
 AliMathBase.cxx:328
 AliMathBase.cxx:329
 AliMathBase.cxx:330
 AliMathBase.cxx:331
 AliMathBase.cxx:332
 AliMathBase.cxx:333
 AliMathBase.cxx:334
 AliMathBase.cxx:335
 AliMathBase.cxx:336
 AliMathBase.cxx:337
 AliMathBase.cxx:338
 AliMathBase.cxx:339
 AliMathBase.cxx:340
 AliMathBase.cxx:341
 AliMathBase.cxx:342
 AliMathBase.cxx:343
 AliMathBase.cxx:344
 AliMathBase.cxx:345
 AliMathBase.cxx:346
 AliMathBase.cxx:347
 AliMathBase.cxx:348
 AliMathBase.cxx:349
 AliMathBase.cxx:350
 AliMathBase.cxx:351
 AliMathBase.cxx:352
 AliMathBase.cxx:353
 AliMathBase.cxx:354
 AliMathBase.cxx:355
 AliMathBase.cxx:356
 AliMathBase.cxx:357
 AliMathBase.cxx:358
 AliMathBase.cxx:359
 AliMathBase.cxx:360
 AliMathBase.cxx:361
 AliMathBase.cxx:362
 AliMathBase.cxx:363
 AliMathBase.cxx:364
 AliMathBase.cxx:365
 AliMathBase.cxx:366
 AliMathBase.cxx:367
 AliMathBase.cxx:368
 AliMathBase.cxx:369
 AliMathBase.cxx:370
 AliMathBase.cxx:371
 AliMathBase.cxx:372
 AliMathBase.cxx:373
 AliMathBase.cxx:374
 AliMathBase.cxx:375
 AliMathBase.cxx:376
 AliMathBase.cxx:377
 AliMathBase.cxx:378
 AliMathBase.cxx:379
 AliMathBase.cxx:380
 AliMathBase.cxx:381
 AliMathBase.cxx:382
 AliMathBase.cxx:383
 AliMathBase.cxx:384
 AliMathBase.cxx:385
 AliMathBase.cxx:386
 AliMathBase.cxx:387
 AliMathBase.cxx:388
 AliMathBase.cxx:389
 AliMathBase.cxx:390
 AliMathBase.cxx:391
 AliMathBase.cxx:392
 AliMathBase.cxx:393
 AliMathBase.cxx:394
 AliMathBase.cxx:395
 AliMathBase.cxx:396
 AliMathBase.cxx:397
 AliMathBase.cxx:398
 AliMathBase.cxx:399
 AliMathBase.cxx:400
 AliMathBase.cxx:401
 AliMathBase.cxx:402
 AliMathBase.cxx:403
 AliMathBase.cxx:404
 AliMathBase.cxx:405
 AliMathBase.cxx:406
 AliMathBase.cxx:407
 AliMathBase.cxx:408
 AliMathBase.cxx:409
 AliMathBase.cxx:410
 AliMathBase.cxx:411
 AliMathBase.cxx:412
 AliMathBase.cxx:413
 AliMathBase.cxx:414
 AliMathBase.cxx:415
 AliMathBase.cxx:416
 AliMathBase.cxx:417
 AliMathBase.cxx:418
 AliMathBase.cxx:419
 AliMathBase.cxx:420
 AliMathBase.cxx:421
 AliMathBase.cxx:422
 AliMathBase.cxx:423
 AliMathBase.cxx:424
 AliMathBase.cxx:425
 AliMathBase.cxx:426
 AliMathBase.cxx:427
 AliMathBase.cxx:428
 AliMathBase.cxx:429
 AliMathBase.cxx:430
 AliMathBase.cxx:431
 AliMathBase.cxx:432
 AliMathBase.cxx:433
 AliMathBase.cxx:434
 AliMathBase.cxx:435
 AliMathBase.cxx:436
 AliMathBase.cxx:437
 AliMathBase.cxx:438
 AliMathBase.cxx:439
 AliMathBase.cxx:440
 AliMathBase.cxx:441
 AliMathBase.cxx:442
 AliMathBase.cxx:443
 AliMathBase.cxx:444
 AliMathBase.cxx:445
 AliMathBase.cxx:446
 AliMathBase.cxx:447
 AliMathBase.cxx:448
 AliMathBase.cxx:449
 AliMathBase.cxx:450
 AliMathBase.cxx:451
 AliMathBase.cxx:452
 AliMathBase.cxx:453
 AliMathBase.cxx:454
 AliMathBase.cxx:455
 AliMathBase.cxx:456
 AliMathBase.cxx:457
 AliMathBase.cxx:458
 AliMathBase.cxx:459
 AliMathBase.cxx:460
 AliMathBase.cxx:461
 AliMathBase.cxx:462
 AliMathBase.cxx:463
 AliMathBase.cxx:464
 AliMathBase.cxx:465
 AliMathBase.cxx:466
 AliMathBase.cxx:467
 AliMathBase.cxx:468
 AliMathBase.cxx:469
 AliMathBase.cxx:470
 AliMathBase.cxx:471
 AliMathBase.cxx:472
 AliMathBase.cxx:473
 AliMathBase.cxx:474
 AliMathBase.cxx:475
 AliMathBase.cxx:476
 AliMathBase.cxx:477
 AliMathBase.cxx:478
 AliMathBase.cxx:479
 AliMathBase.cxx:480
 AliMathBase.cxx:481
 AliMathBase.cxx:482
 AliMathBase.cxx:483
 AliMathBase.cxx:484
 AliMathBase.cxx:485
 AliMathBase.cxx:486
 AliMathBase.cxx:487
 AliMathBase.cxx:488
 AliMathBase.cxx:489
 AliMathBase.cxx:490
 AliMathBase.cxx:491
 AliMathBase.cxx:492
 AliMathBase.cxx:493
 AliMathBase.cxx:494
 AliMathBase.cxx:495
 AliMathBase.cxx:496
 AliMathBase.cxx:497
 AliMathBase.cxx:498
 AliMathBase.cxx:499
 AliMathBase.cxx:500
 AliMathBase.cxx:501
 AliMathBase.cxx:502
 AliMathBase.cxx:503
 AliMathBase.cxx:504
 AliMathBase.cxx:505
 AliMathBase.cxx:506
 AliMathBase.cxx:507
 AliMathBase.cxx:508
 AliMathBase.cxx:509
 AliMathBase.cxx:510
 AliMathBase.cxx:511
 AliMathBase.cxx:512
 AliMathBase.cxx:513
 AliMathBase.cxx:514
 AliMathBase.cxx:515
 AliMathBase.cxx:516
 AliMathBase.cxx:517
 AliMathBase.cxx:518
 AliMathBase.cxx:519
 AliMathBase.cxx:520
 AliMathBase.cxx:521
 AliMathBase.cxx:522
 AliMathBase.cxx:523
 AliMathBase.cxx:524
 AliMathBase.cxx:525
 AliMathBase.cxx:526
 AliMathBase.cxx:527
 AliMathBase.cxx:528
 AliMathBase.cxx:529
 AliMathBase.cxx:530
 AliMathBase.cxx:531
 AliMathBase.cxx:532
 AliMathBase.cxx:533
 AliMathBase.cxx:534
 AliMathBase.cxx:535
 AliMathBase.cxx:536
 AliMathBase.cxx:537
 AliMathBase.cxx:538
 AliMathBase.cxx:539
 AliMathBase.cxx:540
 AliMathBase.cxx:541
 AliMathBase.cxx:542
 AliMathBase.cxx:543
 AliMathBase.cxx:544
 AliMathBase.cxx:545
 AliMathBase.cxx:546
 AliMathBase.cxx:547
 AliMathBase.cxx:548
 AliMathBase.cxx:549
 AliMathBase.cxx:550
 AliMathBase.cxx:551
 AliMathBase.cxx:552
 AliMathBase.cxx:553
 AliMathBase.cxx:554
 AliMathBase.cxx:555
 AliMathBase.cxx:556
 AliMathBase.cxx:557
 AliMathBase.cxx:558
 AliMathBase.cxx:559
 AliMathBase.cxx:560
 AliMathBase.cxx:561
 AliMathBase.cxx:562
 AliMathBase.cxx:563
 AliMathBase.cxx:564
 AliMathBase.cxx:565
 AliMathBase.cxx:566
 AliMathBase.cxx:567
 AliMathBase.cxx:568
 AliMathBase.cxx:569
 AliMathBase.cxx:570
 AliMathBase.cxx:571
 AliMathBase.cxx:572
 AliMathBase.cxx:573
 AliMathBase.cxx:574
 AliMathBase.cxx:575
 AliMathBase.cxx:576
 AliMathBase.cxx:577
 AliMathBase.cxx:578
 AliMathBase.cxx:579
 AliMathBase.cxx:580
 AliMathBase.cxx:581
 AliMathBase.cxx:582
 AliMathBase.cxx:583
 AliMathBase.cxx:584
 AliMathBase.cxx:585
 AliMathBase.cxx:586
 AliMathBase.cxx:587
 AliMathBase.cxx:588
 AliMathBase.cxx:589
 AliMathBase.cxx:590
 AliMathBase.cxx:591
 AliMathBase.cxx:592
 AliMathBase.cxx:593
 AliMathBase.cxx:594
 AliMathBase.cxx:595
 AliMathBase.cxx:596
 AliMathBase.cxx:597
 AliMathBase.cxx:598
 AliMathBase.cxx:599
 AliMathBase.cxx:600
 AliMathBase.cxx:601
 AliMathBase.cxx:602
 AliMathBase.cxx:603
 AliMathBase.cxx:604
 AliMathBase.cxx:605
 AliMathBase.cxx:606
 AliMathBase.cxx:607
 AliMathBase.cxx:608
 AliMathBase.cxx:609
 AliMathBase.cxx:610
 AliMathBase.cxx:611
 AliMathBase.cxx:612
 AliMathBase.cxx:613
 AliMathBase.cxx:614
 AliMathBase.cxx:615
 AliMathBase.cxx:616
 AliMathBase.cxx:617
 AliMathBase.cxx:618
 AliMathBase.cxx:619
 AliMathBase.cxx:620
 AliMathBase.cxx:621
 AliMathBase.cxx:622
 AliMathBase.cxx:623
 AliMathBase.cxx:624
 AliMathBase.cxx:625
 AliMathBase.cxx:626
 AliMathBase.cxx:627
 AliMathBase.cxx:628
 AliMathBase.cxx:629
 AliMathBase.cxx:630
 AliMathBase.cxx:631
 AliMathBase.cxx:632
 AliMathBase.cxx:633
 AliMathBase.cxx:634
 AliMathBase.cxx:635
 AliMathBase.cxx:636
 AliMathBase.cxx:637
 AliMathBase.cxx:638
 AliMathBase.cxx:639
 AliMathBase.cxx:640
 AliMathBase.cxx:641
 AliMathBase.cxx:642
 AliMathBase.cxx:643
 AliMathBase.cxx:644
 AliMathBase.cxx:645
 AliMathBase.cxx:646
 AliMathBase.cxx:647
 AliMathBase.cxx:648
 AliMathBase.cxx:649
 AliMathBase.cxx:650
 AliMathBase.cxx:651
 AliMathBase.cxx:652
 AliMathBase.cxx:653
 AliMathBase.cxx:654
 AliMathBase.cxx:655
 AliMathBase.cxx:656
 AliMathBase.cxx:657
 AliMathBase.cxx:658
 AliMathBase.cxx:659
 AliMathBase.cxx:660
 AliMathBase.cxx:661
 AliMathBase.cxx:662
 AliMathBase.cxx:663
 AliMathBase.cxx:664
 AliMathBase.cxx:665
 AliMathBase.cxx:666
 AliMathBase.cxx:667
 AliMathBase.cxx:668
 AliMathBase.cxx:669
 AliMathBase.cxx:670
 AliMathBase.cxx:671
 AliMathBase.cxx:672
 AliMathBase.cxx:673
 AliMathBase.cxx:674
 AliMathBase.cxx:675
 AliMathBase.cxx:676
 AliMathBase.cxx:677
 AliMathBase.cxx:678
 AliMathBase.cxx:679
 AliMathBase.cxx:680
 AliMathBase.cxx:681
 AliMathBase.cxx:682
 AliMathBase.cxx:683
 AliMathBase.cxx:684
 AliMathBase.cxx:685
 AliMathBase.cxx:686
 AliMathBase.cxx:687
 AliMathBase.cxx:688
 AliMathBase.cxx:689
 AliMathBase.cxx:690
 AliMathBase.cxx:691
 AliMathBase.cxx:692
 AliMathBase.cxx:693
 AliMathBase.cxx:694
 AliMathBase.cxx:695
 AliMathBase.cxx:696
 AliMathBase.cxx:697
 AliMathBase.cxx:698
 AliMathBase.cxx:699
 AliMathBase.cxx:700
 AliMathBase.cxx:701
 AliMathBase.cxx:702
 AliMathBase.cxx:703
 AliMathBase.cxx:704
 AliMathBase.cxx:705
 AliMathBase.cxx:706
 AliMathBase.cxx:707
 AliMathBase.cxx:708
 AliMathBase.cxx:709
 AliMathBase.cxx:710
 AliMathBase.cxx:711
 AliMathBase.cxx:712
 AliMathBase.cxx:713
 AliMathBase.cxx:714
 AliMathBase.cxx:715
 AliMathBase.cxx:716
 AliMathBase.cxx:717
 AliMathBase.cxx:718
 AliMathBase.cxx:719
 AliMathBase.cxx:720
 AliMathBase.cxx:721
 AliMathBase.cxx:722
 AliMathBase.cxx:723
 AliMathBase.cxx:724
 AliMathBase.cxx:725
 AliMathBase.cxx:726
 AliMathBase.cxx:727
 AliMathBase.cxx:728
 AliMathBase.cxx:729
 AliMathBase.cxx:730
 AliMathBase.cxx:731
 AliMathBase.cxx:732
 AliMathBase.cxx:733
 AliMathBase.cxx:734
 AliMathBase.cxx:735
 AliMathBase.cxx:736
 AliMathBase.cxx:737
 AliMathBase.cxx:738
 AliMathBase.cxx:739
 AliMathBase.cxx:740
 AliMathBase.cxx:741
 AliMathBase.cxx:742
 AliMathBase.cxx:743
 AliMathBase.cxx:744
 AliMathBase.cxx:745
 AliMathBase.cxx:746
 AliMathBase.cxx:747
 AliMathBase.cxx:748
 AliMathBase.cxx:749
 AliMathBase.cxx:750
 AliMathBase.cxx:751
 AliMathBase.cxx:752
 AliMathBase.cxx:753
 AliMathBase.cxx:754
 AliMathBase.cxx:755
 AliMathBase.cxx:756
 AliMathBase.cxx:757
 AliMathBase.cxx:758
 AliMathBase.cxx:759
 AliMathBase.cxx:760
 AliMathBase.cxx:761
 AliMathBase.cxx:762
 AliMathBase.cxx:763
 AliMathBase.cxx:764
 AliMathBase.cxx:765
 AliMathBase.cxx:766
 AliMathBase.cxx:767
 AliMathBase.cxx:768
 AliMathBase.cxx:769
 AliMathBase.cxx:770
 AliMathBase.cxx:771
 AliMathBase.cxx:772
 AliMathBase.cxx:773
 AliMathBase.cxx:774
 AliMathBase.cxx:775
 AliMathBase.cxx:776
 AliMathBase.cxx:777
 AliMathBase.cxx:778
 AliMathBase.cxx:779
 AliMathBase.cxx:780
 AliMathBase.cxx:781
 AliMathBase.cxx:782
 AliMathBase.cxx:783
 AliMathBase.cxx:784
 AliMathBase.cxx:785
 AliMathBase.cxx:786
 AliMathBase.cxx:787
 AliMathBase.cxx:788
 AliMathBase.cxx:789
 AliMathBase.cxx:790
 AliMathBase.cxx:791
 AliMathBase.cxx:792
 AliMathBase.cxx:793
 AliMathBase.cxx:794
 AliMathBase.cxx:795
 AliMathBase.cxx:796
 AliMathBase.cxx:797
 AliMathBase.cxx:798
 AliMathBase.cxx:799
 AliMathBase.cxx:800
 AliMathBase.cxx:801
 AliMathBase.cxx:802
 AliMathBase.cxx:803
 AliMathBase.cxx:804
 AliMathBase.cxx:805
 AliMathBase.cxx:806
 AliMathBase.cxx:807
 AliMathBase.cxx:808
 AliMathBase.cxx:809
 AliMathBase.cxx:810
 AliMathBase.cxx:811
 AliMathBase.cxx:812
 AliMathBase.cxx:813
 AliMathBase.cxx:814
 AliMathBase.cxx:815
 AliMathBase.cxx:816
 AliMathBase.cxx:817
 AliMathBase.cxx:818
 AliMathBase.cxx:819
 AliMathBase.cxx:820
 AliMathBase.cxx:821
 AliMathBase.cxx:822
 AliMathBase.cxx:823
 AliMathBase.cxx:824
 AliMathBase.cxx:825
 AliMathBase.cxx:826
 AliMathBase.cxx:827
 AliMathBase.cxx:828
 AliMathBase.cxx:829
 AliMathBase.cxx:830
 AliMathBase.cxx:831
 AliMathBase.cxx:832
 AliMathBase.cxx:833
 AliMathBase.cxx:834