ROOT logo
// ----------------------------------------------------------------------
//                     AliPWGHistoTools
// 
// This class provides some tools which can be useful in the analsis
// of spectra, to fit or transform histograms. See the comments of the
// individual methods for details
//
// Author: M. Floris (CERN)
// ----------------------------------------------------------------------

#include "AliPWGHistoTools.h"
#include "TH1D.h"
#include "TF1.h"
#include "TH1.h"
#include "TMath.h"
#include "TGraphErrors.h"
#include "TVirtualFitter.h"
#include "TMinuit.h"
#include "AliLog.h"
#include "TSpline.h"
#include <iostream>
#include "TGraphAsymmErrors.h"

using namespace std;


ClassImp(AliPWGHistoTools)

TF1 * AliPWGHistoTools::fdNdptForETCalc = 0;

AliPWGHistoTools::AliPWGHistoTools() {
  // ctor
}

AliPWGHistoTools::~AliPWGHistoTools(){
  // dtor
}

TH1 * AliPWGHistoTools::GetdNdmtFromdNdpt(const TH1 * hpt, Double_t mass) {
  // convert the x axis from pt to mt. Assumes you have 1/pt dNdpt in the histo you start with

  Int_t nbins = hpt->GetNbinsX();
  Float_t * xbins = new Float_t[nbins+1];
  for(Int_t ibins = 0; ibins <= nbins; ibins++){
    xbins[ibins] = TMath::Sqrt(hpt->GetBinLowEdge(ibins+1)*hpt->GetBinLowEdge(ibins+1) +
    			       mass *mass) - mass;
    // // xbins[ibins] = TMath::Sqrt(hpt->GetBinLowEdge(ibins+1)*hpt->GetBinLowEdge(ibins+1) +
    // // 			       mass *mass);
    //    cout << ibins << " "<< xbins[ibins]  << endl;

  }

  TH1D * hmt = new TH1D(TString(hpt->GetName())+"_mt",
		      TString(hpt->GetName())+"_mt",
		      nbins, xbins);
  for(Int_t ibins = 1; ibins <= nbins; ibins++){
    hmt->SetBinContent(ibins, hpt->GetBinContent(ibins));
    hmt->SetBinError(ibins,   hpt->GetBinError(ibins));

  }

  hmt->SetXTitle("m_{t} - m_{0} (GeV/c^{2})");
  hmt->SetYTitle("1/m_{t} dN/dm_{t} (a.u.)");
  hmt->SetMarkerStyle(hpt->GetMarkerStyle());
  hmt->SetMarkerColor(hpt->GetMarkerColor());
  hmt->SetLineColor(hpt->GetLineColor());
  
  return hmt;

}

TH1 * AliPWGHistoTools::GetdNdptFromdNdmt(const TH1 * hmt, Double_t mass) {
  // convert the x axis from mt to pt. Assumes you have 1/mt dNdmt in the histo you start with

  Int_t nbins = hmt->GetNbinsX();
  Float_t * xbins = new Float_t[nbins+1];
  for(Int_t ibins = 0; ibins <= nbins; ibins++){
    xbins[ibins] = TMath::Sqrt((hmt->GetBinLowEdge(ibins+1)+mass)*(hmt->GetBinLowEdge(ibins+1)+mass) -
    			       mass *mass);
    xbins[ibins] = Float_t(TMath::Nint(xbins[ibins]*100))/100;
    // // xbins[ibins] = TMath::Sqrt(hmt->GetBinLowEdge(ibins+1)*hmt->GetBinLowEdge(ibins+1) +
    // // 			       mass *mass);
    cout << ibins << " "<< xbins[ibins]  << endl;

  }

  TH1D * hptL = new TH1D(TString(hmt->GetName())+"_pt",
		      TString(hmt->GetName())+"_pt",
		      nbins, xbins);
  for(Int_t ibins = 1; ibins <= nbins; ibins++){
    hptL->SetBinContent(ibins, hmt->GetBinContent(ibins));
    hptL->SetBinError(ibins,   hmt->GetBinError(ibins));

  }

  hptL->SetXTitle("p_{t} (GeV/c)");
  hptL->SetYTitle("1/p_{t} dN/dp_{t} (a.u.)");
  hptL->SetMarkerStyle(hmt->GetMarkerStyle());
  hptL->SetMarkerColor(hmt->GetMarkerColor());
  hptL->SetLineColor(hmt->GetLineColor());

  return hptL;

}


TH1 * AliPWGHistoTools::GetdNdPtFromOneOverPt(const TH1 * h1Pt) {

  // convert an histo from 1/pt dNdpt to dNdpt, using the pt at the center of the bin


  TH1 * hPt = (TH1 *) h1Pt->Clone((TString(h1Pt->GetName()) + "_inv").Data());
  hPt->Reset();

  Int_t nbinx = hPt->GetNbinsX();

  for(Int_t ibinx = 1; ibinx <= nbinx; ibinx++){

    Double_t cont = h1Pt->GetBinContent(ibinx);
    Double_t err  = h1Pt->GetBinError(ibinx);
    
    Double_t pt   = h1Pt->GetBinCenter(ibinx);
    
    if(pt > 0) {
      cont *= pt;
      err  *= pt;
    } else {
      cont = 0;
      err  = 0;
    }

    hPt->SetBinContent(ibinx, cont);
    hPt->SetBinError(ibinx, err);
    
  }

  hPt->SetXTitle("p_{t} (GeV)");
  hPt->SetYTitle("dN/dp_{t} (GeV^{-2})");

  return hPt;    

}




TH1 * AliPWGHistoTools::GetOneOverPtdNdPt(const TH1 * hPt) {

  // convert an histo from dNdpt to 1/pt dNdpt, using the pt at the center of the bin

  TH1 * h1Pt = (TH1 *) hPt->Clone((TString(hPt->GetName()) + "_inv").Data());
  h1Pt->Reset();

  Int_t nbinx = h1Pt->GetNbinsX();

  for(Int_t ibinx = 1; ibinx <= nbinx; ibinx++){

    Double_t cont = hPt->GetBinContent(ibinx);
    Double_t err  = hPt->GetBinError(ibinx);
    
    Double_t pt   = hPt->GetBinCenter(ibinx);
    
    if(pt > 0) {
      cont /= pt;
      err  /= pt;
    } else {
      cont = 0;
      err  = 0;
    }

    h1Pt->SetBinContent(ibinx, cont);
    h1Pt->SetBinError(ibinx, err);
    
  }

  h1Pt->SetXTitle("p_{t} (GeV)");
  h1Pt->SetYTitle("1/p_{t} dN/dp_{t} (GeV^{-2})");

  return h1Pt;    

}


TGraphErrors * AliPWGHistoTools::GetGraphFromHisto(const TH1F * h, Bool_t binWidth) {
  // Convert a histo to a graph
  // if binWidth is true ex is set to the bin width of the histos, otherwise it is set to zero
  Int_t nbin = h->GetNbinsX();

  TGraphErrors * g = new TGraphErrors();
  Int_t ipoint = 0;
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    Double_t xerr = binWidth ? h->GetBinWidth(ibin)/2 : 0;
    if (h->GetBinContent(ibin)) {
      g->SetPoint     (ipoint,   h->GetBinCenter(ibin), h->GetBinContent(ibin));
      g->SetPointError(ipoint,   xerr,  h->GetBinError(ibin));
      ipoint++;
    }
  }
  
  g->SetMarkerStyle(h->GetMarkerStyle());
  g->SetMarkerColor(h->GetMarkerColor());
  g->SetLineColor(h->GetLineColor());
  g->SetLineStyle(h->GetLineStyle());
  g->SetLineWidth(h->GetLineWidth());

  g->SetTitle(h->GetTitle());
  g->SetName(TString("g_")+h->GetName());

  return g;

}

TH1F * AliPWGHistoTools::GetHistoFromGraph(const TGraphErrors * g, const TH1F* hTemplate) {

  // convert a graph to histo with the binning of hTemplate.
  // Warning: the template should be chosen
  // properly: if you have two graph points in the same histo bin this
  // won't work!

  TH1F * h = (TH1F*) hTemplate->Clone(TString("h_")+g->GetName());
  h->Reset();
  Int_t npoint = g->GetN();
  //  g->Print();
  for(Int_t ipoint = 0; ipoint < npoint; ipoint++){
    Float_t x  = g->GetX() [ipoint];
    Float_t y  = g->GetY() [ipoint];
    Float_t ey = 0;
    if(g->InheritsFrom("TGraphErrors"))
       ey = g->GetEY()[ipoint];
    Int_t bin = h->FindBin(x);
    //    cout << "bin: "<< bin << " -> " << x << ", "<< y <<", " << ey << endl;
    
    h->SetBinContent(bin,y);
    h->SetBinError  (bin,ey);
  }
 
  h->SetMarkerStyle(g->GetMarkerStyle());
  h->SetMarkerColor(g->GetMarkerColor());
  h->SetLineColor  (g->GetLineColor());

 
  return h;
}

TGraphErrors * AliPWGHistoTools::ConcatenateGraphs(const TGraphErrors * g1,const TGraphErrors * g2){

  // concatenates two graphs

  Int_t npoint1=g1->GetN();
  Int_t npoint2=g2->GetN();

  TGraphErrors * gClone = (TGraphErrors*) g1->Clone();

//   for(Int_t ipoint = 0; ipoint < npoint1; ipoint++){
//     gClone->SetPointError(ipoint,0,g1->GetEY()[ipoint]);

//   }
  for(Int_t ipoint = 0; ipoint < npoint2; ipoint++){
    gClone->SetPoint(ipoint+npoint1,g2->GetX()[ipoint],g2->GetY()[ipoint]);
    gClone->SetPointError(ipoint+npoint1,g2->GetEX()[ipoint],g2->GetEY()[ipoint]);
    //    gClone->SetPointError(ipoint+npoint1,0,g2->GetEY()[ipoint]);
  }

  gClone->GetHistogram()->GetXaxis()->SetTimeDisplay(1);
  gClone->SetTitle(TString(gClone->GetTitle())+" + "+g2->GetTitle());
  gClone->SetName(TString(gClone->GetName())+"_"+g2->GetName());

  return gClone;
}


TH1F * AliPWGHistoTools::Combine3HistosWithErrors(const TH1 * h1,  const TH1 * h2,  const TH1* h3, 
					    TH1 * he1,  TH1 * he2,  TH1 * he3, 
					    const TH1* htemplate, Int_t statFrom, 
					    Float_t renorm1, Float_t renorm2, Float_t renorm3,
					    TH1 ** hSyst, Bool_t errorFromBinContent) {

  // Combines 3 histos (h1,h2,h3), weighting by the errors provided in
  // he1,he2,he3, supposed to be the independent systematic errors.
  // he1,he2,he3 are also assumed to have the same binning as h1,h2,h3
  // The combined histo must fit the template provided (no check is performed on this)
  // The histogram are supposed to come from the same (nearly) sample
  // of tracks. statFrom tells the stat error of which of the 3 is
  // going to be assigned to the combined
  // Optionally, it is possible to rescale any of the histograms.
  // if hSyst is give, the histo is filled with combined syst error vs pt
  // if errorFromBinContent is true, the weights are taken from the he* content rather than errors
  TH1F * hcomb = (TH1F*) htemplate->Clone(TString("hComb_")+h1->GetName()+"_"+h2->GetName()+h3->GetName());

  // TODO: I should have used an array for h*local...

  // Clone histos locally to rescale them
  TH1F * h1local = (TH1F*) h1->Clone();
  TH1F * h2local = (TH1F*) h2->Clone();
  TH1F * h3local = (TH1F*) h3->Clone();
  h1local->Scale(renorm1);
  h2local->Scale(renorm2);
  h3local->Scale(renorm3);

  const TH1 * hStatError = 0;
  if (statFrom == 0)      hStatError = h1; 
  else if (statFrom == 1) hStatError = h2; 
  else if (statFrom == 2) hStatError = h3; 
  else {
    Printf("AliPWGHistoTools::Combine3HistosWithErrors: wrong value of the statFrom parameter");
    return NULL;
  }
  Printf("AliPWGHistoTools::Combine3HistosWithErrors: improve error on combined");
  // Loop over all bins and take weighted mean of all points
  Int_t nBinComb = hcomb->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nBinComb; ibin++){
    Int_t ibin1 = h1local->FindBin(hcomb->GetBinCenter(ibin));
    Int_t ibin2 = h2local->FindBin(hcomb->GetBinCenter(ibin));
    Int_t ibin3 = h3local->FindBin(hcomb->GetBinCenter(ibin));
    Int_t ibinError = -1; // bin used to get stat error

    if (statFrom == 0)      ibinError = ibin1; 
    else if (statFrom == 1) ibinError = ibin2; 
    else if (statFrom == 2) ibinError = ibin3; 
    else Printf("AliPWGHistoTools::Combine3HistosWithErrors: wrong value of the statFrom parameter");


    Double_t y[3];
    Double_t ye[3];
    y[0]  = h1local->GetBinContent(ibin1);
    y[1]  = h2local->GetBinContent(ibin2);
    y[2]  = h3local->GetBinContent(ibin3);
    if (errorFromBinContent) {
      ye[0] = he1->GetBinContent(he1->FindBin(hcomb->GetBinCenter(ibin)));
      ye[1] = he2->GetBinContent(he2->FindBin(hcomb->GetBinCenter(ibin)));
      ye[2] = he3->GetBinContent(he3->FindBin(hcomb->GetBinCenter(ibin)));
    } else {
      ye[0] = he1->GetBinError(ibin1);
      ye[1] = he2->GetBinError(ibin2);
      ye[2] = he3->GetBinError(ibin3);
    }
    // Set error to 0 if content is 0 (means it was not filled)
    if(!h1local->GetBinContent(ibin1)) ye[0] = 0;
    if(!h2local->GetBinContent(ibin2)) ye[1] = 0;
    if(!h3local->GetBinContent(ibin3)) ye[2] = 0;
    
    // Compute weighted mean
    //    cout << "Bins:  "<< ibin1 << " " << ibin2 << " " << ibin3 << endl;    
    Double_t mean, err;
    WeightedMean(3,y,ye,mean,err);


    // Fill combined
    hcomb->SetBinContent(ibin,mean);
    Double_t statError = 0;
    if (hStatError->GetBinContent(ibinError)) {
      //      cout << "def" << endl;
      statError = hStatError->GetBinError(ibinError)/hStatError->GetBinContent(ibinError)*hcomb->GetBinContent(ibin);
    }
    else if (h1local->GetBinContent(ibin1)) {
      //      cout << "1" << endl;
      statError = h1local->GetBinError(ibin1)/h1local->GetBinContent(ibin1)*hcomb->GetBinContent(ibin);
    }
    else if (h2local->GetBinContent(ibin2)) {
      //      cout << "2" << endl;
      statError = h2local->GetBinError(ibin2)/h2local->GetBinContent(ibin2)*hcomb->GetBinContent(ibin);
    }
    else if (h3local->GetBinContent(ibin3)) {
      //      cout << "3" << endl;
      statError = h3local->GetBinError(ibin3)/h3local->GetBinContent(ibin3)*hcomb->GetBinContent(ibin);
    }
    hcomb->SetBinError  (ibin,statError);
    if(hSyst) (*hSyst)->SetBinContent(ibin,err);
    //    cout << "BIN " << ibin << " " << mean << " " << statError << endl;

  }

  hcomb->SetMarkerStyle(hStatError->GetMarkerStyle());
  hcomb->SetMarkerColor(hStatError->GetMarkerColor());
  hcomb->SetLineColor  (hStatError->GetLineColor());

  hcomb->SetXTitle(hStatError->GetXaxis()->GetTitle());
  hcomb->SetYTitle(hStatError->GetYaxis()->GetTitle());

  delete h1local;
  delete h2local;
  delete h3local;

  return hcomb;
}


void AliPWGHistoTools::GetMeanDataAndExtrapolation(const TH1 * hData, TF1 * fExtrapolation, Double_t &mean, Double_t &error, Float_t min, Float_t max){
  // Computes the mean of the combined data and extrapolation in a
  // given range, use data where they are available and the function
  // where data are not available
  // ERROR from DATA ONLY is returned in this version! 
  //
  Printf("AliPWGHistoTools::GetMeanDataAndExtrapolation: WARNING from data only");
  Float_t minData    = GetLowestNotEmptyBinEdge (hData);
  Int_t minDataBin   = GetLowestNotEmptyBin     (hData);
  Float_t maxData    = GetHighestNotEmptyBinEdge(hData);
  Int_t maxDataBin   = GetHighestNotEmptyBin    (hData);
  Double_t integral  = 0; 
  mean      = 0;
  error     = 0; 
  if (min < minData) {
    // Compute integral exploiting root function to calculate moments, "unnormalizing" them
    mean     += fExtrapolation->Mean(min,minData)*fExtrapolation->Integral(min,minData);
    integral += fExtrapolation->Integral(min,minData);
    cout << " Low "<< mean << " " << integral << endl;
    
  }
  
  if (max > maxData) {
    // Compute integral exploiting root function to calculate moments, "unnormalizing" them
    mean     += fExtrapolation->Mean(maxData,max)*fExtrapolation->Integral(maxData,max);
    integral += fExtrapolation->Integral(maxData,max);
    cout << " Hi "<< mean << " " << integral << endl;
  } 
  Float_t err2 = 0;
  
  for(Int_t ibin = minDataBin; ibin <= maxDataBin; ibin++){
    if(hData->GetBinCenter(ibin) < min) continue;
    if(hData->GetBinCenter(ibin) > max) continue;
    if(hData->GetBinContent(ibin)) {
      mean     = mean + (hData->GetBinCenter(ibin) *  hData->GetBinWidth(ibin)* hData->GetBinContent(ibin));
      err2     = err2 + TMath::Power(hData->GetBinError(ibin) * hData->GetBinCenter(ibin) *  hData->GetBinWidth(ibin),2);
      integral = integral + hData->GetBinContent(ibin) * hData->GetBinWidth(ibin);
    }
    else {
      Double_t locMin = hData->GetBinLowEdge(ibin);
      Double_t locMax = hData->GetBinLowEdge(ibin) + hData->GetBinWidth(ibin);
      mean     += fExtrapolation->Mean(locMin,locMax)*fExtrapolation->Integral(locMin,locMax);
      Double_t deltaIntegral = fExtrapolation->Integral(locMin,locMax);
      integral += deltaIntegral;
      cout << "WARNING: bin " << ibin << " is empty, patching with function (+" << deltaIntegral<<")" << endl;
    }
  }
  cout << " Data "<< mean << " " << integral << endl;
  
  mean = mean/integral;
  error = TMath::Sqrt(err2)/integral;


}

TH1F * AliPWGHistoTools::CombineHistos(const TH1 * h1, TH1 * h2, const TH1* htemplate, Float_t renorm1){
  // Combine two histos. This assumes the histos have the same binning
  // in the overlapping region. It computes the arithmetic mean in the
  // overlapping region and assigns as an error the relative error
  // h2. TO BE IMPROVED

  TH1F * hcomb = (TH1F*) htemplate->Clone(TString(h1->GetName())+"_"+h2->GetName());

  TH1F * h1local = (TH1F*) h1->Clone();
  h1local->Scale(renorm1);
  
  Int_t nBinComb = hcomb->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nBinComb; ibin++){
    Int_t ibin1 = h1local->FindBin(hcomb->GetBinCenter(ibin));
    Int_t ibin2 = h2->FindBin(hcomb->GetBinCenter(ibin));
    
      if (!h1local->GetBinContent(ibin1) && !h2->GetBinContent(ibin2) ) {
	// None has data: go to next bin
	hcomb->SetBinContent(ibin,0);
	hcomb->SetBinError  (ibin,0);	
      } else if(h1local->GetBinContent(ibin1) && !h2->GetBinContent(ibin2)) {
	// take data from h1local:
	hcomb->SetBinContent(ibin,h1local->GetBinContent(ibin1));
	hcomb->SetBinError  (ibin,h1local->GetBinError(ibin1));
      } else if(!h1local->GetBinContent(ibin1) && h2->GetBinContent(ibin2)) {
	// take data from h2:
	hcomb->SetBinContent(ibin,h2->GetBinContent(ibin2));
	hcomb->SetBinError  (ibin,h2->GetBinError(ibin2));
      }  else {
	hcomb->SetBinContent(ibin,(h1local->GetBinContent(ibin1) +h2->GetBinContent(ibin2))/2);
	//	hcomb->SetBinError  (ibin,h1local->GetBinError(ibin1)/h1local->GetBinContent(ibin1)*hcomb->GetBinContent(ibin));
	hcomb->SetBinError  (ibin,h2->GetBinError(ibin2)/h2->GetBinContent(ibin2)*hcomb->GetBinContent(ibin));
      }


  }
  

  hcomb->SetMarkerStyle(h1local->GetMarkerStyle());
  hcomb->SetMarkerColor(h1local->GetMarkerColor());
  hcomb->SetLineColor  (h1local->GetLineColor());

  hcomb->SetXTitle(h1local->GetXaxis()->GetTitle());
  hcomb->SetYTitle(h1local->GetYaxis()->GetTitle());
  delete h1local;
  return hcomb;  

}


void AliPWGHistoTools::GetFromHistoGraphDifferentX(const TH1F * h, TF1 * f, TGraphErrors ** gBarycentre, TGraphErrors ** gXlw) {

  // Computes the baycentre in each bin and the xlw as defined in NIMA
  // 355 - 541 using f. Returs 2 graphs with the same y content of h
  // but with a different x (barycentre and xlw)

  Int_t nbin = h->GetNbinsX();
  
  (*gBarycentre) = new TGraphErrors();
  (*gXlw)        = new TGraphErrors();

  Int_t ipoint = 0;
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    Float_t min = h->GetBinLowEdge(ibin);
    Float_t max = h->GetBinLowEdge(ibin+1);
    Double_t xerr = 0;
    Double_t xbar = f->Mean(min,max);
    // compute xLW
    Double_t temp = 1./(max-min) * f->Integral(min,max);
    Double_t epsilon   = 0.000000001;
    Double_t increment = 0.0000000001;
    Double_t xLW = min;

    while ((f->Eval(xLW)- temp) > epsilon) {
      xLW += increment;
      if(xLW > max) {
	Printf("Cannot find xLW");
	break;
      }
    }
      
    if (h->GetBinContent(ibin)!=0) {
      (*gBarycentre)->SetPoint     (ipoint,   xbar, h->GetBinContent(ibin));
      (*gBarycentre)->SetPointError(ipoint,   xerr, h->GetBinError(ibin));
      (*gXlw)       ->SetPoint     (ipoint,   xLW,  h->GetBinContent(ibin));
      (*gXlw)       ->SetPointError(ipoint,   xerr, h->GetBinError(ibin));
      ipoint++;
    }
  }
  
  (*gBarycentre)->SetMarkerStyle(h->GetMarkerStyle());
  (*gBarycentre)->SetMarkerColor(h->GetMarkerColor());
  (*gBarycentre)->SetLineColor(h->GetLineColor());

  (*gBarycentre)->SetTitle(h->GetTitle());
  (*gBarycentre)->SetName(TString("g_")+h->GetName());

  (*gXlw)->SetMarkerStyle(h->GetMarkerStyle());
  (*gXlw)->SetMarkerColor(h->GetMarkerColor());
  (*gXlw)->SetLineColor(h->GetLineColor());
  (*gXlw)->SetTitle(h->GetTitle());
  (*gXlw)->SetName(TString("g_")+h->GetName());


}


Float_t AliPWGHistoTools::GetMean(TH1F * h, Float_t min, Float_t max, Float_t * error) {

  // Get the mean of histo in a range; root is not reliable in sub
  // ranges with variable binning.  
  Int_t minBin = h->FindBin(min);
  Int_t maxBin = h->FindBin(max-0.00001);

  Float_t mean = 0 ;
  Float_t integral = 0;
  Float_t err2 = 0;
  for(Int_t ibin = minBin; ibin <= maxBin; ibin++){
    mean     = mean + (h->GetBinCenter(ibin) *  h->GetBinWidth(ibin)* h->GetBinContent(ibin));
    err2     = err2 + TMath::Power(h->GetBinError(ibin) * h->GetBinCenter(ibin) *  h->GetBinWidth(ibin),2);
    integral = integral + h->GetBinContent(ibin) * h->GetBinWidth(ibin);
  }
  
  Float_t value = mean/integral;  
  if (error) (*error) = TMath::Sqrt(err2);
  return value;


}

void AliPWGHistoTools::GetMean(TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) {

  // Get the mean of function in a range; If normPar is >= 0, it means
  // that the function is defined such that par[normPar] is its
  // integral.  In this case the error on meanpt can be calculated
  // correctly. Otherwise, the function is normalized in get moment,
  // but the error is not computed correctly.

  return GetMoment("fMean", TString("x*")+func->GetExpFormula(), func, mean, error, min, max, normPar);

}

void AliPWGHistoTools::GetMeanSquare(TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) {

  // Get the mean2 of function in a range;  If normPar is >= 0, it means
  // that the function is defined such that par[normPar] is its
  // integral.  In this case the error on meanpt can be calculated
  // correctly. Otherwise, the function is normalized in get moment,
  // but the error is not computed correctly.

  return GetMoment("fMean2", TString("x*x*")+func->GetExpFormula(), func, mean, error, min, max, normPar);


}

void AliPWGHistoTools::GetMoment(TString name, TString var, TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) {

  // returns the integral of a function derived from func by prefixing some operation.
  // the derived function MUST have the same parameter in the same order
  // Used as a base method for mean and mean2
  //  If normPar is >= 0, it means that the function is defined such
  // that par[normPar] is its integral.  In this case the error on
  // meanpt can be calculated correctly. Otherwise, the function is
  // normalized using its numerical integral, but the error is not computed
  // correctly. 

  // TODO:
  // - improve to propagate error also in the case you need the
  //   numerical integrals (will be rather slow)
  // - this version assumes that func is defined using a
  //   TFormula. Generalize to the case of a C++ function

  if (normPar<0) Printf("AliPWGHistoTools::GetMoment: Warning: If normalization is required, the error may bot be correct");
  if (!strcmp(func->GetExpFormula(),"")) Printf("AliPWGHistoTools::GetMoment: Warning: Empty formula in the base function");
  Int_t npar = func->GetNpar();

  // Definition changes according to the value of normPar
  TF1 * f = normPar < 0 ? 
    new TF1(name, var) :	              // not normalized
    new TF1(name, var+Form("/[%d]",normPar)); // normalized with par normPar

  // integr is used to normalize if no parameter is provided
  Double_t integr  = normPar < 0 ? func->Integral(min,max) : 1;
  
  // The parameter of the function used to compute the mean should be
  // the same as the parent function: fixed if needed and they should
  // also have the same errors.

  //  cout << "npar :" << npar << endl;
  
  for(Int_t ipar = 0; ipar < npar; ipar++){
    Double_t parmin, parmax;
    Double_t value = func->GetParameter(ipar);
    f->SetParameter(ipar,value);
    func->GetParLimits(ipar, parmin, parmax);
    if ( parmin == parmax )   {
      //      if ( parmin || (parmin == 1 && !value) ) { // not sure we I check parmin == 1 here. 
      if ( parmin || (TMath::Abs(parmin-1)<0.000001 && !value) ) { // not sure we I check parmin == 1 here. Changed like this because of coding conventions. Does it still work? FIXME
	f->FixParameter(ipar,func->GetParameter(ipar));
	//	cout << "Fixing " << ipar << "("<<value<<","<<parmin<<","<<parmax<<")"<<endl;
      }       
      else {
	f->SetParError (ipar,func->GetParError(ipar) );
	//	cout << "Setting Err" << ipar << "("<<func->GetParError(ipar)<<")"<<endl;      
      }
    }
    else {
      f->SetParError (ipar,func->GetParError(ipar) );
      //      cout << "Setting Err" << ipar << "("<<func->GetParError(ipar)<<")"<<endl;      
    }
  }

//   f->Print();
//   cout << "----" << endl;
//   func->Print();

  mean  = normPar < 0 ? f->Integral     (min,max)/integr : f->Integral     (min,max);
  error = normPar < 0 ? f->IntegralError(min,max)/integr : f->IntegralError(min,max);
//   cout << "Mean " << mean <<"+-"<< error<< endl;
//   cout << "Integral Error "  << error << endl;
  
}



Bool_t AliPWGHistoTools::Fit (TH1 * h1, TF1* func, Float_t min, Float_t max) { 
  
  // Fits h1 with func, preforms several checks on the quality of the
  // fit and tries to improve it. If the fit is not good enough, it
  // returs false.

  Double_t amin; Double_t edm; Double_t errdef; Int_t nvpar; Int_t nparx;
  TVirtualFitter *fitter;
  cout << "--- Fitting : " << h1->GetName() << " ["<< h1->GetTitle() <<"] ---"<< endl;

  h1-> Fit(func,"IME0","",min,max);      
  Int_t fitResult = h1-> Fit(func,"IME0","",min,max);      
//   h1-> Fit(func,"0","",min,max);      
//   Int_t fitResult = h1-> Fit(func,"0IE","",min,max);      
  

// From TH1:
// The fitStatus is 0 if the fit is OK (i.e no error occurred).  The
// value of the fit status code is negative in case of an error not
// connected with the minimization procedure, for example when a wrong
// function is used.  Otherwise the return value is the one returned
// from the minimization procedure.  When TMinuit (default case) or
// Minuit2 are used as minimizer the status returned is : fitStatus =
// migradResult + 10*minosResult + 100*hesseResult +
// 1000*improveResult.  TMinuit will return 0 (for migrad, minos,
// hesse or improve) in case of success and 4 in case of error (see
// the documentation of TMinuit::mnexcm). So for example, for an error
// only in Minos but not in Migrad a fitStatus of 40 will be returned.
// Minuit2 will return also 0 in case of success and different values
// in migrad minos or hesse depending on the error. See in this case
// the documentation of Minuit2Minimizer::Minimize for the
// migradResult, Minuit2Minimizer::GetMinosError for the minosResult
// and Minuit2Minimizer::Hesse for the hesseResult.  If other
// minimizers are used see their specific documentation for the status
// code returned.  For example in the case of Fumili, for the status
// returned see TFumili::Minimize.
 

  if( gMinuit->fLimset ) {
    Printf("ERROR: AliPWGHistoTools: Parameters at limits");
    return kFALSE;
  } 


  ///
  fitter = TVirtualFitter::GetFitter();   
  Int_t  fitStat = fitter->GetStats(amin, edm, errdef, nvpar, nparx);  

  if( ( (fitStat < 3  && gMinuit->fCstatu != "UNCHANGED ")|| (edm > 1e6) || (fitResult !=0 && fitResult < 4000) ) && 
      TString(gMinuit->fCstatu) != "SUCCESSFUL"  &&
      TString(gMinuit->fCstatu) != "CONVERGED "  ) {
    if(fitStat < 3 && gMinuit->fCstatu != "UNCHANGED ") {
      Printf("WARNING: AliPWGHistoTools: Cannot properly compute errors");
      if (fitStat == 0) Printf(" not calculated at all");
      if (fitStat == 1) Printf(" approximation only, not accurate");
      if (fitStat == 2) Printf(" full matrix, but forced positive-definite");
    }

    if (edm > 1e6) {

      Printf("WARNING: AliPWGHistoTools: Huge EDM  (%f): Fit probably failed!", edm);
    }
    if (fitResult != 0) {
      Printf("WARNING: AliPWGHistoTools: Fit Result (%d)", fitResult);
    }
      
    Printf ("AliPWGHistoTools: Trying Again with Strategy = 2");

    gMinuit->Command("SET STRATEGY 2"); // more effort
    fitResult = 0;
    fitResult = h1-> Fit(func,"0","",min,max);      
    fitResult = h1-> Fit(func,"IME0","",min,max);      
    fitResult = h1-> Fit(func,"IME0","",min,max);      
      
    fitter = TVirtualFitter::GetFitter();   
  
    fitStat = fitter->GetStats(amin, edm, errdef, nvpar, nparx);  

    if(fitStat < 3 && gMinuit->fCstatu != "UNCHANGED ") {
      Printf("ERROR: AliPWGHistoTools: Cannot properly compute errors");
      if (fitStat == 0) Printf(" not calculated at all");
      if (fitStat == 1) Printf(" approximation only, not accurate");
      if (fitStat == 2) Printf(" full matrix, but forced positive-definite");
      cout << "[" <<gMinuit->fCstatu<<"]" << endl;
      return kFALSE;
    }

    if (edm > 1e6) {
      Printf("ERROR: AliPWGHistoTools: Huge EDM  (%f): Fit probably failed!", edm);

      return kFALSE;
    }
    if (fitResult != 0) {
      Printf("ERROR: AliPWGHistoTools: Fit Result (%d)", fitResult);

      return kFALSE;
    }
      
    gMinuit->Command("SET STRATEGY 1"); // back to normal value

  }

  cout << "---- FIT OK ----" << endl;
  
  return kTRUE;
  
}

Int_t AliPWGHistoTools::GetLowestNotEmptyBin(const TH1*h) {

  // Return the index of the lowest non empty bin in the histo h

  Int_t nbin = h->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    if(h->GetBinContent(ibin)>0) return ibin;
  }
  
  return -1;

}

Int_t AliPWGHistoTools::GetHighestNotEmptyBin(const TH1*h) {

  // Return the index of the highest non empty bin in the histo h

  Int_t nbin = h->GetNbinsX();
  for(Int_t ibin = nbin; ibin > 0; ibin--){
    if(h->GetBinContent(ibin)>0) return ibin;
  }
  
  return -1;

}

void AliPWGHistoTools::GetResiduals(const TGraphErrors * gdata, const TF1 * func, TH1F ** hres, TGraphErrors ** gres) {

  // Returns a graph of residuals vs point and the res/err distribution

  Int_t npoint = gdata->GetN();

  (*gres) =new TGraphErrors();
  (*hres) = new TH1F(TString("hres_")+gdata->GetName()+"-"+func->GetName(),
                  TString("hres_")+gdata->GetName()+"-"+func->GetName(),
                  20,-5,5);


  for(Int_t ipoint = 0; ipoint < npoint; ipoint++){
    Float_t x   = gdata->GetX()[ipoint];
    Float_t res = (gdata->GetY()[ipoint] - func->Eval(x))/func->Eval(x);
    Float_t err = gdata->GetEY()[ipoint]/func->Eval(x);
    (*hres)->Fill(res/err);
    (*gres)->SetPoint(ipoint, x, res/err);
    //    (*gres)->SetPointError(ipoint, 0, err);
    
  }
  
  (*gres)->SetMarkerStyle(gdata->GetMarkerStyle());
  (*gres)->SetMarkerColor(gdata->GetMarkerColor());
  (*gres)->SetLineColor  (gdata->GetLineColor());
  (*gres)->GetHistogram()->GetYaxis()->SetTitle("(data-function)/function");
  (*hres)->SetMarkerStyle(gdata->GetMarkerStyle());
  (*hres)->SetMarkerColor(gdata->GetMarkerColor());
  (*hres)->SetLineColor  (gdata->GetLineColor());



}

void AliPWGHistoTools::GetResiduals(const TH1F* hdata, const TF1 * func, TH1F ** hres, TH1F ** hresVsBin) {

  // Returns an histo of residuals bin by bin and the res/err distribution

  if (!func) {
    Printf("AliPWGHistoTools::GetResiduals: No function provided");
    return;
  }
  if (!hdata) {
    Printf("AliPWGHistoTools::GetResiduals: No data provided");
    return;
  }

  (*hresVsBin) = (TH1F*) hdata->Clone(TString("hres_")+hdata->GetName());
  (*hresVsBin)->Reset();
  (*hres) = new TH1F(TString("hres_")+hdata->GetName()+"-"+func->GetName(),
		     TString("hres_")+hdata->GetName()+"-"+func->GetName(),
		     20,-5,5);

  Int_t nbin = hdata->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    if(!hdata->GetBinContent(ibin)) continue;
    Float_t res = (hdata->GetBinContent(ibin) - func->Eval(hdata->GetBinCenter(ibin)) ) / 
      func->Eval(hdata->GetBinCenter(ibin));
    Float_t err = hdata->GetBinError  (ibin) /  func->Eval(hdata->GetBinCenter(ibin));
    (*hresVsBin)->SetBinContent(ibin,res);
    (*hresVsBin)->SetBinError  (ibin,err);
    (*hres)->Fill(res/err);
    
  }
  
  (*hresVsBin)->SetMarkerStyle(hdata->GetMarkerStyle());
  (*hresVsBin)->SetMarkerColor(hdata->GetMarkerColor());
  (*hresVsBin)->SetLineColor  (hdata->GetLineColor()  );
  (*hresVsBin)->GetYaxis()->SetTitle("(data-function)/function");
  (*hres)->SetMarkerStyle(hdata->GetMarkerStyle());
  (*hres)->SetMarkerColor(hdata->GetMarkerColor());
  (*hres)->SetLineColor  (hdata->GetLineColor()  );

}

void AliPWGHistoTools::GetYield(TH1* h,  TF1 * f, Double_t &yield, Double_t &yieldError, Float_t min, Float_t max,
			  Double_t *partialYields, Double_t *partialYieldsErrors){

  // Returns the yield extracted from the data in the histo where
  // there are points and from the fit to extrapolate, in the given
  // range.

  // Partial yields are also returned if the corresponding pointers are non null

  Int_t bin1 = h->FindBin(min);
  Int_t bin2 = h->FindBin(max);
  Float_t bin1Edge = GetLowestNotEmptyBinEdge (h);
  Float_t bin2Edge = GetHighestNotEmptyBinEdge(h);
  Int_t lowestNE  = GetLowestNotEmptyBin (h);
  Int_t highestNE = GetHighestNotEmptyBin(h);

  Double_t integralFromHistoError = 0;
  Double_t integralFromHisto = 0;
  //  Double_t integralFromHisto = DoIntegral(h,bin1,bin2,-1,-1,-1,-1,integralFromHistoError,"width",1);

  for(Int_t ibin = TMath::Max(bin1,lowestNE); ibin <= TMath::Min(bin2,highestNE); ibin++){
    if(h->GetBinContent(ibin)) {
      integralFromHisto = integralFromHisto + h->GetBinContent(ibin) * h->GetBinWidth(ibin);
      integralFromHistoError    = integralFromHistoError    + h->GetBinError(ibin)*h->GetBinError(ibin) * h->GetBinWidth(ibin);
    }
    else {
      Double_t locMin = h->GetBinLowEdge(ibin);
      Double_t locMax = h->GetBinLowEdge(ibin) + h->GetBinWidth(ibin);
      Double_t deltaIntegral = f->Integral(locMin,locMax);
      integralFromHisto += deltaIntegral;
      integralFromHistoError = integralFromHistoError + f->IntegralError(locMin,locMax)*f->IntegralError(locMin,locMax);
      cout << "WARNING: bin " << ibin << " is empty, patching with function (+" << deltaIntegral<<")" << endl;
    }
  }
  
  integralFromHistoError = TMath::Sqrt(integralFromHistoError);// ok, this is a bit dumb. Just in case the code above is changed again.
  
  
  Double_t integralBelow      = min < bin1Edge ? f->Integral(min,bin1Edge) : 0;
  Double_t integralBelowError = min < bin1Edge ? f->IntegralError(min,bin1Edge) : 0;
  Double_t integralAbove      = max > bin2Edge ? f->Integral(bin2Edge,max) : 0;
  Double_t integralAboveError = max > bin2Edge ? f->IntegralError(bin2Edge,max) : 0;

//   cout << "GetYield INFO" << endl;
//   cout << " " << bin1Edge << " " << bin2Edge << endl;  
//   cout << " " << integralFromHisto      << " " << integralBelow      << " " << integralAbove      << endl;
//   cout << " " << integralFromHistoError << " " << integralBelowError << " " << integralAboveError << endl;
  
  if(partialYields) {
    partialYields[0] = integralFromHisto;
    partialYields[1] = integralBelow;
    partialYields[2] = integralAbove;
  }
  if(partialYieldsErrors) {
    partialYieldsErrors[0] = integralFromHistoError;
    partialYieldsErrors[1] = integralBelowError;
    partialYieldsErrors[2] = integralAboveError;
  }
  yield      = integralFromHisto+integralBelow+integralAbove;
  yieldError = TMath::Sqrt(integralFromHistoError*integralFromHistoError+
			   integralBelowError*integralBelowError+
			   integralAboveError*integralAboveError);

}

TGraphErrors * AliPWGHistoTools::DivideGraphByFunc(const TGraphErrors * g, const TF1 * f, Bool_t invert){ 

  // Divides g/f. If invert == true => f/g

  TGraphErrors * gRatio = new TGraphErrors();
  Int_t npoint = g->GetN();
  for(Int_t ipoint = 0; ipoint < npoint; ipoint++){
    Double_t x = g->GetX()[ipoint];
    Double_t ratio  = invert ? f->Eval(x)/g->GetY()[ipoint] :g->GetY()[ipoint]/f->Eval(x);
    gRatio->SetPoint     (ipoint, x, ratio);
    if(f->Eval(x) && strcmp(g->ClassName(),"TGraphAsymmErrors")) gRatio->SetPointError(ipoint, 0, g->GetEY()[ipoint]/f->Eval(x));
    //    cout << x << " " << g->GetY()[ipoint] << " " << f->Eval(x) << endl;
    
  }
  gRatio->SetMarkerStyle(20);
  //gRatio->Print();
  return gRatio;

}

TGraphErrors * AliPWGHistoTools::Divide2Graphs(const TGraph * g1, const TGraph * g2, Int_t strategy){ 

  // Divides g1/g2, 
  // 2 strategies are possible: 
  // - strategy = 0 looks for point with very close centers (also propagates error in this case)
  // - strategy = 1 Loop over all points from g1 (xi,yi), interpolates the TGraphs with a spline and computes the ratio as yi / SplineG2(xi)

  TGraphErrors * gRatio = new TGraphErrors();

  if(strategy == 0) {
    Int_t ipoint=0;
    Int_t npoint1 = g1->GetN();
    Int_t npoint2 = g2->GetN();
    for(Int_t ipoint1 = 0; ipoint1 < npoint1; ipoint1++){
      Double_t x1 = g1->GetX()[ipoint1];
      for(Int_t ipoint2 = 0; ipoint2 < npoint2; ipoint2++){
	Double_t x2 = g2->GetX()[ipoint2];
	if((TMath::Abs(x1-x2)/(x1+x2)*2)<0.01) {
	  Double_t ratio   = g2->GetY()[ipoint2]  ? g1->GetY()[ipoint1]/g2->GetY()[ipoint2] : 0;
	  Double_t eratio  = 0;
	  if(g1->InheritsFrom("TGraphErrors") && g2->InheritsFrom("TGraphErrors")) {
	    eratio = g2->GetY()[ipoint2]  ? 
	      TMath::Sqrt(g1->GetEY()[ipoint1]*g1->GetEY()[ipoint1]/g1->GetY()[ipoint1]/g1->GetY()[ipoint1] + 
			  g2->GetEY()[ipoint2]/g2->GetY()[ipoint2]/g2->GetY()[ipoint2] ) * ratio
	      : 0;
	  }
	  gRatio->SetPoint     (ipoint, x1, ratio);
	  gRatio->SetPointError(ipoint, 0, eratio);
	  ipoint++;
	  cout << ipoint << " [" << x1 << "] " <<  g1->GetY()[ipoint1] << "/" << g2->GetY()[ipoint2] << " = " << ratio <<"+-"<<eratio<< endl;
	  
	  //    cout << x << " " << g->GetY()[ipoint] << " " << f->Eval(x) << endl;
	}
	
      }
    }
  }
  else if (strategy == 1) {
    TSpline3 * splG2 = new TSpline3("fGraphSplG2",g2);
    Int_t npoint1 = g1->GetN();
    for(Int_t ipoint1 = 0; ipoint1 < npoint1; ipoint1++){
      Double_t x1    = g1->GetX()[ipoint1];
      Double_t y1    = g1->GetY()[ipoint1];
      Double_t ratio = y1/splG2->Eval(x1);
      gRatio->SetPoint(ipoint1, x1, ratio);

    }
    delete splG2;
  }
  else Printf("AliPWGHistoTools::Divide2Graphs(): Invalid strategy [%d]", strategy);

    gRatio->SetMarkerStyle(20);
  //gRatio->Print();
  return gRatio;

}
TGraphErrors * AliPWGHistoTools::Add2Graphs(const TGraphErrors * g1, const TGraphErrors * g2){ 

  // Sums g1/g2, looks for point with very close centers
  Int_t ipoint=0;
  TGraphErrors * gSum = new TGraphErrors();
  Int_t npoint1 = g1->GetN();
  Int_t npoint2 = g2->GetN();
  for(Int_t ipoint1 = 0; ipoint1 < npoint1; ipoint1++){
    Double_t x1 = g1->GetX()[ipoint1];
    for(Int_t ipoint2 = 0; ipoint2 < npoint2; ipoint2++){
      Double_t x2 = g2->GetX()[ipoint2];
      if((TMath::Abs(x1-x2)/(x1+x2)*2)<0.01) {
	Double_t sum   = g1->GetY()[ipoint1]+g2->GetY()[ipoint2];
	Double_t esum  = 0;
	if(g1->InheritsFrom("TGraphErrors") && g2->InheritsFrom("TGraphErrors")) {
	  esum = 
	    TMath::Sqrt(g1->GetEY()[ipoint1]*g1->GetEY()[ipoint1] + 
			g2->GetEY()[ipoint2]*g2->GetEY()[ipoint2] ) ;
	    
	}
	gSum->SetPoint     (ipoint, x1, sum);
	gSum->SetPointError(ipoint, 0, esum);
	ipoint++;
	cout << ipoint << " [" << x1 << "] " <<  g1->GetY()[ipoint1] << "+" << g2->GetY()[ipoint2] << " = " << sum <<"+-"<<esum<< endl;
	
    //    cout << x << " " << g->GetY()[ipoint] << " " << f->Eval(x) << endl;
      }
    
    }
  }
  gSum->SetMarkerStyle(20);
  //gSum->Print();
  return gSum;

}

TGraphErrors * AliPWGHistoTools::DivideGraphByHisto(const TGraphErrors * g, TH1 * h, Bool_t invert){ 

  // Divides g/h. If invert == true => h/g

  Bool_t skipError = kFALSE;
  if(!strcmp(g->ClassName(),"TGraph")) skipError = kTRUE;
  if(!strcmp(g->ClassName(),"TGraphAsymmErrors")) skipError = kTRUE;
  if(skipError) {
    Printf("WARNING: Skipping graph errors");
  }
  TGraphErrors * gRatio = new TGraphErrors();
  Int_t npoint = g->GetN();
  for(Int_t ipoint = 0; ipoint < npoint; ipoint++){
    Double_t xj  = g->GetX()[ipoint];
    Double_t yj  = g->GetY()[ipoint];
    Double_t yje = skipError ? 0 : g->GetEY()[ipoint];

    Int_t binData = h->FindBin(xj);
    Double_t yd   = h->GetBinContent(binData);
    Double_t yde  = h->GetBinError(binData);
    Double_t xd   = h->GetBinCenter(binData);
    

     
    if (!yd) continue;
    
    if (TMath::Abs(xd-xj)/TMath::Abs(xd) > 0.01){
      Printf( "WARNING: bin center (%f)  and x graph (%f) are more than 1 %% away, skipping",xd,xj );
      continue;
      
    }

    Double_t ratio = invert ? yd/yj : yj/yd;

    gRatio->SetPoint(ipoint, xj, ratio);
    gRatio->SetPointError(ipoint, 0, TMath::Sqrt(yde*yde/yd/yd + yje*yje/yj/yj)*ratio);
    //    gRatio->SetPointError(ipoint, 0, yje/yj * ratio);
  }

  return gRatio;


}

TH1F * AliPWGHistoTools::DivideHistoByFunc(TH1F * h, TF1 * f, Bool_t invert){ 

  // Divides h/f. If invert == true => f/g
  // Performs the integral of f on the bin range to perform the ratio
  // Returns a histo with the same binnig as h

  // Prepare histo for ratio
  TH1F * hRatio = (TH1F*) h->Clone(TString("hRatio_")+h->GetName()+"_"+f->GetName());
  hRatio->Reset();
  // Set y title
  if(!invert) hRatio->SetYTitle(TString(h->GetName())+"/"+f->GetName());
  else        hRatio->SetYTitle(TString(f->GetName())+"/"+h->GetName());

  // Loop over all bins
  Int_t nbin = hRatio->GetNbinsX();

  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    Double_t yhisto = h->GetBinContent(ibin);
    Double_t yerror = h->GetBinError(ibin);
    Double_t xmin   = h->GetBinLowEdge(ibin);
    Double_t xmax   = h->GetBinLowEdge(ibin+1);
    Double_t yfunc  = f->Integral(xmin,xmax)/(xmax-xmin);
    Double_t ratio = yfunc&&yhisto ? (invert ? yfunc/yhisto : yhisto/yfunc) : 0 ;
    Double_t error = yfunc ? yerror/yfunc : 0  ;
    hRatio->SetBinContent(ibin,ratio);
    hRatio->SetBinError  (ibin,error);
  }

  return hRatio;

}

void AliPWGHistoTools::WeightedMean(Int_t npoints, const Double_t *x, const Double_t *xerr, Double_t &mean, Double_t &meanerr){

  // Performs the weighted mean of npoints numbers in x with errors in xerr

  mean = 0;
  meanerr = 0;

  Double_t sumweight = 0;

  for (Int_t ipoint = 0; ipoint < npoints; ipoint++){
    
    Double_t xerr2 = xerr[ipoint]*xerr[ipoint];
    if(xerr2>0){
      //      cout << "xe2 " << xerr2 << endl;
      Double_t weight = 1. / xerr2;
      sumweight += weight;
      mean += weight * x[ipoint];
    }// else cout << " Skipping " << ipoint << endl;
    
  }


  if(sumweight){
    mean /= sumweight;
    meanerr = TMath::Sqrt(1./ sumweight);
  }
  else {
    //    cout << " No sumweight" << endl;
    mean = 0;
    meanerr = 0;
  }

  
}

TH1 * AliPWGHistoTools::GetRelativeError(TH1 * h){
  // Returns an histogram with the same binning as h, filled with the relative error bin by bin
  TH1 * hrel = (TH1*) h->Clone(TString(h->GetName())+"_rel");
  hrel->Reset();
  Int_t nbin = hrel->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    if(h->GetBinError(ibin)){
      hrel->SetBinContent(ibin,h->GetBinError(ibin)/h->GetBinContent(ibin));
    } else {
      hrel->SetBinContent(ibin,0);
    }
    hrel->SetBinError(ibin,0);
  }
  
  return hrel;
}


void AliPWGHistoTools::GetValueAndError(TH1 * hdest, const TH1 * hvalue, const TH1 * herror, Bool_t isPercentError) {
  
  // Put into source, bin-by-bin, the values from hvalue and the
  // errors from content from herror. 
  // Used mainly to combine histos of systemati errors with their spectra
  // Set isPercentError to kTRUE if the error is given in % 

  if(hdest == NULL){ 
    Printf("AliPWGHistoTools::GetValueAndError Errror: hdest is null");
    return;
  }


  Int_t nbin  = hdest->GetNbinsX();
  Int_t nBinSourceVal = hvalue->GetNbinsX();
  Int_t nBinSourceErr = herror->GetNbinsX();
  
  for(Int_t iBinDest = 1; iBinDest <= nbin; iBinDest++){
    Float_t lowPtDest=hdest->GetBinLowEdge(iBinDest);
    Float_t binWidDest=hdest->GetBinWidth(iBinDest);
    // Loop over Source bins and find overlapping bins to Dest
    // First value then error
    // Value
    Bool_t foundValue = kFALSE;
    for(Int_t iBinSourceVal=1; iBinSourceVal<=nBinSourceVal; iBinSourceVal++){
      Float_t lowPtSource=  hvalue->GetBinLowEdge(iBinSourceVal) ;
      Float_t binWidSource= hvalue->GetBinWidth(iBinSourceVal);
      if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){
	Double_t content = hvalue->GetBinContent(iBinSourceVal);
	hdest->SetBinContent(iBinDest, content);
	foundValue = kTRUE;
	break;
      }
    }
    // if (!foundValue){
    //   Printf("AliPWGHistoTools::GetValueAndError: Error: cannot find matching value source bin for destination %d",iBinDest);
    // }

    // Error
    Bool_t foundError = kFALSE;
    for(Int_t iBinSourceErr=1; iBinSourceErr<=nBinSourceErr; iBinSourceErr++){
      Float_t lowPtSource=  herror->GetBinLowEdge(iBinSourceErr) ;
      Float_t binWidSource= herror->GetBinWidth(iBinSourceErr);
      if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){
	Double_t error = herror->GetBinContent(iBinSourceErr);
	//	cout << "-> " << iBinDest << " " << error << " " << hdest->GetBinContent(iBinDest) << endl;
	
	hdest->SetBinError(iBinDest, isPercentError ? error * hdest->GetBinContent(iBinDest) : error);
	foundError=kTRUE;
	break;
      }      
    }
    // if (!foundError ){
    //   Printf("AliPWGHistoTools::GetValueAndError: Error: cannot find matching error source bin for destination %d",iBinDest);
    // }
  }
  

}

void AliPWGHistoTools::AddHisto(TH1 * hdest, const TH1* hsource, Bool_t getMirrorBins) {

  // Adds hsource to hdest bin by bin, even if they have a different
  // binning If getMirrorBins is true, it takes the negative bins
  // (Needed because sometimes the TPC uses the positive axis for
  // negative particles and the possitive axis for positive
  // particles).


  if (hdest == NULL) {
    Printf("Error: hdest is NULL\n");
    return;
  } 
  if (hsource == NULL) {
    Printf("Error: hsource is NULL\n");
    return;
  } 

  Int_t nBinSource = hsource->GetNbinsX();
  Int_t nBinDest = hdest->GetNbinsX();

  // Loop over destination bins, 
  for(Int_t iBinDest=1; iBinDest<=nBinDest; iBinDest++){
    Float_t lowPtDest=hdest->GetBinLowEdge(iBinDest);
    Float_t binWidDest=hdest->GetBinWidth(iBinDest);
    // Loop over Source bins and find overlapping bins to Dest
    Bool_t found = kFALSE;
    for(Int_t iBinSource=1; iBinSource<=nBinSource; iBinSource++){      
      Float_t lowPtSource= getMirrorBins ? -hsource->GetBinLowEdge(iBinSource)+hsource->GetBinWidth(iBinSource) : hsource->GetBinLowEdge(iBinSource) ;
      Float_t binWidSource= hsource->GetBinWidth(iBinSource)  ;
      if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){
	Float_t dest=hdest->GetBinContent(iBinDest);
	Float_t source=hsource->GetBinContent(iBinSource);
	Float_t edest=hdest->GetBinError(iBinDest);
	Float_t esource=hsource->GetBinError(iBinSource);
	Double_t cont=dest+source;
	Double_t econt=TMath::Sqrt(edest*edest+esource*esource);
	hdest->SetBinContent(iBinDest,cont);
	hdest->SetBinError  (iBinDest,econt);
	found = kTRUE;
	
	break;
      }
    }
    // if (!found){
    //   Printf("Error: cannot find matching source bin for destination %d",iBinDest);
    // }
  }


}

void AliPWGHistoTools::GetHistoCombinedErrors(TH1 * hdest, const TH1 * h1) {

  // Combine the errors of hdest with the errors of h1, summing in
  // quadrature. Results are put in hdest. Histograms are assumed to
  // have the same binning

  Int_t nbin = hdest->GetNbinsX();
  for(Int_t ibin = 1; ibin <= nbin; ibin++){
    Double_t e1 = hdest->GetBinError(ibin);
    Double_t e2 = h1->GetBinError(ibin);
    hdest->SetBinError(ibin, TMath::Sqrt(e1*e1+e2*e2));
  }
  
  
}

TH1F * AliPWGHistoTools::DivideHistosDifferentBins(const TH1F* h1, const TH1F* h2) {
  // Divides 2 histos even if they have a different binning. Finds
  // overlapping bins and divides them
  
  // 1. clone histo
  TH1F * hRatio = new TH1F(*h1);
  Int_t nBinsH1=h1->GetNbinsX();
  Int_t nBinsH2=h2->GetNbinsX();
  // Loop over H1 bins, 
  for(Int_t iBin=1; iBin<=nBinsH1; iBin++){
    hRatio->SetBinContent(iBin,0.);
    hRatio->SetBinContent(iBin,0.);
    Float_t lowPtH1=h1->GetBinLowEdge(iBin);
    Float_t binWidH1=h1->GetBinWidth(iBin);
    // Loop over H2 bins and find overlapping bins to H1
    for(Int_t jBin=1; jBin<=nBinsH2; jBin++){
      Float_t lowPtH2=h2->GetBinLowEdge(jBin);
      Float_t binWidH2=h2->GetBinWidth(jBin);
      if(TMath::Abs(lowPtH1-lowPtH2)<0.001 && TMath::Abs(binWidH2-binWidH1)<0.001){
	Float_t numer=h1->GetBinContent(iBin);
	Float_t denom=h2->GetBinContent(jBin);
	Float_t enumer=h1->GetBinError(iBin);
	Float_t edenom=h2->GetBinError(jBin);
	Double_t ratio=0.;
	Double_t eratio=0.;
	if(numer>0. && denom>0.){
	  ratio=numer/denom;
	  eratio=ratio*TMath::Sqrt((enumer/numer)*(enumer/numer)+(edenom/denom)*(edenom/denom));
	}
	hRatio->SetBinContent(iBin,ratio);
	hRatio->SetBinError(iBin,eratio);
	break;
      }
    }
  }
  return hRatio;
}

Double_t AliPWGHistoTools::DoIntegral(TH1* h, Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2, Double_t & error ,
				Option_t *option, Bool_t doError) 
{
   // function to compute integral and optionally the error  between the limits
   // specified by the bin number values working for all histograms (1D, 2D and 3D)
  // copied from TH! to fix a bug still present in 5-27-06b
   Int_t nbinsx = h->GetNbinsX();
   if (binx1 < 0) binx1 = 0;
   if (binx2 > nbinsx+1 || binx2 < binx1) binx2 = nbinsx+1;
   if (h->GetDimension() > 1) {
      Int_t nbinsy = h->GetNbinsY();
      if (biny1 < 0) biny1 = 0;
      if (biny2 > nbinsy+1 || biny2 < biny1) biny2 = nbinsy+1;
   } else {
      biny1 = 0; biny2 = 0;
   }
   if (h->GetDimension() > 2) {
      Int_t nbinsz = h->GetNbinsZ();
      if (binz1 < 0) binz1 = 0;
      if (binz2 > nbinsz+1 || binz2 < binz1) binz2 = nbinsz+1;
   } else {
      binz1 = 0; binz2 = 0;
   }

   //   - Loop on bins in specified range
   TString opt = option;
   opt.ToLower();
   Bool_t width   = kFALSE;
   if (opt.Contains("width")) width = kTRUE;


   Double_t dx = 1.;
   Double_t dy = 1.;
   Double_t dz = 1.;
   Double_t integral = 0;
   Double_t igerr2 = 0;
   for (Int_t binx = binx1; binx <= binx2; ++binx) {
     if (width) dx = h->GetXaxis()->GetBinWidth(binx);
     for (Int_t biny = biny1; biny <= biny2; ++biny) {
       if (width) dy = h->GetYaxis()->GetBinWidth(biny);
       for (Int_t binz = binz1; binz <= binz2; ++binz) {
	 if (width) dz = h->GetZaxis()->GetBinWidth(binz);
	 Int_t bin  = h->GetBin(binx, biny, binz);
	 if (width) integral += h->GetBinContent(bin)*dx*dy*dz;
	 else       integral += h->GetBinContent(bin);
	 if (doError) {
	   if (width)  igerr2 += h->GetBinError(bin)*h->GetBinError(bin)*dx*dy*dz*dx*dy*dz;
	   else        igerr2 += h->GetBinError(bin)*h->GetBinError(bin);
	 }
	 //	 cout << h->GetBinContent(bin) << " " <<  h->GetBinError(bin) << " " << dx*dy*dz << " "  << integral << " +- " << igerr2 << endl;
	 
       }
     }
   }
   
   if (doError) error = TMath::Sqrt(igerr2);
   return integral;
}

Double_t AliPWGHistoTools::dMtdptFunction(Double_t *x, Double_t *p) {

  // Computes the dmt/dptdeta function using the dN/dpt function
  // This is a protected function used internally by GetdMtdy to integrate dN/dpt function using mt as a weight
  // The mass of the particle is given as p[0]
  Double_t pt   = x[0];
  Double_t mass = p[0]; 
  Double_t mt = TMath::Sqrt(pt*pt + mass*mass);
  Double_t jacobian = pt/mt;
  if(!fdNdptForETCalc){
    Printf("AliPWGHistoTools::dMtdptFunction: ERROR: fdNdptForETCalc not defined");
    return 0;
  }
  Double_t dNdpt = fdNdptForETCalc->Eval(pt);
  return dNdpt*mt*jacobian; // FIXME: do I have to normalize somehow?
  
}

Double_t AliPWGHistoTools::GetdMtdEta(TH1 *hData, TF1 * fExtrapolation, Double_t mass) {
  // Computes dMtdEta integrating dN/dptdy with the proper weights and jacobian.
  Printf("WARNING ALIBWTOOLS::GetdMtdEta: ONLY USING FUNCTION FOR THE TIME BEING, hData");
  if(!hData) {
    Printf("hData not set");
  }

  // Assign the fiunction used internally by dMtdptFunction
  fdNdptForETCalc = fExtrapolation;
  // Create the function to be integrated
  TF1 * funcdMtdPt = new TF1 ("funcdMtdPt", dMtdptFunction, 0.0, 20, 1);
  funcdMtdPt->SetParameter(0,mass);
  // Integrate it
  Double_t dMtdEta = funcdMtdPt->Integral(0,100);
  // Clean up
  fdNdptForETCalc=0;
  delete funcdMtdPt;
  //return 
  return dMtdEta;

}




void AliPWGHistoTools::ScaleGraph(TGraph * g1, Double_t scale) {

  // Scale a grh by the factor scale. If g1 is a TGraph(Asym)Error, errors are also scaled

  if(!g1) return;
  Bool_t isErrors  = TString(g1->ClassName()) == "TGraphErrors";
  Bool_t isAErrors = TString(g1->ClassName()) == "TGraphAsymmErrors";

  Int_t npoint = g1->GetN();
  for (Int_t ipoint = 0; ipoint<npoint; ipoint++) {
    Double_t x = g1->GetX()[ipoint];
    Double_t y = g1->GetY()[ipoint];
    Double_t val = y*scale;
    g1->SetPoint(ipoint, x, val);
    //    if it is  aclass with errors, also scale those
    if (isAErrors) {
      Double_t errxlow  = g1->GetEXlow() [ipoint];      
      Double_t errxhigh = g1->GetEXhigh()[ipoint];      
      Double_t errylow  = g1->GetEYlow() [ipoint]*scale;      
      Double_t erryhigh = g1->GetEYhigh()[ipoint]*scale;      
      ((TGraphAsymmErrors*)g1)->SetPointError(ipoint, errxlow, errxhigh, errylow, erryhigh);     
    }
    else if (isErrors) {
      Double_t erry = g1->GetEY()[ipoint]*scale;
      ((TGraphErrors*)g1)->SetPointError(ipoint, g1->GetEX()[ipoint], erry);
    }
  }


}
 AliPWGHistoTools.cxx:1
 AliPWGHistoTools.cxx:2
 AliPWGHistoTools.cxx:3
 AliPWGHistoTools.cxx:4
 AliPWGHistoTools.cxx:5
 AliPWGHistoTools.cxx:6
 AliPWGHistoTools.cxx:7
 AliPWGHistoTools.cxx:8
 AliPWGHistoTools.cxx:9
 AliPWGHistoTools.cxx:10
 AliPWGHistoTools.cxx:11
 AliPWGHistoTools.cxx:12
 AliPWGHistoTools.cxx:13
 AliPWGHistoTools.cxx:14
 AliPWGHistoTools.cxx:15
 AliPWGHistoTools.cxx:16
 AliPWGHistoTools.cxx:17
 AliPWGHistoTools.cxx:18
 AliPWGHistoTools.cxx:19
 AliPWGHistoTools.cxx:20
 AliPWGHistoTools.cxx:21
 AliPWGHistoTools.cxx:22
 AliPWGHistoTools.cxx:23
 AliPWGHistoTools.cxx:24
 AliPWGHistoTools.cxx:25
 AliPWGHistoTools.cxx:26
 AliPWGHistoTools.cxx:27
 AliPWGHistoTools.cxx:28
 AliPWGHistoTools.cxx:29
 AliPWGHistoTools.cxx:30
 AliPWGHistoTools.cxx:31
 AliPWGHistoTools.cxx:32
 AliPWGHistoTools.cxx:33
 AliPWGHistoTools.cxx:34
 AliPWGHistoTools.cxx:35
 AliPWGHistoTools.cxx:36
 AliPWGHistoTools.cxx:37
 AliPWGHistoTools.cxx:38
 AliPWGHistoTools.cxx:39
 AliPWGHistoTools.cxx:40
 AliPWGHistoTools.cxx:41
 AliPWGHistoTools.cxx:42
 AliPWGHistoTools.cxx:43
 AliPWGHistoTools.cxx:44
 AliPWGHistoTools.cxx:45
 AliPWGHistoTools.cxx:46
 AliPWGHistoTools.cxx:47
 AliPWGHistoTools.cxx:48
 AliPWGHistoTools.cxx:49
 AliPWGHistoTools.cxx:50
 AliPWGHistoTools.cxx:51
 AliPWGHistoTools.cxx:52
 AliPWGHistoTools.cxx:53
 AliPWGHistoTools.cxx:54
 AliPWGHistoTools.cxx:55
 AliPWGHistoTools.cxx:56
 AliPWGHistoTools.cxx:57
 AliPWGHistoTools.cxx:58
 AliPWGHistoTools.cxx:59
 AliPWGHistoTools.cxx:60
 AliPWGHistoTools.cxx:61
 AliPWGHistoTools.cxx:62
 AliPWGHistoTools.cxx:63
 AliPWGHistoTools.cxx:64
 AliPWGHistoTools.cxx:65
 AliPWGHistoTools.cxx:66
 AliPWGHistoTools.cxx:67
 AliPWGHistoTools.cxx:68
 AliPWGHistoTools.cxx:69
 AliPWGHistoTools.cxx:70
 AliPWGHistoTools.cxx:71
 AliPWGHistoTools.cxx:72
 AliPWGHistoTools.cxx:73
 AliPWGHistoTools.cxx:74
 AliPWGHistoTools.cxx:75
 AliPWGHistoTools.cxx:76
 AliPWGHistoTools.cxx:77
 AliPWGHistoTools.cxx:78
 AliPWGHistoTools.cxx:79
 AliPWGHistoTools.cxx:80
 AliPWGHistoTools.cxx:81
 AliPWGHistoTools.cxx:82
 AliPWGHistoTools.cxx:83
 AliPWGHistoTools.cxx:84
 AliPWGHistoTools.cxx:85
 AliPWGHistoTools.cxx:86
 AliPWGHistoTools.cxx:87
 AliPWGHistoTools.cxx:88
 AliPWGHistoTools.cxx:89
 AliPWGHistoTools.cxx:90
 AliPWGHistoTools.cxx:91
 AliPWGHistoTools.cxx:92
 AliPWGHistoTools.cxx:93
 AliPWGHistoTools.cxx:94
 AliPWGHistoTools.cxx:95
 AliPWGHistoTools.cxx:96
 AliPWGHistoTools.cxx:97
 AliPWGHistoTools.cxx:98
 AliPWGHistoTools.cxx:99
 AliPWGHistoTools.cxx:100
 AliPWGHistoTools.cxx:101
 AliPWGHistoTools.cxx:102
 AliPWGHistoTools.cxx:103
 AliPWGHistoTools.cxx:104
 AliPWGHistoTools.cxx:105
 AliPWGHistoTools.cxx:106
 AliPWGHistoTools.cxx:107
 AliPWGHistoTools.cxx:108
 AliPWGHistoTools.cxx:109
 AliPWGHistoTools.cxx:110
 AliPWGHistoTools.cxx:111
 AliPWGHistoTools.cxx:112
 AliPWGHistoTools.cxx:113
 AliPWGHistoTools.cxx:114
 AliPWGHistoTools.cxx:115
 AliPWGHistoTools.cxx:116
 AliPWGHistoTools.cxx:117
 AliPWGHistoTools.cxx:118
 AliPWGHistoTools.cxx:119
 AliPWGHistoTools.cxx:120
 AliPWGHistoTools.cxx:121
 AliPWGHistoTools.cxx:122
 AliPWGHistoTools.cxx:123
 AliPWGHistoTools.cxx:124
 AliPWGHistoTools.cxx:125
 AliPWGHistoTools.cxx:126
 AliPWGHistoTools.cxx:127
 AliPWGHistoTools.cxx:128
 AliPWGHistoTools.cxx:129
 AliPWGHistoTools.cxx:130
 AliPWGHistoTools.cxx:131
 AliPWGHistoTools.cxx:132
 AliPWGHistoTools.cxx:133
 AliPWGHistoTools.cxx:134
 AliPWGHistoTools.cxx:135
 AliPWGHistoTools.cxx:136
 AliPWGHistoTools.cxx:137
 AliPWGHistoTools.cxx:138
 AliPWGHistoTools.cxx:139
 AliPWGHistoTools.cxx:140
 AliPWGHistoTools.cxx:141
 AliPWGHistoTools.cxx:142
 AliPWGHistoTools.cxx:143
 AliPWGHistoTools.cxx:144
 AliPWGHistoTools.cxx:145
 AliPWGHistoTools.cxx:146
 AliPWGHistoTools.cxx:147
 AliPWGHistoTools.cxx:148
 AliPWGHistoTools.cxx:149
 AliPWGHistoTools.cxx:150
 AliPWGHistoTools.cxx:151
 AliPWGHistoTools.cxx:152
 AliPWGHistoTools.cxx:153
 AliPWGHistoTools.cxx:154
 AliPWGHistoTools.cxx:155
 AliPWGHistoTools.cxx:156
 AliPWGHistoTools.cxx:157
 AliPWGHistoTools.cxx:158
 AliPWGHistoTools.cxx:159
 AliPWGHistoTools.cxx:160
 AliPWGHistoTools.cxx:161
 AliPWGHistoTools.cxx:162
 AliPWGHistoTools.cxx:163
 AliPWGHistoTools.cxx:164
 AliPWGHistoTools.cxx:165
 AliPWGHistoTools.cxx:166
 AliPWGHistoTools.cxx:167
 AliPWGHistoTools.cxx:168
 AliPWGHistoTools.cxx:169
 AliPWGHistoTools.cxx:170
 AliPWGHistoTools.cxx:171
 AliPWGHistoTools.cxx:172
 AliPWGHistoTools.cxx:173
 AliPWGHistoTools.cxx:174
 AliPWGHistoTools.cxx:175
 AliPWGHistoTools.cxx:176
 AliPWGHistoTools.cxx:177
 AliPWGHistoTools.cxx:178
 AliPWGHistoTools.cxx:179
 AliPWGHistoTools.cxx:180
 AliPWGHistoTools.cxx:181
 AliPWGHistoTools.cxx:182
 AliPWGHistoTools.cxx:183
 AliPWGHistoTools.cxx:184
 AliPWGHistoTools.cxx:185
 AliPWGHistoTools.cxx:186
 AliPWGHistoTools.cxx:187
 AliPWGHistoTools.cxx:188
 AliPWGHistoTools.cxx:189
 AliPWGHistoTools.cxx:190
 AliPWGHistoTools.cxx:191
 AliPWGHistoTools.cxx:192
 AliPWGHistoTools.cxx:193
 AliPWGHistoTools.cxx:194
 AliPWGHistoTools.cxx:195
 AliPWGHistoTools.cxx:196
 AliPWGHistoTools.cxx:197
 AliPWGHistoTools.cxx:198
 AliPWGHistoTools.cxx:199
 AliPWGHistoTools.cxx:200
 AliPWGHistoTools.cxx:201
 AliPWGHistoTools.cxx:202
 AliPWGHistoTools.cxx:203
 AliPWGHistoTools.cxx:204
 AliPWGHistoTools.cxx:205
 AliPWGHistoTools.cxx:206
 AliPWGHistoTools.cxx:207
 AliPWGHistoTools.cxx:208
 AliPWGHistoTools.cxx:209
 AliPWGHistoTools.cxx:210
 AliPWGHistoTools.cxx:211
 AliPWGHistoTools.cxx:212
 AliPWGHistoTools.cxx:213
 AliPWGHistoTools.cxx:214
 AliPWGHistoTools.cxx:215
 AliPWGHistoTools.cxx:216
 AliPWGHistoTools.cxx:217
 AliPWGHistoTools.cxx:218
 AliPWGHistoTools.cxx:219
 AliPWGHistoTools.cxx:220
 AliPWGHistoTools.cxx:221
 AliPWGHistoTools.cxx:222
 AliPWGHistoTools.cxx:223
 AliPWGHistoTools.cxx:224
 AliPWGHistoTools.cxx:225
 AliPWGHistoTools.cxx:226
 AliPWGHistoTools.cxx:227
 AliPWGHistoTools.cxx:228
 AliPWGHistoTools.cxx:229
 AliPWGHistoTools.cxx:230
 AliPWGHistoTools.cxx:231
 AliPWGHistoTools.cxx:232
 AliPWGHistoTools.cxx:233
 AliPWGHistoTools.cxx:234
 AliPWGHistoTools.cxx:235
 AliPWGHistoTools.cxx:236
 AliPWGHistoTools.cxx:237
 AliPWGHistoTools.cxx:238
 AliPWGHistoTools.cxx:239
 AliPWGHistoTools.cxx:240
 AliPWGHistoTools.cxx:241
 AliPWGHistoTools.cxx:242
 AliPWGHistoTools.cxx:243
 AliPWGHistoTools.cxx:244
 AliPWGHistoTools.cxx:245
 AliPWGHistoTools.cxx:246
 AliPWGHistoTools.cxx:247
 AliPWGHistoTools.cxx:248
 AliPWGHistoTools.cxx:249
 AliPWGHistoTools.cxx:250
 AliPWGHistoTools.cxx:251
 AliPWGHistoTools.cxx:252
 AliPWGHistoTools.cxx:253
 AliPWGHistoTools.cxx:254
 AliPWGHistoTools.cxx:255
 AliPWGHistoTools.cxx:256
 AliPWGHistoTools.cxx:257
 AliPWGHistoTools.cxx:258
 AliPWGHistoTools.cxx:259
 AliPWGHistoTools.cxx:260
 AliPWGHistoTools.cxx:261
 AliPWGHistoTools.cxx:262
 AliPWGHistoTools.cxx:263
 AliPWGHistoTools.cxx:264
 AliPWGHistoTools.cxx:265
 AliPWGHistoTools.cxx:266
 AliPWGHistoTools.cxx:267
 AliPWGHistoTools.cxx:268
 AliPWGHistoTools.cxx:269
 AliPWGHistoTools.cxx:270
 AliPWGHistoTools.cxx:271
 AliPWGHistoTools.cxx:272
 AliPWGHistoTools.cxx:273
 AliPWGHistoTools.cxx:274
 AliPWGHistoTools.cxx:275
 AliPWGHistoTools.cxx:276
 AliPWGHistoTools.cxx:277
 AliPWGHistoTools.cxx:278
 AliPWGHistoTools.cxx:279
 AliPWGHistoTools.cxx:280
 AliPWGHistoTools.cxx:281
 AliPWGHistoTools.cxx:282
 AliPWGHistoTools.cxx:283
 AliPWGHistoTools.cxx:284
 AliPWGHistoTools.cxx:285
 AliPWGHistoTools.cxx:286
 AliPWGHistoTools.cxx:287
 AliPWGHistoTools.cxx:288
 AliPWGHistoTools.cxx:289
 AliPWGHistoTools.cxx:290
 AliPWGHistoTools.cxx:291
 AliPWGHistoTools.cxx:292
 AliPWGHistoTools.cxx:293
 AliPWGHistoTools.cxx:294
 AliPWGHistoTools.cxx:295
 AliPWGHistoTools.cxx:296
 AliPWGHistoTools.cxx:297
 AliPWGHistoTools.cxx:298
 AliPWGHistoTools.cxx:299
 AliPWGHistoTools.cxx:300
 AliPWGHistoTools.cxx:301
 AliPWGHistoTools.cxx:302
 AliPWGHistoTools.cxx:303
 AliPWGHistoTools.cxx:304
 AliPWGHistoTools.cxx:305
 AliPWGHistoTools.cxx:306
 AliPWGHistoTools.cxx:307
 AliPWGHistoTools.cxx:308
 AliPWGHistoTools.cxx:309
 AliPWGHistoTools.cxx:310
 AliPWGHistoTools.cxx:311
 AliPWGHistoTools.cxx:312
 AliPWGHistoTools.cxx:313
 AliPWGHistoTools.cxx:314
 AliPWGHistoTools.cxx:315
 AliPWGHistoTools.cxx:316
 AliPWGHistoTools.cxx:317
 AliPWGHistoTools.cxx:318
 AliPWGHistoTools.cxx:319
 AliPWGHistoTools.cxx:320
 AliPWGHistoTools.cxx:321
 AliPWGHistoTools.cxx:322
 AliPWGHistoTools.cxx:323
 AliPWGHistoTools.cxx:324
 AliPWGHistoTools.cxx:325
 AliPWGHistoTools.cxx:326
 AliPWGHistoTools.cxx:327
 AliPWGHistoTools.cxx:328
 AliPWGHistoTools.cxx:329
 AliPWGHistoTools.cxx:330
 AliPWGHistoTools.cxx:331
 AliPWGHistoTools.cxx:332
 AliPWGHistoTools.cxx:333
 AliPWGHistoTools.cxx:334
 AliPWGHistoTools.cxx:335
 AliPWGHistoTools.cxx:336
 AliPWGHistoTools.cxx:337
 AliPWGHistoTools.cxx:338
 AliPWGHistoTools.cxx:339
 AliPWGHistoTools.cxx:340
 AliPWGHistoTools.cxx:341
 AliPWGHistoTools.cxx:342
 AliPWGHistoTools.cxx:343
 AliPWGHistoTools.cxx:344
 AliPWGHistoTools.cxx:345
 AliPWGHistoTools.cxx:346
 AliPWGHistoTools.cxx:347
 AliPWGHistoTools.cxx:348
 AliPWGHistoTools.cxx:349
 AliPWGHistoTools.cxx:350
 AliPWGHistoTools.cxx:351
 AliPWGHistoTools.cxx:352
 AliPWGHistoTools.cxx:353
 AliPWGHistoTools.cxx:354
 AliPWGHistoTools.cxx:355
 AliPWGHistoTools.cxx:356
 AliPWGHistoTools.cxx:357
 AliPWGHistoTools.cxx:358
 AliPWGHistoTools.cxx:359
 AliPWGHistoTools.cxx:360
 AliPWGHistoTools.cxx:361
 AliPWGHistoTools.cxx:362
 AliPWGHistoTools.cxx:363
 AliPWGHistoTools.cxx:364
 AliPWGHistoTools.cxx:365
 AliPWGHistoTools.cxx:366
 AliPWGHistoTools.cxx:367
 AliPWGHistoTools.cxx:368
 AliPWGHistoTools.cxx:369
 AliPWGHistoTools.cxx:370
 AliPWGHistoTools.cxx:371
 AliPWGHistoTools.cxx:372
 AliPWGHistoTools.cxx:373
 AliPWGHistoTools.cxx:374
 AliPWGHistoTools.cxx:375
 AliPWGHistoTools.cxx:376
 AliPWGHistoTools.cxx:377
 AliPWGHistoTools.cxx:378
 AliPWGHistoTools.cxx:379
 AliPWGHistoTools.cxx:380
 AliPWGHistoTools.cxx:381
 AliPWGHistoTools.cxx:382
 AliPWGHistoTools.cxx:383
 AliPWGHistoTools.cxx:384
 AliPWGHistoTools.cxx:385
 AliPWGHistoTools.cxx:386
 AliPWGHistoTools.cxx:387
 AliPWGHistoTools.cxx:388
 AliPWGHistoTools.cxx:389
 AliPWGHistoTools.cxx:390
 AliPWGHistoTools.cxx:391
 AliPWGHistoTools.cxx:392
 AliPWGHistoTools.cxx:393
 AliPWGHistoTools.cxx:394
 AliPWGHistoTools.cxx:395
 AliPWGHistoTools.cxx:396
 AliPWGHistoTools.cxx:397
 AliPWGHistoTools.cxx:398
 AliPWGHistoTools.cxx:399
 AliPWGHistoTools.cxx:400
 AliPWGHistoTools.cxx:401
 AliPWGHistoTools.cxx:402
 AliPWGHistoTools.cxx:403
 AliPWGHistoTools.cxx:404
 AliPWGHistoTools.cxx:405
 AliPWGHistoTools.cxx:406
 AliPWGHistoTools.cxx:407
 AliPWGHistoTools.cxx:408
 AliPWGHistoTools.cxx:409
 AliPWGHistoTools.cxx:410
 AliPWGHistoTools.cxx:411
 AliPWGHistoTools.cxx:412
 AliPWGHistoTools.cxx:413
 AliPWGHistoTools.cxx:414
 AliPWGHistoTools.cxx:415
 AliPWGHistoTools.cxx:416
 AliPWGHistoTools.cxx:417
 AliPWGHistoTools.cxx:418
 AliPWGHistoTools.cxx:419
 AliPWGHistoTools.cxx:420
 AliPWGHistoTools.cxx:421
 AliPWGHistoTools.cxx:422
 AliPWGHistoTools.cxx:423
 AliPWGHistoTools.cxx:424
 AliPWGHistoTools.cxx:425
 AliPWGHistoTools.cxx:426
 AliPWGHistoTools.cxx:427
 AliPWGHistoTools.cxx:428
 AliPWGHistoTools.cxx:429
 AliPWGHistoTools.cxx:430
 AliPWGHistoTools.cxx:431
 AliPWGHistoTools.cxx:432
 AliPWGHistoTools.cxx:433
 AliPWGHistoTools.cxx:434
 AliPWGHistoTools.cxx:435
 AliPWGHistoTools.cxx:436
 AliPWGHistoTools.cxx:437
 AliPWGHistoTools.cxx:438
 AliPWGHistoTools.cxx:439
 AliPWGHistoTools.cxx:440
 AliPWGHistoTools.cxx:441
 AliPWGHistoTools.cxx:442
 AliPWGHistoTools.cxx:443
 AliPWGHistoTools.cxx:444
 AliPWGHistoTools.cxx:445
 AliPWGHistoTools.cxx:446
 AliPWGHistoTools.cxx:447
 AliPWGHistoTools.cxx:448
 AliPWGHistoTools.cxx:449
 AliPWGHistoTools.cxx:450
 AliPWGHistoTools.cxx:451
 AliPWGHistoTools.cxx:452
 AliPWGHistoTools.cxx:453
 AliPWGHistoTools.cxx:454
 AliPWGHistoTools.cxx:455
 AliPWGHistoTools.cxx:456
 AliPWGHistoTools.cxx:457
 AliPWGHistoTools.cxx:458
 AliPWGHistoTools.cxx:459
 AliPWGHistoTools.cxx:460
 AliPWGHistoTools.cxx:461
 AliPWGHistoTools.cxx:462
 AliPWGHistoTools.cxx:463
 AliPWGHistoTools.cxx:464
 AliPWGHistoTools.cxx:465
 AliPWGHistoTools.cxx:466
 AliPWGHistoTools.cxx:467
 AliPWGHistoTools.cxx:468
 AliPWGHistoTools.cxx:469
 AliPWGHistoTools.cxx:470
 AliPWGHistoTools.cxx:471
 AliPWGHistoTools.cxx:472
 AliPWGHistoTools.cxx:473
 AliPWGHistoTools.cxx:474
 AliPWGHistoTools.cxx:475
 AliPWGHistoTools.cxx:476
 AliPWGHistoTools.cxx:477
 AliPWGHistoTools.cxx:478
 AliPWGHistoTools.cxx:479
 AliPWGHistoTools.cxx:480
 AliPWGHistoTools.cxx:481
 AliPWGHistoTools.cxx:482
 AliPWGHistoTools.cxx:483
 AliPWGHistoTools.cxx:484
 AliPWGHistoTools.cxx:485
 AliPWGHistoTools.cxx:486
 AliPWGHistoTools.cxx:487
 AliPWGHistoTools.cxx:488
 AliPWGHistoTools.cxx:489
 AliPWGHistoTools.cxx:490
 AliPWGHistoTools.cxx:491
 AliPWGHistoTools.cxx:492
 AliPWGHistoTools.cxx:493
 AliPWGHistoTools.cxx:494
 AliPWGHistoTools.cxx:495
 AliPWGHistoTools.cxx:496
 AliPWGHistoTools.cxx:497
 AliPWGHistoTools.cxx:498
 AliPWGHistoTools.cxx:499
 AliPWGHistoTools.cxx:500
 AliPWGHistoTools.cxx:501
 AliPWGHistoTools.cxx:502
 AliPWGHistoTools.cxx:503
 AliPWGHistoTools.cxx:504
 AliPWGHistoTools.cxx:505
 AliPWGHistoTools.cxx:506
 AliPWGHistoTools.cxx:507
 AliPWGHistoTools.cxx:508
 AliPWGHistoTools.cxx:509
 AliPWGHistoTools.cxx:510
 AliPWGHistoTools.cxx:511
 AliPWGHistoTools.cxx:512
 AliPWGHistoTools.cxx:513
 AliPWGHistoTools.cxx:514
 AliPWGHistoTools.cxx:515
 AliPWGHistoTools.cxx:516
 AliPWGHistoTools.cxx:517
 AliPWGHistoTools.cxx:518
 AliPWGHistoTools.cxx:519
 AliPWGHistoTools.cxx:520
 AliPWGHistoTools.cxx:521
 AliPWGHistoTools.cxx:522
 AliPWGHistoTools.cxx:523
 AliPWGHistoTools.cxx:524
 AliPWGHistoTools.cxx:525
 AliPWGHistoTools.cxx:526
 AliPWGHistoTools.cxx:527
 AliPWGHistoTools.cxx:528
 AliPWGHistoTools.cxx:529
 AliPWGHistoTools.cxx:530
 AliPWGHistoTools.cxx:531
 AliPWGHistoTools.cxx:532
 AliPWGHistoTools.cxx:533
 AliPWGHistoTools.cxx:534
 AliPWGHistoTools.cxx:535
 AliPWGHistoTools.cxx:536
 AliPWGHistoTools.cxx:537
 AliPWGHistoTools.cxx:538
 AliPWGHistoTools.cxx:539
 AliPWGHistoTools.cxx:540
 AliPWGHistoTools.cxx:541
 AliPWGHistoTools.cxx:542
 AliPWGHistoTools.cxx:543
 AliPWGHistoTools.cxx:544
 AliPWGHistoTools.cxx:545
 AliPWGHistoTools.cxx:546
 AliPWGHistoTools.cxx:547
 AliPWGHistoTools.cxx:548
 AliPWGHistoTools.cxx:549
 AliPWGHistoTools.cxx:550
 AliPWGHistoTools.cxx:551
 AliPWGHistoTools.cxx:552
 AliPWGHistoTools.cxx:553
 AliPWGHistoTools.cxx:554
 AliPWGHistoTools.cxx:555
 AliPWGHistoTools.cxx:556
 AliPWGHistoTools.cxx:557
 AliPWGHistoTools.cxx:558
 AliPWGHistoTools.cxx:559
 AliPWGHistoTools.cxx:560
 AliPWGHistoTools.cxx:561
 AliPWGHistoTools.cxx:562
 AliPWGHistoTools.cxx:563
 AliPWGHistoTools.cxx:564
 AliPWGHistoTools.cxx:565
 AliPWGHistoTools.cxx:566
 AliPWGHistoTools.cxx:567
 AliPWGHistoTools.cxx:568
 AliPWGHistoTools.cxx:569
 AliPWGHistoTools.cxx:570
 AliPWGHistoTools.cxx:571
 AliPWGHistoTools.cxx:572
 AliPWGHistoTools.cxx:573
 AliPWGHistoTools.cxx:574
 AliPWGHistoTools.cxx:575
 AliPWGHistoTools.cxx:576
 AliPWGHistoTools.cxx:577
 AliPWGHistoTools.cxx:578
 AliPWGHistoTools.cxx:579
 AliPWGHistoTools.cxx:580
 AliPWGHistoTools.cxx:581
 AliPWGHistoTools.cxx:582
 AliPWGHistoTools.cxx:583
 AliPWGHistoTools.cxx:584
 AliPWGHistoTools.cxx:585
 AliPWGHistoTools.cxx:586
 AliPWGHistoTools.cxx:587
 AliPWGHistoTools.cxx:588
 AliPWGHistoTools.cxx:589
 AliPWGHistoTools.cxx:590
 AliPWGHistoTools.cxx:591
 AliPWGHistoTools.cxx:592
 AliPWGHistoTools.cxx:593
 AliPWGHistoTools.cxx:594
 AliPWGHistoTools.cxx:595
 AliPWGHistoTools.cxx:596
 AliPWGHistoTools.cxx:597
 AliPWGHistoTools.cxx:598
 AliPWGHistoTools.cxx:599
 AliPWGHistoTools.cxx:600
 AliPWGHistoTools.cxx:601
 AliPWGHistoTools.cxx:602
 AliPWGHistoTools.cxx:603
 AliPWGHistoTools.cxx:604
 AliPWGHistoTools.cxx:605
 AliPWGHistoTools.cxx:606
 AliPWGHistoTools.cxx:607
 AliPWGHistoTools.cxx:608
 AliPWGHistoTools.cxx:609
 AliPWGHistoTools.cxx:610
 AliPWGHistoTools.cxx:611
 AliPWGHistoTools.cxx:612
 AliPWGHistoTools.cxx:613
 AliPWGHistoTools.cxx:614
 AliPWGHistoTools.cxx:615
 AliPWGHistoTools.cxx:616
 AliPWGHistoTools.cxx:617
 AliPWGHistoTools.cxx:618
 AliPWGHistoTools.cxx:619
 AliPWGHistoTools.cxx:620
 AliPWGHistoTools.cxx:621
 AliPWGHistoTools.cxx:622
 AliPWGHistoTools.cxx:623
 AliPWGHistoTools.cxx:624
 AliPWGHistoTools.cxx:625
 AliPWGHistoTools.cxx:626
 AliPWGHistoTools.cxx:627
 AliPWGHistoTools.cxx:628
 AliPWGHistoTools.cxx:629
 AliPWGHistoTools.cxx:630
 AliPWGHistoTools.cxx:631
 AliPWGHistoTools.cxx:632
 AliPWGHistoTools.cxx:633
 AliPWGHistoTools.cxx:634
 AliPWGHistoTools.cxx:635
 AliPWGHistoTools.cxx:636
 AliPWGHistoTools.cxx:637
 AliPWGHistoTools.cxx:638
 AliPWGHistoTools.cxx:639
 AliPWGHistoTools.cxx:640
 AliPWGHistoTools.cxx:641
 AliPWGHistoTools.cxx:642
 AliPWGHistoTools.cxx:643
 AliPWGHistoTools.cxx:644
 AliPWGHistoTools.cxx:645
 AliPWGHistoTools.cxx:646
 AliPWGHistoTools.cxx:647
 AliPWGHistoTools.cxx:648
 AliPWGHistoTools.cxx:649
 AliPWGHistoTools.cxx:650
 AliPWGHistoTools.cxx:651
 AliPWGHistoTools.cxx:652
 AliPWGHistoTools.cxx:653
 AliPWGHistoTools.cxx:654
 AliPWGHistoTools.cxx:655
 AliPWGHistoTools.cxx:656
 AliPWGHistoTools.cxx:657
 AliPWGHistoTools.cxx:658
 AliPWGHistoTools.cxx:659
 AliPWGHistoTools.cxx:660
 AliPWGHistoTools.cxx:661
 AliPWGHistoTools.cxx:662
 AliPWGHistoTools.cxx:663
 AliPWGHistoTools.cxx:664
 AliPWGHistoTools.cxx:665
 AliPWGHistoTools.cxx:666
 AliPWGHistoTools.cxx:667
 AliPWGHistoTools.cxx:668
 AliPWGHistoTools.cxx:669
 AliPWGHistoTools.cxx:670
 AliPWGHistoTools.cxx:671
 AliPWGHistoTools.cxx:672
 AliPWGHistoTools.cxx:673
 AliPWGHistoTools.cxx:674
 AliPWGHistoTools.cxx:675
 AliPWGHistoTools.cxx:676
 AliPWGHistoTools.cxx:677
 AliPWGHistoTools.cxx:678
 AliPWGHistoTools.cxx:679
 AliPWGHistoTools.cxx:680
 AliPWGHistoTools.cxx:681
 AliPWGHistoTools.cxx:682
 AliPWGHistoTools.cxx:683
 AliPWGHistoTools.cxx:684
 AliPWGHistoTools.cxx:685
 AliPWGHistoTools.cxx:686
 AliPWGHistoTools.cxx:687
 AliPWGHistoTools.cxx:688
 AliPWGHistoTools.cxx:689
 AliPWGHistoTools.cxx:690
 AliPWGHistoTools.cxx:691
 AliPWGHistoTools.cxx:692
 AliPWGHistoTools.cxx:693
 AliPWGHistoTools.cxx:694
 AliPWGHistoTools.cxx:695
 AliPWGHistoTools.cxx:696
 AliPWGHistoTools.cxx:697
 AliPWGHistoTools.cxx:698
 AliPWGHistoTools.cxx:699
 AliPWGHistoTools.cxx:700
 AliPWGHistoTools.cxx:701
 AliPWGHistoTools.cxx:702
 AliPWGHistoTools.cxx:703
 AliPWGHistoTools.cxx:704
 AliPWGHistoTools.cxx:705
 AliPWGHistoTools.cxx:706
 AliPWGHistoTools.cxx:707
 AliPWGHistoTools.cxx:708
 AliPWGHistoTools.cxx:709
 AliPWGHistoTools.cxx:710
 AliPWGHistoTools.cxx:711
 AliPWGHistoTools.cxx:712
 AliPWGHistoTools.cxx:713
 AliPWGHistoTools.cxx:714
 AliPWGHistoTools.cxx:715
 AliPWGHistoTools.cxx:716
 AliPWGHistoTools.cxx:717
 AliPWGHistoTools.cxx:718
 AliPWGHistoTools.cxx:719
 AliPWGHistoTools.cxx:720
 AliPWGHistoTools.cxx:721
 AliPWGHistoTools.cxx:722
 AliPWGHistoTools.cxx:723
 AliPWGHistoTools.cxx:724
 AliPWGHistoTools.cxx:725
 AliPWGHistoTools.cxx:726
 AliPWGHistoTools.cxx:727
 AliPWGHistoTools.cxx:728
 AliPWGHistoTools.cxx:729
 AliPWGHistoTools.cxx:730
 AliPWGHistoTools.cxx:731
 AliPWGHistoTools.cxx:732
 AliPWGHistoTools.cxx:733
 AliPWGHistoTools.cxx:734
 AliPWGHistoTools.cxx:735
 AliPWGHistoTools.cxx:736
 AliPWGHistoTools.cxx:737
 AliPWGHistoTools.cxx:738
 AliPWGHistoTools.cxx:739
 AliPWGHistoTools.cxx:740
 AliPWGHistoTools.cxx:741
 AliPWGHistoTools.cxx:742
 AliPWGHistoTools.cxx:743
 AliPWGHistoTools.cxx:744
 AliPWGHistoTools.cxx:745
 AliPWGHistoTools.cxx:746
 AliPWGHistoTools.cxx:747
 AliPWGHistoTools.cxx:748
 AliPWGHistoTools.cxx:749
 AliPWGHistoTools.cxx:750
 AliPWGHistoTools.cxx:751
 AliPWGHistoTools.cxx:752
 AliPWGHistoTools.cxx:753
 AliPWGHistoTools.cxx:754
 AliPWGHistoTools.cxx:755
 AliPWGHistoTools.cxx:756
 AliPWGHistoTools.cxx:757
 AliPWGHistoTools.cxx:758
 AliPWGHistoTools.cxx:759
 AliPWGHistoTools.cxx:760
 AliPWGHistoTools.cxx:761
 AliPWGHistoTools.cxx:762
 AliPWGHistoTools.cxx:763
 AliPWGHistoTools.cxx:764
 AliPWGHistoTools.cxx:765
 AliPWGHistoTools.cxx:766
 AliPWGHistoTools.cxx:767
 AliPWGHistoTools.cxx:768
 AliPWGHistoTools.cxx:769
 AliPWGHistoTools.cxx:770
 AliPWGHistoTools.cxx:771
 AliPWGHistoTools.cxx:772
 AliPWGHistoTools.cxx:773
 AliPWGHistoTools.cxx:774
 AliPWGHistoTools.cxx:775
 AliPWGHistoTools.cxx:776
 AliPWGHistoTools.cxx:777
 AliPWGHistoTools.cxx:778
 AliPWGHistoTools.cxx:779
 AliPWGHistoTools.cxx:780
 AliPWGHistoTools.cxx:781
 AliPWGHistoTools.cxx:782
 AliPWGHistoTools.cxx:783
 AliPWGHistoTools.cxx:784
 AliPWGHistoTools.cxx:785
 AliPWGHistoTools.cxx:786
 AliPWGHistoTools.cxx:787
 AliPWGHistoTools.cxx:788
 AliPWGHistoTools.cxx:789
 AliPWGHistoTools.cxx:790
 AliPWGHistoTools.cxx:791
 AliPWGHistoTools.cxx:792
 AliPWGHistoTools.cxx:793
 AliPWGHistoTools.cxx:794
 AliPWGHistoTools.cxx:795
 AliPWGHistoTools.cxx:796
 AliPWGHistoTools.cxx:797
 AliPWGHistoTools.cxx:798
 AliPWGHistoTools.cxx:799
 AliPWGHistoTools.cxx:800
 AliPWGHistoTools.cxx:801
 AliPWGHistoTools.cxx:802
 AliPWGHistoTools.cxx:803
 AliPWGHistoTools.cxx:804
 AliPWGHistoTools.cxx:805
 AliPWGHistoTools.cxx:806
 AliPWGHistoTools.cxx:807
 AliPWGHistoTools.cxx:808
 AliPWGHistoTools.cxx:809
 AliPWGHistoTools.cxx:810
 AliPWGHistoTools.cxx:811
 AliPWGHistoTools.cxx:812
 AliPWGHistoTools.cxx:813
 AliPWGHistoTools.cxx:814
 AliPWGHistoTools.cxx:815
 AliPWGHistoTools.cxx:816
 AliPWGHistoTools.cxx:817
 AliPWGHistoTools.cxx:818
 AliPWGHistoTools.cxx:819
 AliPWGHistoTools.cxx:820
 AliPWGHistoTools.cxx:821
 AliPWGHistoTools.cxx:822
 AliPWGHistoTools.cxx:823
 AliPWGHistoTools.cxx:824
 AliPWGHistoTools.cxx:825
 AliPWGHistoTools.cxx:826
 AliPWGHistoTools.cxx:827
 AliPWGHistoTools.cxx:828
 AliPWGHistoTools.cxx:829
 AliPWGHistoTools.cxx:830
 AliPWGHistoTools.cxx:831
 AliPWGHistoTools.cxx:832
 AliPWGHistoTools.cxx:833
 AliPWGHistoTools.cxx:834
 AliPWGHistoTools.cxx:835
 AliPWGHistoTools.cxx:836
 AliPWGHistoTools.cxx:837
 AliPWGHistoTools.cxx:838
 AliPWGHistoTools.cxx:839
 AliPWGHistoTools.cxx:840
 AliPWGHistoTools.cxx:841
 AliPWGHistoTools.cxx:842
 AliPWGHistoTools.cxx:843
 AliPWGHistoTools.cxx:844
 AliPWGHistoTools.cxx:845
 AliPWGHistoTools.cxx:846
 AliPWGHistoTools.cxx:847
 AliPWGHistoTools.cxx:848
 AliPWGHistoTools.cxx:849
 AliPWGHistoTools.cxx:850
 AliPWGHistoTools.cxx:851
 AliPWGHistoTools.cxx:852
 AliPWGHistoTools.cxx:853
 AliPWGHistoTools.cxx:854
 AliPWGHistoTools.cxx:855
 AliPWGHistoTools.cxx:856
 AliPWGHistoTools.cxx:857
 AliPWGHistoTools.cxx:858
 AliPWGHistoTools.cxx:859
 AliPWGHistoTools.cxx:860
 AliPWGHistoTools.cxx:861
 AliPWGHistoTools.cxx:862
 AliPWGHistoTools.cxx:863
 AliPWGHistoTools.cxx:864
 AliPWGHistoTools.cxx:865
 AliPWGHistoTools.cxx:866
 AliPWGHistoTools.cxx:867
 AliPWGHistoTools.cxx:868
 AliPWGHistoTools.cxx:869
 AliPWGHistoTools.cxx:870
 AliPWGHistoTools.cxx:871
 AliPWGHistoTools.cxx:872
 AliPWGHistoTools.cxx:873
 AliPWGHistoTools.cxx:874
 AliPWGHistoTools.cxx:875
 AliPWGHistoTools.cxx:876
 AliPWGHistoTools.cxx:877
 AliPWGHistoTools.cxx:878
 AliPWGHistoTools.cxx:879
 AliPWGHistoTools.cxx:880
 AliPWGHistoTools.cxx:881
 AliPWGHistoTools.cxx:882
 AliPWGHistoTools.cxx:883
 AliPWGHistoTools.cxx:884
 AliPWGHistoTools.cxx:885
 AliPWGHistoTools.cxx:886
 AliPWGHistoTools.cxx:887
 AliPWGHistoTools.cxx:888
 AliPWGHistoTools.cxx:889
 AliPWGHistoTools.cxx:890
 AliPWGHistoTools.cxx:891
 AliPWGHistoTools.cxx:892
 AliPWGHistoTools.cxx:893
 AliPWGHistoTools.cxx:894
 AliPWGHistoTools.cxx:895
 AliPWGHistoTools.cxx:896
 AliPWGHistoTools.cxx:897
 AliPWGHistoTools.cxx:898
 AliPWGHistoTools.cxx:899
 AliPWGHistoTools.cxx:900
 AliPWGHistoTools.cxx:901
 AliPWGHistoTools.cxx:902
 AliPWGHistoTools.cxx:903
 AliPWGHistoTools.cxx:904
 AliPWGHistoTools.cxx:905
 AliPWGHistoTools.cxx:906
 AliPWGHistoTools.cxx:907
 AliPWGHistoTools.cxx:908
 AliPWGHistoTools.cxx:909
 AliPWGHistoTools.cxx:910
 AliPWGHistoTools.cxx:911
 AliPWGHistoTools.cxx:912
 AliPWGHistoTools.cxx:913
 AliPWGHistoTools.cxx:914
 AliPWGHistoTools.cxx:915
 AliPWGHistoTools.cxx:916
 AliPWGHistoTools.cxx:917
 AliPWGHistoTools.cxx:918
 AliPWGHistoTools.cxx:919
 AliPWGHistoTools.cxx:920
 AliPWGHistoTools.cxx:921
 AliPWGHistoTools.cxx:922
 AliPWGHistoTools.cxx:923
 AliPWGHistoTools.cxx:924
 AliPWGHistoTools.cxx:925
 AliPWGHistoTools.cxx:926
 AliPWGHistoTools.cxx:927
 AliPWGHistoTools.cxx:928
 AliPWGHistoTools.cxx:929
 AliPWGHistoTools.cxx:930
 AliPWGHistoTools.cxx:931
 AliPWGHistoTools.cxx:932
 AliPWGHistoTools.cxx:933
 AliPWGHistoTools.cxx:934
 AliPWGHistoTools.cxx:935
 AliPWGHistoTools.cxx:936
 AliPWGHistoTools.cxx:937
 AliPWGHistoTools.cxx:938
 AliPWGHistoTools.cxx:939
 AliPWGHistoTools.cxx:940
 AliPWGHistoTools.cxx:941
 AliPWGHistoTools.cxx:942
 AliPWGHistoTools.cxx:943
 AliPWGHistoTools.cxx:944
 AliPWGHistoTools.cxx:945
 AliPWGHistoTools.cxx:946
 AliPWGHistoTools.cxx:947
 AliPWGHistoTools.cxx:948
 AliPWGHistoTools.cxx:949
 AliPWGHistoTools.cxx:950
 AliPWGHistoTools.cxx:951
 AliPWGHistoTools.cxx:952
 AliPWGHistoTools.cxx:953
 AliPWGHistoTools.cxx:954
 AliPWGHistoTools.cxx:955
 AliPWGHistoTools.cxx:956
 AliPWGHistoTools.cxx:957
 AliPWGHistoTools.cxx:958
 AliPWGHistoTools.cxx:959
 AliPWGHistoTools.cxx:960
 AliPWGHistoTools.cxx:961
 AliPWGHistoTools.cxx:962
 AliPWGHistoTools.cxx:963
 AliPWGHistoTools.cxx:964
 AliPWGHistoTools.cxx:965
 AliPWGHistoTools.cxx:966
 AliPWGHistoTools.cxx:967
 AliPWGHistoTools.cxx:968
 AliPWGHistoTools.cxx:969
 AliPWGHistoTools.cxx:970
 AliPWGHistoTools.cxx:971
 AliPWGHistoTools.cxx:972
 AliPWGHistoTools.cxx:973
 AliPWGHistoTools.cxx:974
 AliPWGHistoTools.cxx:975
 AliPWGHistoTools.cxx:976
 AliPWGHistoTools.cxx:977
 AliPWGHistoTools.cxx:978
 AliPWGHistoTools.cxx:979
 AliPWGHistoTools.cxx:980
 AliPWGHistoTools.cxx:981
 AliPWGHistoTools.cxx:982
 AliPWGHistoTools.cxx:983
 AliPWGHistoTools.cxx:984
 AliPWGHistoTools.cxx:985
 AliPWGHistoTools.cxx:986
 AliPWGHistoTools.cxx:987
 AliPWGHistoTools.cxx:988
 AliPWGHistoTools.cxx:989
 AliPWGHistoTools.cxx:990
 AliPWGHistoTools.cxx:991
 AliPWGHistoTools.cxx:992
 AliPWGHistoTools.cxx:993
 AliPWGHistoTools.cxx:994
 AliPWGHistoTools.cxx:995
 AliPWGHistoTools.cxx:996
 AliPWGHistoTools.cxx:997
 AliPWGHistoTools.cxx:998
 AliPWGHistoTools.cxx:999
 AliPWGHistoTools.cxx:1000
 AliPWGHistoTools.cxx:1001
 AliPWGHistoTools.cxx:1002
 AliPWGHistoTools.cxx:1003
 AliPWGHistoTools.cxx:1004
 AliPWGHistoTools.cxx:1005
 AliPWGHistoTools.cxx:1006
 AliPWGHistoTools.cxx:1007
 AliPWGHistoTools.cxx:1008
 AliPWGHistoTools.cxx:1009
 AliPWGHistoTools.cxx:1010
 AliPWGHistoTools.cxx:1011
 AliPWGHistoTools.cxx:1012
 AliPWGHistoTools.cxx:1013
 AliPWGHistoTools.cxx:1014
 AliPWGHistoTools.cxx:1015
 AliPWGHistoTools.cxx:1016
 AliPWGHistoTools.cxx:1017
 AliPWGHistoTools.cxx:1018
 AliPWGHistoTools.cxx:1019
 AliPWGHistoTools.cxx:1020
 AliPWGHistoTools.cxx:1021
 AliPWGHistoTools.cxx:1022
 AliPWGHistoTools.cxx:1023
 AliPWGHistoTools.cxx:1024
 AliPWGHistoTools.cxx:1025
 AliPWGHistoTools.cxx:1026
 AliPWGHistoTools.cxx:1027
 AliPWGHistoTools.cxx:1028
 AliPWGHistoTools.cxx:1029
 AliPWGHistoTools.cxx:1030
 AliPWGHistoTools.cxx:1031
 AliPWGHistoTools.cxx:1032
 AliPWGHistoTools.cxx:1033
 AliPWGHistoTools.cxx:1034
 AliPWGHistoTools.cxx:1035
 AliPWGHistoTools.cxx:1036
 AliPWGHistoTools.cxx:1037
 AliPWGHistoTools.cxx:1038
 AliPWGHistoTools.cxx:1039
 AliPWGHistoTools.cxx:1040
 AliPWGHistoTools.cxx:1041
 AliPWGHistoTools.cxx:1042
 AliPWGHistoTools.cxx:1043
 AliPWGHistoTools.cxx:1044
 AliPWGHistoTools.cxx:1045
 AliPWGHistoTools.cxx:1046
 AliPWGHistoTools.cxx:1047
 AliPWGHistoTools.cxx:1048
 AliPWGHistoTools.cxx:1049
 AliPWGHistoTools.cxx:1050
 AliPWGHistoTools.cxx:1051
 AliPWGHistoTools.cxx:1052
 AliPWGHistoTools.cxx:1053
 AliPWGHistoTools.cxx:1054
 AliPWGHistoTools.cxx:1055
 AliPWGHistoTools.cxx:1056
 AliPWGHistoTools.cxx:1057
 AliPWGHistoTools.cxx:1058
 AliPWGHistoTools.cxx:1059
 AliPWGHistoTools.cxx:1060
 AliPWGHistoTools.cxx:1061
 AliPWGHistoTools.cxx:1062
 AliPWGHistoTools.cxx:1063
 AliPWGHistoTools.cxx:1064
 AliPWGHistoTools.cxx:1065
 AliPWGHistoTools.cxx:1066
 AliPWGHistoTools.cxx:1067
 AliPWGHistoTools.cxx:1068
 AliPWGHistoTools.cxx:1069
 AliPWGHistoTools.cxx:1070
 AliPWGHistoTools.cxx:1071
 AliPWGHistoTools.cxx:1072
 AliPWGHistoTools.cxx:1073
 AliPWGHistoTools.cxx:1074
 AliPWGHistoTools.cxx:1075
 AliPWGHistoTools.cxx:1076
 AliPWGHistoTools.cxx:1077
 AliPWGHistoTools.cxx:1078
 AliPWGHistoTools.cxx:1079
 AliPWGHistoTools.cxx:1080
 AliPWGHistoTools.cxx:1081
 AliPWGHistoTools.cxx:1082
 AliPWGHistoTools.cxx:1083
 AliPWGHistoTools.cxx:1084
 AliPWGHistoTools.cxx:1085
 AliPWGHistoTools.cxx:1086
 AliPWGHistoTools.cxx:1087
 AliPWGHistoTools.cxx:1088
 AliPWGHistoTools.cxx:1089
 AliPWGHistoTools.cxx:1090
 AliPWGHistoTools.cxx:1091
 AliPWGHistoTools.cxx:1092
 AliPWGHistoTools.cxx:1093
 AliPWGHistoTools.cxx:1094
 AliPWGHistoTools.cxx:1095
 AliPWGHistoTools.cxx:1096
 AliPWGHistoTools.cxx:1097
 AliPWGHistoTools.cxx:1098
 AliPWGHistoTools.cxx:1099
 AliPWGHistoTools.cxx:1100
 AliPWGHistoTools.cxx:1101
 AliPWGHistoTools.cxx:1102
 AliPWGHistoTools.cxx:1103
 AliPWGHistoTools.cxx:1104
 AliPWGHistoTools.cxx:1105
 AliPWGHistoTools.cxx:1106
 AliPWGHistoTools.cxx:1107
 AliPWGHistoTools.cxx:1108
 AliPWGHistoTools.cxx:1109
 AliPWGHistoTools.cxx:1110
 AliPWGHistoTools.cxx:1111
 AliPWGHistoTools.cxx:1112
 AliPWGHistoTools.cxx:1113
 AliPWGHistoTools.cxx:1114
 AliPWGHistoTools.cxx:1115
 AliPWGHistoTools.cxx:1116
 AliPWGHistoTools.cxx:1117
 AliPWGHistoTools.cxx:1118
 AliPWGHistoTools.cxx:1119
 AliPWGHistoTools.cxx:1120
 AliPWGHistoTools.cxx:1121
 AliPWGHistoTools.cxx:1122
 AliPWGHistoTools.cxx:1123
 AliPWGHistoTools.cxx:1124
 AliPWGHistoTools.cxx:1125
 AliPWGHistoTools.cxx:1126
 AliPWGHistoTools.cxx:1127
 AliPWGHistoTools.cxx:1128
 AliPWGHistoTools.cxx:1129
 AliPWGHistoTools.cxx:1130
 AliPWGHistoTools.cxx:1131
 AliPWGHistoTools.cxx:1132
 AliPWGHistoTools.cxx:1133
 AliPWGHistoTools.cxx:1134
 AliPWGHistoTools.cxx:1135
 AliPWGHistoTools.cxx:1136
 AliPWGHistoTools.cxx:1137
 AliPWGHistoTools.cxx:1138
 AliPWGHistoTools.cxx:1139
 AliPWGHistoTools.cxx:1140
 AliPWGHistoTools.cxx:1141
 AliPWGHistoTools.cxx:1142
 AliPWGHistoTools.cxx:1143
 AliPWGHistoTools.cxx:1144
 AliPWGHistoTools.cxx:1145
 AliPWGHistoTools.cxx:1146
 AliPWGHistoTools.cxx:1147
 AliPWGHistoTools.cxx:1148
 AliPWGHistoTools.cxx:1149
 AliPWGHistoTools.cxx:1150
 AliPWGHistoTools.cxx:1151
 AliPWGHistoTools.cxx:1152
 AliPWGHistoTools.cxx:1153
 AliPWGHistoTools.cxx:1154
 AliPWGHistoTools.cxx:1155
 AliPWGHistoTools.cxx:1156
 AliPWGHistoTools.cxx:1157
 AliPWGHistoTools.cxx:1158
 AliPWGHistoTools.cxx:1159
 AliPWGHistoTools.cxx:1160
 AliPWGHistoTools.cxx:1161
 AliPWGHistoTools.cxx:1162
 AliPWGHistoTools.cxx:1163
 AliPWGHistoTools.cxx:1164
 AliPWGHistoTools.cxx:1165
 AliPWGHistoTools.cxx:1166
 AliPWGHistoTools.cxx:1167
 AliPWGHistoTools.cxx:1168
 AliPWGHistoTools.cxx:1169
 AliPWGHistoTools.cxx:1170
 AliPWGHistoTools.cxx:1171
 AliPWGHistoTools.cxx:1172
 AliPWGHistoTools.cxx:1173
 AliPWGHistoTools.cxx:1174
 AliPWGHistoTools.cxx:1175
 AliPWGHistoTools.cxx:1176
 AliPWGHistoTools.cxx:1177
 AliPWGHistoTools.cxx:1178
 AliPWGHistoTools.cxx:1179
 AliPWGHistoTools.cxx:1180
 AliPWGHistoTools.cxx:1181
 AliPWGHistoTools.cxx:1182
 AliPWGHistoTools.cxx:1183
 AliPWGHistoTools.cxx:1184
 AliPWGHistoTools.cxx:1185
 AliPWGHistoTools.cxx:1186
 AliPWGHistoTools.cxx:1187
 AliPWGHistoTools.cxx:1188
 AliPWGHistoTools.cxx:1189
 AliPWGHistoTools.cxx:1190
 AliPWGHistoTools.cxx:1191
 AliPWGHistoTools.cxx:1192
 AliPWGHistoTools.cxx:1193
 AliPWGHistoTools.cxx:1194
 AliPWGHistoTools.cxx:1195
 AliPWGHistoTools.cxx:1196
 AliPWGHistoTools.cxx:1197
 AliPWGHistoTools.cxx:1198
 AliPWGHistoTools.cxx:1199
 AliPWGHistoTools.cxx:1200
 AliPWGHistoTools.cxx:1201
 AliPWGHistoTools.cxx:1202
 AliPWGHistoTools.cxx:1203
 AliPWGHistoTools.cxx:1204
 AliPWGHistoTools.cxx:1205
 AliPWGHistoTools.cxx:1206
 AliPWGHistoTools.cxx:1207
 AliPWGHistoTools.cxx:1208
 AliPWGHistoTools.cxx:1209
 AliPWGHistoTools.cxx:1210
 AliPWGHistoTools.cxx:1211
 AliPWGHistoTools.cxx:1212
 AliPWGHistoTools.cxx:1213
 AliPWGHistoTools.cxx:1214
 AliPWGHistoTools.cxx:1215
 AliPWGHistoTools.cxx:1216
 AliPWGHistoTools.cxx:1217
 AliPWGHistoTools.cxx:1218
 AliPWGHistoTools.cxx:1219
 AliPWGHistoTools.cxx:1220
 AliPWGHistoTools.cxx:1221
 AliPWGHistoTools.cxx:1222
 AliPWGHistoTools.cxx:1223
 AliPWGHistoTools.cxx:1224
 AliPWGHistoTools.cxx:1225
 AliPWGHistoTools.cxx:1226
 AliPWGHistoTools.cxx:1227
 AliPWGHistoTools.cxx:1228
 AliPWGHistoTools.cxx:1229
 AliPWGHistoTools.cxx:1230
 AliPWGHistoTools.cxx:1231
 AliPWGHistoTools.cxx:1232
 AliPWGHistoTools.cxx:1233
 AliPWGHistoTools.cxx:1234
 AliPWGHistoTools.cxx:1235
 AliPWGHistoTools.cxx:1236
 AliPWGHistoTools.cxx:1237
 AliPWGHistoTools.cxx:1238
 AliPWGHistoTools.cxx:1239
 AliPWGHistoTools.cxx:1240
 AliPWGHistoTools.cxx:1241
 AliPWGHistoTools.cxx:1242
 AliPWGHistoTools.cxx:1243
 AliPWGHistoTools.cxx:1244
 AliPWGHistoTools.cxx:1245
 AliPWGHistoTools.cxx:1246
 AliPWGHistoTools.cxx:1247
 AliPWGHistoTools.cxx:1248
 AliPWGHistoTools.cxx:1249
 AliPWGHistoTools.cxx:1250
 AliPWGHistoTools.cxx:1251
 AliPWGHistoTools.cxx:1252
 AliPWGHistoTools.cxx:1253
 AliPWGHistoTools.cxx:1254
 AliPWGHistoTools.cxx:1255
 AliPWGHistoTools.cxx:1256
 AliPWGHistoTools.cxx:1257
 AliPWGHistoTools.cxx:1258
 AliPWGHistoTools.cxx:1259
 AliPWGHistoTools.cxx:1260
 AliPWGHistoTools.cxx:1261
 AliPWGHistoTools.cxx:1262
 AliPWGHistoTools.cxx:1263
 AliPWGHistoTools.cxx:1264
 AliPWGHistoTools.cxx:1265
 AliPWGHistoTools.cxx:1266
 AliPWGHistoTools.cxx:1267
 AliPWGHistoTools.cxx:1268
 AliPWGHistoTools.cxx:1269
 AliPWGHistoTools.cxx:1270
 AliPWGHistoTools.cxx:1271
 AliPWGHistoTools.cxx:1272
 AliPWGHistoTools.cxx:1273
 AliPWGHistoTools.cxx:1274
 AliPWGHistoTools.cxx:1275
 AliPWGHistoTools.cxx:1276
 AliPWGHistoTools.cxx:1277
 AliPWGHistoTools.cxx:1278
 AliPWGHistoTools.cxx:1279
 AliPWGHistoTools.cxx:1280
 AliPWGHistoTools.cxx:1281
 AliPWGHistoTools.cxx:1282
 AliPWGHistoTools.cxx:1283
 AliPWGHistoTools.cxx:1284
 AliPWGHistoTools.cxx:1285
 AliPWGHistoTools.cxx:1286
 AliPWGHistoTools.cxx:1287
 AliPWGHistoTools.cxx:1288
 AliPWGHistoTools.cxx:1289
 AliPWGHistoTools.cxx:1290
 AliPWGHistoTools.cxx:1291
 AliPWGHistoTools.cxx:1292
 AliPWGHistoTools.cxx:1293
 AliPWGHistoTools.cxx:1294
 AliPWGHistoTools.cxx:1295
 AliPWGHistoTools.cxx:1296
 AliPWGHistoTools.cxx:1297
 AliPWGHistoTools.cxx:1298
 AliPWGHistoTools.cxx:1299
 AliPWGHistoTools.cxx:1300
 AliPWGHistoTools.cxx:1301
 AliPWGHistoTools.cxx:1302
 AliPWGHistoTools.cxx:1303
 AliPWGHistoTools.cxx:1304
 AliPWGHistoTools.cxx:1305
 AliPWGHistoTools.cxx:1306
 AliPWGHistoTools.cxx:1307
 AliPWGHistoTools.cxx:1308
 AliPWGHistoTools.cxx:1309
 AliPWGHistoTools.cxx:1310
 AliPWGHistoTools.cxx:1311
 AliPWGHistoTools.cxx:1312
 AliPWGHistoTools.cxx:1313
 AliPWGHistoTools.cxx:1314
 AliPWGHistoTools.cxx:1315
 AliPWGHistoTools.cxx:1316
 AliPWGHistoTools.cxx:1317
 AliPWGHistoTools.cxx:1318
 AliPWGHistoTools.cxx:1319
 AliPWGHistoTools.cxx:1320
 AliPWGHistoTools.cxx:1321
 AliPWGHistoTools.cxx:1322
 AliPWGHistoTools.cxx:1323
 AliPWGHistoTools.cxx:1324
 AliPWGHistoTools.cxx:1325
 AliPWGHistoTools.cxx:1326
 AliPWGHistoTools.cxx:1327
 AliPWGHistoTools.cxx:1328
 AliPWGHistoTools.cxx:1329
 AliPWGHistoTools.cxx:1330
 AliPWGHistoTools.cxx:1331
 AliPWGHistoTools.cxx:1332
 AliPWGHistoTools.cxx:1333
 AliPWGHistoTools.cxx:1334
 AliPWGHistoTools.cxx:1335
 AliPWGHistoTools.cxx:1336
 AliPWGHistoTools.cxx:1337
 AliPWGHistoTools.cxx:1338
 AliPWGHistoTools.cxx:1339
 AliPWGHistoTools.cxx:1340
 AliPWGHistoTools.cxx:1341
 AliPWGHistoTools.cxx:1342
 AliPWGHistoTools.cxx:1343
 AliPWGHistoTools.cxx:1344
 AliPWGHistoTools.cxx:1345
 AliPWGHistoTools.cxx:1346
 AliPWGHistoTools.cxx:1347
 AliPWGHistoTools.cxx:1348
 AliPWGHistoTools.cxx:1349
 AliPWGHistoTools.cxx:1350
 AliPWGHistoTools.cxx:1351
 AliPWGHistoTools.cxx:1352
 AliPWGHistoTools.cxx:1353
 AliPWGHistoTools.cxx:1354
 AliPWGHistoTools.cxx:1355
 AliPWGHistoTools.cxx:1356
 AliPWGHistoTools.cxx:1357
 AliPWGHistoTools.cxx:1358
 AliPWGHistoTools.cxx:1359
 AliPWGHistoTools.cxx:1360
 AliPWGHistoTools.cxx:1361
 AliPWGHistoTools.cxx:1362
 AliPWGHistoTools.cxx:1363
 AliPWGHistoTools.cxx:1364
 AliPWGHistoTools.cxx:1365
 AliPWGHistoTools.cxx:1366
 AliPWGHistoTools.cxx:1367
 AliPWGHistoTools.cxx:1368
 AliPWGHistoTools.cxx:1369
 AliPWGHistoTools.cxx:1370
 AliPWGHistoTools.cxx:1371
 AliPWGHistoTools.cxx:1372
 AliPWGHistoTools.cxx:1373
 AliPWGHistoTools.cxx:1374
 AliPWGHistoTools.cxx:1375
 AliPWGHistoTools.cxx:1376
 AliPWGHistoTools.cxx:1377
 AliPWGHistoTools.cxx:1378
 AliPWGHistoTools.cxx:1379
 AliPWGHistoTools.cxx:1380
 AliPWGHistoTools.cxx:1381
 AliPWGHistoTools.cxx:1382
 AliPWGHistoTools.cxx:1383
 AliPWGHistoTools.cxx:1384
 AliPWGHistoTools.cxx:1385
 AliPWGHistoTools.cxx:1386
 AliPWGHistoTools.cxx:1387
 AliPWGHistoTools.cxx:1388
 AliPWGHistoTools.cxx:1389
 AliPWGHistoTools.cxx:1390
 AliPWGHistoTools.cxx:1391
 AliPWGHistoTools.cxx:1392
 AliPWGHistoTools.cxx:1393
 AliPWGHistoTools.cxx:1394
 AliPWGHistoTools.cxx:1395
 AliPWGHistoTools.cxx:1396
 AliPWGHistoTools.cxx:1397
 AliPWGHistoTools.cxx:1398
 AliPWGHistoTools.cxx:1399
 AliPWGHistoTools.cxx:1400
 AliPWGHistoTools.cxx:1401
 AliPWGHistoTools.cxx:1402
 AliPWGHistoTools.cxx:1403
 AliPWGHistoTools.cxx:1404
 AliPWGHistoTools.cxx:1405
 AliPWGHistoTools.cxx:1406
 AliPWGHistoTools.cxx:1407
 AliPWGHistoTools.cxx:1408
 AliPWGHistoTools.cxx:1409
 AliPWGHistoTools.cxx:1410
 AliPWGHistoTools.cxx:1411
 AliPWGHistoTools.cxx:1412
 AliPWGHistoTools.cxx:1413
 AliPWGHistoTools.cxx:1414
 AliPWGHistoTools.cxx:1415
 AliPWGHistoTools.cxx:1416
 AliPWGHistoTools.cxx:1417
 AliPWGHistoTools.cxx:1418
 AliPWGHistoTools.cxx:1419
 AliPWGHistoTools.cxx:1420
 AliPWGHistoTools.cxx:1421
 AliPWGHistoTools.cxx:1422
 AliPWGHistoTools.cxx:1423
 AliPWGHistoTools.cxx:1424
 AliPWGHistoTools.cxx:1425
 AliPWGHistoTools.cxx:1426
 AliPWGHistoTools.cxx:1427
 AliPWGHistoTools.cxx:1428
 AliPWGHistoTools.cxx:1429
 AliPWGHistoTools.cxx:1430
 AliPWGHistoTools.cxx:1431
 AliPWGHistoTools.cxx:1432
 AliPWGHistoTools.cxx:1433
 AliPWGHistoTools.cxx:1434
 AliPWGHistoTools.cxx:1435
 AliPWGHistoTools.cxx:1436
 AliPWGHistoTools.cxx:1437
 AliPWGHistoTools.cxx:1438
 AliPWGHistoTools.cxx:1439
 AliPWGHistoTools.cxx:1440
 AliPWGHistoTools.cxx:1441
 AliPWGHistoTools.cxx:1442
 AliPWGHistoTools.cxx:1443
 AliPWGHistoTools.cxx:1444
 AliPWGHistoTools.cxx:1445
 AliPWGHistoTools.cxx:1446
 AliPWGHistoTools.cxx:1447
 AliPWGHistoTools.cxx:1448
 AliPWGHistoTools.cxx:1449
 AliPWGHistoTools.cxx:1450
 AliPWGHistoTools.cxx:1451
 AliPWGHistoTools.cxx:1452
 AliPWGHistoTools.cxx:1453
 AliPWGHistoTools.cxx:1454
 AliPWGHistoTools.cxx:1455
 AliPWGHistoTools.cxx:1456
 AliPWGHistoTools.cxx:1457
 AliPWGHistoTools.cxx:1458
 AliPWGHistoTools.cxx:1459
 AliPWGHistoTools.cxx:1460
 AliPWGHistoTools.cxx:1461
 AliPWGHistoTools.cxx:1462
 AliPWGHistoTools.cxx:1463
 AliPWGHistoTools.cxx:1464
 AliPWGHistoTools.cxx:1465
 AliPWGHistoTools.cxx:1466
 AliPWGHistoTools.cxx:1467
 AliPWGHistoTools.cxx:1468
 AliPWGHistoTools.cxx:1469
 AliPWGHistoTools.cxx:1470
 AliPWGHistoTools.cxx:1471
 AliPWGHistoTools.cxx:1472
 AliPWGHistoTools.cxx:1473
 AliPWGHistoTools.cxx:1474
 AliPWGHistoTools.cxx:1475
 AliPWGHistoTools.cxx:1476
 AliPWGHistoTools.cxx:1477
 AliPWGHistoTools.cxx:1478
 AliPWGHistoTools.cxx:1479
 AliPWGHistoTools.cxx:1480