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

/* $Id: AliUnfolding.cxx 31168 2009-02-23 15:18:45Z jgrosseo $ */

// This class allows 1-dimensional unfolding.
// Methods that are implemented are chi2 minimization and bayesian unfolding.
//
//  Author: Jan.Fiete.Grosse-Oetringhaus@cern.ch

#include "AliUnfolding.h"
#include <TH1F.h>
#include <TH2F.h>
#include <TVirtualFitter.h>
#include <TMath.h>
#include <TCanvas.h>
#include <TF1.h>
#include <TExec.h>
#include "Riostream.h"
#include "TROOT.h"

using namespace std; //required for resolving the 'cout' symbol

TMatrixD* AliUnfolding::fgCorrelationMatrix = 0;
TMatrixD* AliUnfolding::fgCorrelationMatrixSquared = 0;
TMatrixD* AliUnfolding::fgCorrelationCovarianceMatrix = 0;
TVectorD* AliUnfolding::fgCurrentESDVector = 0;
TVectorD* AliUnfolding::fgEntropyAPriori = 0;
TVectorD* AliUnfolding::fgEfficiency = 0;

TAxis* AliUnfolding::fgUnfoldedAxis = 0;
TAxis* AliUnfolding::fgMeasuredAxis = 0;

TF1* AliUnfolding::fgFitFunction = 0;

AliUnfolding::MethodType AliUnfolding::fgMethodType = AliUnfolding::kInvalid;
Int_t AliUnfolding::fgMaxInput  = -1;  // bins in measured histogram
Int_t AliUnfolding::fgMaxParams = -1;  // bins in unfolded histogram = number of fit params
Float_t AliUnfolding::fgOverflowBinLimit = -1;

AliUnfolding::RegularizationType AliUnfolding::fgRegularizationType = AliUnfolding::kPol1;
Float_t AliUnfolding::fgRegularizationWeight = 10000;
Int_t AliUnfolding::fgSkipBinsBegin = 0;
Float_t AliUnfolding::fgMinuitStepSize = 0.1;                 // (usually not needed to be changed) step size in minimization
Float_t AliUnfolding::fgMinuitPrecision = 1e-6;               // minuit precision
Int_t   AliUnfolding::fgMinuitMaxIterations = 1000000;            // minuit maximum number of iterations
Double_t AliUnfolding::fgMinuitStrategy = 1.;                 // minuit strategy
Bool_t AliUnfolding::fgMinimumInitialValue = kFALSE;          // set all initial values at least to the smallest value among the initial values
Float_t AliUnfolding::fgMinimumInitialValueFix = -1;
Bool_t AliUnfolding::fgNormalizeInput = kFALSE;               // normalize input spectrum
Float_t AliUnfolding::fgNotFoundEvents = 0;
Bool_t AliUnfolding::fgSkipBin0InChi2 = kFALSE;

Float_t AliUnfolding::fgBayesianSmoothing  = 1;           // smoothing parameter (0 = no smoothing)
Int_t   AliUnfolding::fgBayesianIterations = 10;          // number of iterations in Bayesian method

Bool_t AliUnfolding::fgDebug = kFALSE;

Int_t AliUnfolding::fgCallCount = 0;

Int_t AliUnfolding::fgPowern = 5;

Double_t AliUnfolding::fChi2FromFit = 0.;
Double_t AliUnfolding::fPenaltyVal  = 0.;
Double_t AliUnfolding::fAvgResidual = 0.;

Int_t AliUnfolding::fgPrintChi2Details = 0;

TCanvas *AliUnfolding::fgCanvas = 0;
TH1 *AliUnfolding::fghUnfolded = 0;     
TH2 *AliUnfolding::fghCorrelation = 0;  
TH1 *AliUnfolding::fghEfficiency = 0;   
TH1 *AliUnfolding::fghMeasured = 0;     

ClassImp(AliUnfolding)

//____________________________________________________________________
void AliUnfolding::SetUnfoldingMethod(MethodType methodType)
{
  // set unfolding method
  fgMethodType = methodType; 
  
  const char* name = 0;
  switch (methodType)
  {
    case kInvalid: name = "INVALID"; break;
    case kChi2Minimization: name = "Chi2 Minimization"; break;
    case kBayesian: name = "Bayesian unfolding"; break;
    case kFunction: name = "Functional fit"; break;
  }
  Printf("AliUnfolding::SetUnfoldingMethod: %s enabled.", name);
}

//____________________________________________________________________
void AliUnfolding::SetCreateOverflowBin(Float_t overflowBinLimit) 
{ 
  // enable the creation of a overflow bin that includes all statistics below the given limit
  
  fgOverflowBinLimit = overflowBinLimit; 
  
  Printf("AliUnfolding::SetCreateOverflowBin: overflow bin limit set to %f", overflowBinLimit);
}

//____________________________________________________________________
void AliUnfolding::SetSkipBinsBegin(Int_t nBins)
{
  // set number of skipped bins in regularization
  
  fgSkipBinsBegin = nBins;
  
  Printf("AliUnfolding::SetSkipBinsBegin: skipping %d bins at the beginning of the spectrum in the regularization.", fgSkipBinsBegin);
}

//____________________________________________________________________
void AliUnfolding::SetNbins(Int_t nMeasured, Int_t nUnfolded) 
{ 
  // set number of bins in the input (measured) distribution and in the unfolded distribution
  fgMaxInput = nMeasured; 
  fgMaxParams = nUnfolded; 
  
  if (fgCorrelationMatrix)
  {
    delete fgCorrelationMatrix;
    fgCorrelationMatrix = 0;
  }
  if (fgCorrelationMatrixSquared)
  {
    fgCorrelationMatrixSquared = 0;
    delete fgCorrelationMatrixSquared;
  }
  if (fgCorrelationCovarianceMatrix)
  {
    delete fgCorrelationCovarianceMatrix;
    fgCorrelationCovarianceMatrix = 0;
  }
  if (fgCurrentESDVector)
  {
    delete fgCurrentESDVector;
    fgCurrentESDVector = 0;
  }
  if (fgEntropyAPriori)
  {
    delete fgEntropyAPriori;
    fgEntropyAPriori = 0;
  }
  if (fgEfficiency)
  {
    delete fgEfficiency;
    fgEfficiency = 0;
  }
  if (fgUnfoldedAxis)
  {
    delete fgUnfoldedAxis;
    fgUnfoldedAxis = 0;
  }
  if (fgMeasuredAxis)
  {
    delete fgMeasuredAxis;
    fgMeasuredAxis = 0;
  }
  
  Printf("AliUnfolding::SetNbins: Set %d measured bins and %d unfolded bins", nMeasured, nUnfolded);
}

//____________________________________________________________________
void AliUnfolding::SetChi2Regularization(RegularizationType type, Float_t weight)
{
  //
  // sets the parameters for chi2 minimization
  //

  fgRegularizationType = type;
  fgRegularizationWeight = weight;

  Printf("AliUnfolding::SetChi2Regularization --> Regularization set to %d with weight %f", (Int_t) type, weight);
}

//____________________________________________________________________
void AliUnfolding::SetBayesianParameters(Float_t smoothing, Int_t nIterations)
{
  //
  // sets the parameters for Bayesian unfolding
  //

  fgBayesianSmoothing = smoothing;
  fgBayesianIterations = nIterations;

  Printf("AliUnfolding::SetBayesianParameters --> Parameters set to %d iterations with smoothing %f", fgBayesianIterations, fgBayesianSmoothing);
}

//____________________________________________________________________
void AliUnfolding::SetFunction(TF1* function)
{
  // set function for unfolding with a fit function
  
  fgFitFunction = function;
  
  Printf("AliUnfolding::SetFunction: Set fit function with %d parameters.", function->GetNpar());
}

//____________________________________________________________________
Int_t AliUnfolding::Unfold(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result, Bool_t check)
{
  // unfolds with unfolding method fgMethodType
  //
  // parameters:
  //  correlation: response matrix as measured vs. generated
  //  efficiency:  (optional) efficiency that is applied on the unfolded spectrum, i.e. it has to be in unfolded variables. If 0 no efficiency is applied.
  //  measured:    the measured spectrum
  //  initialConditions: (optional) initial conditions for the unfolding. if 0 the measured spectrum is used as initial conditions.
  //  result:      target for the unfolded result
  //  check:       depends on the unfolding method, see comments in specific functions
  //
  //  return code: see UnfoldWithMinuit/UnfoldWithBayesian/UnfoldWithFunction

  if (fgMaxInput == -1)
  {
    Printf("AliUnfolding::Unfold: WARNING. Number of measured bins not set with SetNbins. Using number of bins in measured distribution");
    fgMaxInput = measured->GetNbinsX();
  }
  if (fgMaxParams == -1)
  {
    Printf("AliUnfolding::Unfold: WARNING. Number of unfolded bins not set with SetNbins. Using number of bins in measured distribution");
    fgMaxParams = measured->GetNbinsX();
  }

  if (fgOverflowBinLimit > 0)
    CreateOverflowBin(correlation, measured);
    
  switch (fgMethodType)
  {
    case kInvalid:
    {
      Printf("AliUnfolding::Unfold: ERROR: Unfolding method not set. Use SetUnfoldingMethod. Exiting...");
      return -1;
    }
    case kChi2Minimization:
      return UnfoldWithMinuit(correlation, efficiency, measured, initialConditions, result, check);
    case kBayesian:
      return UnfoldWithBayesian(correlation, efficiency, measured, initialConditions, result);
    case kFunction:
      return UnfoldWithFunction(correlation, efficiency, measured, initialConditions, result);
  }



  return -1;
}

//____________________________________________________________________
void AliUnfolding::SetStaticVariables(TH2* correlation, TH1* measured, TH1* efficiency)
{
  // fill static variables needed for minuit fit

  if (!fgCorrelationMatrix)
    fgCorrelationMatrix = new TMatrixD(fgMaxInput, fgMaxParams);
  if (!fgCorrelationMatrixSquared)
    fgCorrelationMatrixSquared = new TMatrixD(fgMaxInput, fgMaxParams);
  if (!fgCorrelationCovarianceMatrix)
    fgCorrelationCovarianceMatrix = new TMatrixD(fgMaxInput, fgMaxInput);
  if (!fgCurrentESDVector)
    fgCurrentESDVector = new TVectorD(fgMaxInput);
  if (!fgEntropyAPriori)
    fgEntropyAPriori = new TVectorD(fgMaxParams);
  if (!fgEfficiency)
    fgEfficiency = new TVectorD(fgMaxParams);
  if (fgUnfoldedAxis)
  {
    delete fgUnfoldedAxis;
    fgUnfoldedAxis = 0;
  }
  if (!fgUnfoldedAxis) 
    fgUnfoldedAxis = new TAxis(*(correlation->GetXaxis()));
  if (fgMeasuredAxis)
  {
    delete fgMeasuredAxis;
    fgMeasuredAxis = 0;
  }
  if (!fgMeasuredAxis) 
    fgMeasuredAxis = new TAxis(*(correlation->GetYaxis()));

  fgCorrelationMatrix->Zero();
  fgCorrelationCovarianceMatrix->Zero();
  fgCurrentESDVector->Zero();
  fgEntropyAPriori->Zero();

  // normalize correction for given nPart
  for (Int_t i=1; i<=correlation->GetNbinsX(); ++i)
  {
    Double_t sum = correlation->Integral(i, i, 1, correlation->GetNbinsY());
    if (sum <= 0)
      continue;
    Float_t maxValue = 0;
    //    Int_t maxBin = -1;
    for (Int_t j=1; j<=correlation->GetNbinsY(); ++j)
    {
      // find most probably value
      if (maxValue < correlation->GetBinContent(i, j))
      {
        maxValue = correlation->GetBinContent(i, j);
	//        maxBin = j;
      }

      // npart sum to 1
      correlation->SetBinContent(i, j, correlation->GetBinContent(i, j) / sum);// * correlation->GetXaxis()->GetBinWidth(i));
      correlation->SetBinError(i, j, correlation->GetBinError(i, j) / sum);

      if (i <= fgMaxParams && j <= fgMaxInput)
      {
        (*fgCorrelationMatrix)(j-1, i-1) = correlation->GetBinContent(i, j);
        (*fgCorrelationMatrixSquared)(j-1, i-1) = correlation->GetBinContent(i, j) * correlation->GetBinContent(i, j);
      }
    }

    //printf("MPV for Ntrue = %f is %f\n", fCurrentCorrelation->GetXaxis()->GetBinCenter(i), fCurrentCorrelation->GetYaxis()->GetBinCenter(maxBin));
  }
    
  //normalize measured
  Float_t smallestError = 1;
  if (fgNormalizeInput)
  {
    Float_t sumMeasured = measured->Integral();
    measured->Scale(1.0 / sumMeasured);
    smallestError /= sumMeasured;
  }
  
  for (Int_t i=0; i<fgMaxInput; ++i)
  {
    (*fgCurrentESDVector)[i] = measured->GetBinContent(i+1);
    if (measured->GetBinError(i+1) > 0)
    {
      (*fgCorrelationCovarianceMatrix)(i, i) = (Double_t) 1e-6 / measured->GetBinError(i+1) / measured->GetBinError(i+1);
    }
    else // in this case put error of 1, otherwise 0 bins are not added to the chi2...
      (*fgCorrelationCovarianceMatrix)(i, i) = (Double_t) 1e-6 / smallestError / smallestError;

    if ((*fgCorrelationCovarianceMatrix)(i, i) > 1e7)
      (*fgCorrelationCovarianceMatrix)(i, i) = 0;
    //Printf("%d, %e", i, (*fgCorrelationCovarianceMatrix)(i, i));
  }

  // efficiency is expected to match bin width of result
  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    (*fgEfficiency)(i) = efficiency->GetBinContent(i+1);
  }

  if (correlation->GetNbinsX() != fgMaxParams || correlation->GetNbinsY() != fgMaxInput)
    cout << "Response histo has incorrect dimensions; expect (" << fgMaxParams << ", " << fgMaxInput << "), got (" << correlation->GetNbinsX() << ", " << correlation->GetNbinsY() << ")" << endl;

}

//____________________________________________________________________
Int_t AliUnfolding::UnfoldWithMinuit(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result, Bool_t check)
{
  //
  // implementation of unfolding (internal function)
  //
  // unfolds <measured> using response from <correlation> and effiency <efficiency>
  // output is in <result>
  // <initialConditions> set the initial values for the minimization, if 0 <measured> is used
  //   negative values in initialConditions mean that the given parameter is fixed to the absolute of the value
  // if <check> is true no unfolding is made, instead only the chi2 without unfolding is printed
  //
  // returns minuit status (0 = success), (-1 when check was set)
  //

  SetStaticVariables(correlation, measured, efficiency);
  
  // Initialize TMinuit via generic fitter interface
  Int_t params = fgMaxParams;
  if (fgNotFoundEvents > 0)
    params++;
    
  TVirtualFitter *minuit = TVirtualFitter::Fitter(0, params);
  Double_t arglist[100];
  //  minuit->SetDefaultFitter("Minuit2");

  // disable any output (-1), unfortuantly we do not see warnings anymore then. Have to find another way...
  arglist[0] = 0;
  minuit->ExecuteCommand("SET PRINT", arglist, 1);

  // however, enable warnings
  //minuit->ExecuteCommand("SET WAR", arglist, 0);

  // set minimization function
  minuit->SetFCN(Chi2Function);

  // set precision
  minuit->SetPrecision(fgMinuitPrecision);

  minuit->SetMaxIterations(fgMinuitMaxIterations);

  minuit->ExecuteCommand("SET STRATEGY",&fgMinuitStrategy,1);

  for (Int_t i=0; i<fgMaxParams; i++)
    (*fgEntropyAPriori)[i] = 1;

  // set initial conditions as a-priori distribution for MRX regularization
  /*
  for (Int_t i=0; i<fgMaxParams; i++)
    if (initialConditions && initialConditions->GetBinContent(i+1) > 0)
      (*fgEntropyAPriori)[i] = initialConditions->GetBinContent(i+1);
  */

  if (!initialConditions) {
    initialConditions = measured;
  } else {
    Printf("AliUnfolding::UnfoldWithMinuit: Using different initial conditions...");
    //new TCanvas; initialConditions->DrawCopy();
    if (fgNormalizeInput)
      initialConditions->Scale(1.0 / initialConditions->Integral());
  }

  // extract minimum value from initial conditions (if we set a value to 0 it will stay 0)
  Float_t minValue = 1e35;
  if (fgMinimumInitialValueFix < 0)
  {
    for (Int_t i=0; i<fgMaxParams; ++i)
    {
      Int_t bin = initialConditions->GetXaxis()->FindBin(result->GetXaxis()->GetBinCenter(i+1));
      if (initialConditions->GetBinContent(bin) > 0)
	minValue = TMath::Min(minValue, (Float_t) initialConditions->GetBinContent(bin));
    }
  }
  else
    minValue = fgMinimumInitialValueFix;
  
  Double_t* results = new Double_t[fgMaxParams+1];
  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    Int_t bin = initialConditions->GetXaxis()->FindBin(result->GetXaxis()->GetBinCenter(i+1));
    results[i] = initialConditions->GetBinContent(bin);

    Bool_t fix = kFALSE;
    if (results[i] < 0)
    {
      fix = kTRUE;
      results[i] = -results[i];
    }
 
    if (!fix && fgMinimumInitialValue && results[i] < minValue)
      results[i] = minValue;
      
    // minuit sees squared values to prevent it from going negative...
    results[i] = TMath::Sqrt(results[i]);

    minuit->SetParameter(i, Form("param%d", i), results[i], (fix) ? 0 : fgMinuitStepSize, 0, 0);
  }
  if (fgNotFoundEvents > 0)
  {
    results[fgMaxParams] = efficiency->GetBinContent(1);
    minuit->SetParameter(fgMaxParams, "vtx0", results[fgMaxParams], fgMinuitStepSize / 100, 0.01, 0.80);
  }
  
  Int_t dummy = 0;
  Double_t chi2 = 0;
  Chi2Function(dummy, 0, chi2, results, 0);
  printf("AliUnfolding::UnfoldWithMinuit: Chi2 of initial parameters is = %f\n", chi2);

  if (check)
  {
    DrawGuess(results);
    delete[] results;
    return -1;
  }

  // first param is number of iterations, second is precision....
  arglist[0] = (float)fgMinuitMaxIterations;
  //  arglist[1] = 1e-5;
  //  minuit->ExecuteCommand("SET PRINT", arglist, 3);
  //  minuit->ExecuteCommand("SCAN", arglist, 0);
  Int_t status = minuit->ExecuteCommand("MIGRAD", arglist, 1);
  Printf("AliUnfolding::UnfoldWithMinuit: MINUIT status is %d", status);
  //printf("!!!!!!!!!!!!!! MIGRAD finished: Starting MINOS !!!!!!!!!!!!!!");
  //minuit->ExecuteCommand("MINOS", arglist, 0);

  if (fgNotFoundEvents > 0)
  {
    results[fgMaxParams] = minuit->GetParameter(fgMaxParams);
    Printf("Efficiency for bin 0 changed from %f to %f", efficiency->GetBinContent(1), results[fgMaxParams]);
    efficiency->SetBinContent(1, results[fgMaxParams]);
  }
  
  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    results[i] = minuit->GetParameter(i);
    Double_t value = results[i] * results[i];
   // error is : 2 * (relError on results[i]) * (value) = 2 * (minuit->GetParError(i) / minuit->GetParameter(i)) * (minuit->GetParameter(i) * minuit->GetParameter(i))
    Double_t error = 0;
    if (TMath::IsNaN(minuit->GetParError(i)))
      Printf("WARNING: Parameter %d error is nan", i);
    else 
      error = 2 * minuit->GetParError(i) * results[i];
    
    if (efficiency)
    {	
      //printf("value before efficiency correction: %f\n",value);
      if (efficiency->GetBinContent(i+1) > 0)
      {
	value /= efficiency->GetBinContent(i+1);
	error /= efficiency->GetBinContent(i+1);
      }
      else
      {
	value = 0;
	error = 0;
      }
    }
    //printf("value after efficiency correction: %f +/- %f\n",value,error);
    result->SetBinContent(i+1, value);
    result->SetBinError(i+1, error);
  }

  Int_t tmpCallCount = fgCallCount;
  fgCallCount = 0; // needs to be 0 so that the Chi2Function prints its output
  Chi2Function(dummy, 0, chi2, results, 0);
  
  Printf("AliUnfolding::UnfoldWithMinuit: iterations %d. Chi2 of final parameters is = %f", tmpCallCount, chi2);
  
  delete[] results;

  return status;
}

//____________________________________________________________________
Int_t AliUnfolding::UnfoldWithBayesian(TH2* correlation, TH1* aEfficiency, TH1* measured, TH1* initialConditions, TH1* aResult)
{
  //
  // unfolds a spectrum using the Bayesian method
  //
  
  if (measured->Integral() <= 0)
  {
    Printf("AliUnfolding::UnfoldWithBayesian: ERROR: The measured spectrum is empty");
    return -1;
  }

  const Int_t kStartBin = 0;

  Int_t kMaxM = fgMaxInput;  //<= fCurrentCorrelation->GetNbinsY(); // max measured axis
  Int_t kMaxT = fgMaxParams; //<= fCurrentCorrelation->GetNbinsX(); // max true axis

  // convergence limit: kMaxT * 0.001^2 = kMaxT * 1e-6 (e.g. 250 bins --> 2.5 e-4)
  const Double_t kConvergenceLimit = kMaxT * 1e-6;

  // store information in arrays, to increase processing speed (~ factor 5)
  Double_t* measuredCopy = new Double_t[kMaxM];
  Double_t* measuredError = new Double_t[kMaxM];
  Double_t* prior = new Double_t[kMaxT];
  Double_t* result = new Double_t[kMaxT];
  Double_t* efficiency = new Double_t[kMaxT];
  Double_t* binWidths = new Double_t[kMaxT];
  Double_t* normResponse = new Double_t[kMaxT]; 

  Double_t** response = new Double_t*[kMaxT];
  Double_t** inverseResponse = new Double_t*[kMaxT];
  for (Int_t i=0; i<kMaxT; i++)
  {
    response[i] = new Double_t[kMaxM];
    inverseResponse[i] = new Double_t[kMaxM];
    normResponse[i] = 0;
  }

  // for normalization
  Float_t measuredIntegral = measured->Integral();
  for (Int_t m=0; m<kMaxM; m++)
  {
    measuredCopy[m] = measured->GetBinContent(m+1) / measuredIntegral;
    measuredError[m] = measured->GetBinError(m+1) / measuredIntegral;

    for (Int_t t=0; t<kMaxT; t++)
    {
      response[t][m] = correlation->GetBinContent(t+1, m+1);
      inverseResponse[t][m] = 0;
      normResponse[t] += correlation->GetBinContent(t+1, m+1);
    }
  }

  for (Int_t t=0; t<kMaxT; t++)
  {
    if (aEfficiency)
    {
      efficiency[t] = aEfficiency->GetBinContent(t+1);
    }
    else
      efficiency[t] = 1;
      
    prior[t] = measuredCopy[t];
    result[t] = 0;
    binWidths[t] = aResult->GetXaxis()->GetBinWidth(t+1);

    for (Int_t m=0; m<kMaxM; m++) { // Normalise response matrix
      if (normResponse[t] != 0) 
	response[t][m] /= normResponse[t];
      else
        Printf("AliUnfolding::UnfoldWithBayesian: Empty row,column in response matrix, for truth bin %d",t);
    }
  }

  // pick prior distribution
  if (initialConditions)
  {
    printf("Using different starting conditions...\n");
    // for normalization
    Float_t inputDistIntegral = initialConditions->Integral();
    for (Int_t i=0; i<kMaxT; i++)
      prior[i] = initialConditions->GetBinContent(i+1) / inputDistIntegral;
  }

  //TH1F* convergence = new TH1F("convergence", "convergence", 200, 0.5, 200.5);
  
  //new TCanvas;
  // unfold...
  for (Int_t i=0; i<fgBayesianIterations || fgBayesianIterations < 0; i++)
  {
    if (fgDebug)
      Printf("AliUnfolding::UnfoldWithBayesian: iteration %i", i);

    // calculate IR from Bayes theorem
    // IR_ji = R_ij * prior_i / sum_k(R_kj * prior_k)

    Double_t chi2Measured = 0;
    for (Int_t m=0; m<kMaxM; m++)
    {
      Float_t norm = 0;
      for (Int_t t = kStartBin; t<kMaxT; t++)
        norm += response[t][m] * prior[t] * efficiency[t];

      // calc. chi2: (measured - response * prior) / error
      if (measuredError[m] > 0)
      {
        Double_t value = (measuredCopy[m] - norm) / measuredError[m];
        chi2Measured += value * value;
      }

      if (norm > 0)
      {
        for (Int_t t = kStartBin; t<kMaxT; t++)
          inverseResponse[t][m] = response[t][m] * prior[t] * efficiency[t] / norm;
      }
      else
      {
        for (Int_t t = kStartBin; t<kMaxT; t++)
          inverseResponse[t][m] = 0;
      }
    }
    //Printf("chi2Measured of the last prior is %e", chi2Measured);

    for (Int_t t = kStartBin; t<kMaxT; t++)
    {
      Float_t value = 0;
      for (Int_t m=0; m<kMaxM; m++)
        value += inverseResponse[t][m] * measuredCopy[m];

      if (efficiency[t] > 0)
        result[t] = value / efficiency[t];
      else
        result[t] = 0;
    }
    
    /* 
    // draw intermediate result
    for (Int_t t=0; t<kMaxT; t++)
    {
      aResult->SetBinContent(t+1, result[t]);
    }
    aResult->SetMarkerStyle(24+i);
    aResult->SetMarkerColor(2);
    aResult->DrawCopy((i == 0) ? "P" : "PSAME");
    */
 
    Double_t chi2LastIter = 0;
    // regularization (simple smoothing)
    for (Int_t t=kStartBin; t<kMaxT; t++)
    {
      Float_t newValue = 0;
      
      // 0 bin excluded from smoothing
      if (t > kStartBin+2 && t<kMaxT-1)
      {
        Float_t average = (result[t-1] / binWidths[t-1] + result[t] / binWidths[t] + result[t+1] / binWidths[t+1]) / 3 * binWidths[t];

        // weight the average with the regularization parameter
        newValue = (1 - fgBayesianSmoothing) * result[t] + fgBayesianSmoothing * average;
      }
      else
        newValue = result[t];

      // calculate chi2 (change from last iteration)
      if (prior[t] > 1e-5)
      {
        Double_t diff = (prior[t] - newValue) / prior[t];
        chi2LastIter += diff * diff;
      }

      prior[t] = newValue;
    }
    //printf("Chi2 of %d iteration = %e\n", i, chi2LastIter);
    //convergence->Fill(i+1, chi2LastIter);

    if (fgBayesianIterations < 0 && chi2LastIter < kConvergenceLimit)
    {
      Printf("AliUnfolding::UnfoldWithBayesian: Stopped Bayesian unfolding after %d iterations at chi2(change since last iteration) of %e; chi2Measured of the last prior is %e", i, chi2LastIter, chi2Measured);
      break;
    }
  } // end of iterations

  //new TCanvas; convergence->DrawCopy(); gPad->SetLogy();
  //delete convergence;

  Float_t factor = 1;
  if (!fgNormalizeInput)
    factor = measuredIntegral;
  for (Int_t t=0; t<kMaxT; t++)
    aResult->SetBinContent(t+1, result[t] * factor);

  delete[] measuredCopy;
  delete[] measuredError;
  delete[] prior;
  delete[] result;
  delete[] efficiency;
  delete[] binWidths;
  delete[] normResponse;

  for (Int_t i=0; i<kMaxT; i++)
  {
    delete[] response[i];
    delete[] inverseResponse[i];
  }
  delete[] response;
  delete[] inverseResponse;
  
  return 0;

  // ********
  // Calculate the covariance matrix, all arguments are taken from NIM,A362,487-498,1995

  /*printf("Calculating covariance matrix. This may take some time...\n");

  // check if this is the right one...
  TH1* sumHist = GetMultiplicityMC(inputRange, eventType)->ProjectionY("sumHist", 1, GetMultiplicityMC(inputRange, eventType)->GetNbinsX());

  Int_t xBins = hInverseResponseBayes->GetNbinsX();
  Int_t yBins = hInverseResponseBayes->GetNbinsY();

  // calculate "unfolding matrix" Mij
  Float_t matrixM[251][251];
  for (Int_t i=1; i<=xBins; i++)
  {
    for (Int_t j=1; j<=yBins; j++)
    {
      if (fCurrentEfficiency->GetBinContent(i) > 0)
        matrixM[i-1][j-1] = hInverseResponseBayes->GetBinContent(i, j) / fCurrentEfficiency->GetBinContent(i);
      else
        matrixM[i-1][j-1] = 0;
    }
  }

  Float_t* vectorn = new Float_t[yBins];
  for (Int_t j=1; j<=yBins; j++)
    vectorn[j-1] = fCurrentESD->GetBinContent(j);

  // first part of covariance matrix, depends on input distribution n(E)
  Float_t cov1[251][251];

  Float_t nEvents = fCurrentESD->Integral(); // N

  xBins = 20;
  yBins = 20;

  for (Int_t k=0; k<xBins; k++)
  {
    printf("In Cov1: %d\n", k);
    for (Int_t l=0; l<yBins; l++)
    {
      cov1[k][l] = 0;

      // sum_j Mkj Mlj n(Ej) * (1 - n(Ej) / N)
      for (Int_t j=0; j<yBins; j++)
        cov1[k][l] += matrixM[k][j] * matrixM[l][j] * vectorn[j]
          * (1.0 - vectorn[j] / nEvents);

      // - sum_i,j (i != j) Mki Mlj n(Ei) n(Ej) / N
      for (Int_t i=0; i<yBins; i++)
        for (Int_t j=0; j<yBins; j++)
        {
          if (i == j)
            continue;
          cov1[k][l] -= matrixM[k][i] * matrixM[l][j] * vectorn[i]
            * vectorn[j] / nEvents;
         }
    }
  }

  printf("Cov1 finished\n");

  TH2F* cov = (TH2F*) hInverseResponseBayes->Clone("cov");
  cov->Reset();

  for (Int_t i=1; i<=xBins; i++)
    for (Int_t j=1; j<=yBins; j++)
      cov->SetBinContent(i, j, cov1[i-1][j-1]);

  new TCanvas;
  cov->Draw("COLZ");

  // second part of covariance matrix, depends on response matrix
  Float_t cov2[251][251];

  // Cov[P(Er|Cu), P(Es|Cu)] term
  Float_t covTerm[100][100][100];
  for (Int_t r=0; r<yBins; r++)
    for (Int_t u=0; u<xBins; u++)
      for (Int_t s=0; s<yBins; s++)
      {
        if (r == s)
          covTerm[r][u][s] = 1.0 / sumHist->GetBinContent(u+1) * hResponse->GetBinContent(u+1, r+1)
            * (1.0 - hResponse->GetBinContent(u+1, r+1));
        else
          covTerm[r][u][s] = - 1.0 / sumHist->GetBinContent(u+1) * hResponse->GetBinContent(u+1, r+1)
            * hResponse->GetBinContent(u+1, s+1);
      }

  for (Int_t k=0; k<xBins; k++)
    for (Int_t l=0; l<yBins; l++)
    {
      cov2[k][l] = 0;
      printf("In Cov2: %d %d\n", k, l);
      for (Int_t i=0; i<yBins; i++)
        for (Int_t j=0; j<yBins; j++)
        {
          //printf("In Cov2: %d %d %d %d\n", k, l, i, j);
          // calculate Cov(Mki, Mlj) = sum{ru},{su} ...
          Float_t tmpCov = 0;
          for (Int_t r=0; r<yBins; r++)
            for (Int_t u=0; u<xBins; u++)
              for (Int_t s=0; s<yBins; s++)
              {
                if (hResponse->GetBinContent(u+1, r+1) == 0 || hResponse->GetBinContent(u+1, s+1) == 0
                  || hResponse->GetBinContent(u+1, i+1) == 0)
                  continue;

                tmpCov += BayesCovarianceDerivate(matrixM, hResponse, fCurrentEfficiency, k, i, r, u)
                  * BayesCovarianceDerivate(matrixM, hResponse, fCurrentEfficiency, l, j, s, u)
                  * covTerm[r][u][s];
              }

          cov2[k][l] += fCurrentESD->GetBinContent(i+1) * fCurrentESD->GetBinContent(j+1) * tmpCov;
        }
    }

  printf("Cov2 finished\n");

  for (Int_t i=1; i<=xBins; i++)
    for (Int_t j=1; j<=yBins; j++)
      cov->SetBinContent(i, j, cov1[i-1][j-1] + cov2[i-1][j-1]);

  new TCanvas;
  cov->Draw("COLZ");*/
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationPol0(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // prefers constant function (pol0)
  //
  // Does not take into account efficiency
  Double_t chi2 = 0;

  for (Int_t i=1+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    Double_t right  = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
    Double_t left   = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);

    if (left != 0)
    {
      Double_t diff = (right - left);
      chi2 += diff * diff / left / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
    }
  }

  return chi2 / 100.0;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationPol1(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // prefers linear function (pol1)
  //
  // Does not take into account efficiency
  Double_t chi2 = 0;

  for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    if (params[i-1] == 0)
      continue;

    Double_t right  = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
    Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
    Double_t left   = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);

    Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
    Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);

    //Double_t diff = (der1 - der2) / middle;
    //chi2 += diff * diff;
    chi2 += (der1 - der2) * (der1 - der2) / middle * fgUnfoldedAxis->GetBinWidth(i);
  }

  return chi2;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationLog(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // prefers logarithmic function (log)
  //
  // Does not take into account efficiency

  Double_t chi2 = 0;

  for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    if (params[i-1] == 0 || params[i] == 0 || params[i-2] == 0)
     continue;

    Double_t right  = log(params[i] / fgUnfoldedAxis->GetBinWidth(i+1));
    Double_t middle = log(params[i-1] / fgUnfoldedAxis->GetBinWidth(i));
    Double_t left   = log(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1));
    
    Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
    Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
    
    //Double_t error = 1. / params[i] + 4. / params[i-1] + 1. / params[i-2];

    //if (fgCallCount == 0)
    //  Printf("%d %f %f", i, (der1 - der2) * (der1 - der2), error);
    chi2 += (der1 - der2) * (der1 - der2);// / error;
  }

  return chi2;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationTotalCurvature(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // minimizes the total curvature (from Unfolding Methods In High-Energy Physics Experiments,
  // V. Blobel (Hamburg U.) . DESY 84/118, Dec 1984. 40pp.
  //
  // Does not take into account efficiency

  Double_t chi2 = 0;

  for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    Double_t right  = params[i];
    Double_t middle = params[i-1];
    Double_t left   = params[i-2];

    Double_t der1 = (right - middle);
    Double_t der2 = (middle - left);

    Double_t diff = (der1 - der2);

    chi2 += diff * diff;
  }

  return chi2 * 1e4;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationEntropy(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // calculates entropy, from
  // The method of reduced cross-entropy (M. Schmelling 1993)
  //
  // Does not take into account efficiency

  Double_t paramSum = 0;
  
  for (Int_t i=fgSkipBinsBegin; i<fgMaxParams; ++i)
    paramSum += params[i];

  Double_t chi2 = 0;
  for (Int_t i=fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    Double_t tmp = params[i] / paramSum;
    //Double_t tmp = params[i];
    if (tmp > 0 && (*fgEntropyAPriori)[i] > 0)
    {
      chi2 += tmp * TMath::Log(tmp / (*fgEntropyAPriori)[i]);
    }
    else
      chi2 += 100;
  }

  return -chi2;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationRatio(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  //
  // Does not take into account efficiency

  Double_t chi2 = 0;

  for (Int_t i=5+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    if (params[i-1] == 0 || params[i] == 0)
      continue;

    Double_t right  = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
    Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
    Double_t left   = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);
    Double_t left2   = params[i-3] / fgUnfoldedAxis->GetBinWidth(i-2);
    Double_t left3   = params[i-4] / fgUnfoldedAxis->GetBinWidth(i-3);
    Double_t left4   = params[i-5] / fgUnfoldedAxis->GetBinWidth(i-4);

    //Double_t diff = left / middle - middle / right;
    //Double_t diff = 2 * left / middle - middle / right - left2 / left;
    Double_t diff = 4 * left2 / left - middle / right - left / middle - left3 / left2 - left4 / left3;
    
    chi2 += diff * diff;// / middle;
  }

  return chi2;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationPowerLaw(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // prefers power law with n = -5
  //
  // Does not take into account efficiency

  Double_t chi2 = 0;

  Double_t right = 0.;
  Double_t middle = 0.;
  Double_t left = 0.;

  for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
      continue;

    if (fgUnfoldedAxis->GetBinWidth(i+1) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i-1) < 1e-8)
      continue;
    
    middle = TMath::Power(params[i-1] / fgUnfoldedAxis->GetBinWidth(i),fgPowern);

    if(middle>0) {
      right  = TMath::Power(params[i] / fgUnfoldedAxis->GetBinWidth(i+1),fgPowern)/middle;

      left   = TMath::Power(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1),fgPowern)/middle;
      
      middle = 1.;
      
      Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
      Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
    
      chi2 += (der1 - der2) * (der1 - der2)/ ( fgUnfoldedAxis->GetBinWidth(i+1)/2. + fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)/2.)/( fgUnfoldedAxis->GetBinWidth(i+1)/2. + fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)/2.);// / error;
      //   printf("i: %d chi2 = %f\n",i,chi2);
    }

  }

  return chi2;
}

//____________________________________________________________________
Double_t AliUnfolding::RegularizationLogLog(TVectorD& params)
{
  // homogenity term for minuit fitting
  // pure function of the parameters
  // prefers a powerlaw (linear on a log-log scale)
  //
  // The calculation takes into account the efficiencies

  Double_t chi2 = 0;

  for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    if (params[i-1] == 0 || params[i] == 0 || params[i-2] == 0)
     continue;
    if ((*fgEfficiency)(i-1) == 0 || (*fgEfficiency)(i) == 0 || (*fgEfficiency)(i-2) == 0)
     continue;

    Double_t right  = log(params[i] / (*fgEfficiency)(i) / fgUnfoldedAxis->GetBinWidth(i+1));
    Double_t middle = log(params[i-1] / (*fgEfficiency)(i-1) / fgUnfoldedAxis->GetBinWidth(i));
    Double_t left   = log(params[i-2] / (*fgEfficiency)(i-2) / fgUnfoldedAxis->GetBinWidth(i-1));
    
    Double_t der1 = (right - middle) / ( log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i)) );
    Double_t der2 = (middle - left) /( log(fgUnfoldedAxis->GetBinCenter(i)) - log(fgUnfoldedAxis->GetBinCenter(i-1)) );
    
    double tmp = (log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i-1)))/2.;
    Double_t dder = (der1-der2) / tmp;

    chi2 += dder * dder;
  }

  return chi2;
}



//____________________________________________________________________
void AliUnfolding::Chi2Function(Int_t&, Double_t*, Double_t& chi2, Double_t *params, Int_t)
{
  //
  // fit function for minuit
  // does: (m - Ad)W(m - Ad) where m = measured, A correlation matrix, d = guess, W = covariance matrix
  //
  
  // TODO use static members for the variables here to speed up processing (no construction/deconstruction)

  // d = guess
  TVectorD paramsVector(fgMaxParams);
  for (Int_t i=0; i<fgMaxParams; ++i)
    paramsVector[i] = params[i] * params[i];

  // calculate penalty factor
  Double_t penaltyVal = 0;

  switch (fgRegularizationType)
  {
    case kNone:       break;
    case kPol0:       penaltyVal = RegularizationPol0(paramsVector); break;
    case kPol1:       penaltyVal = RegularizationPol1(paramsVector); break;
    case kCurvature:  penaltyVal = RegularizationTotalCurvature(paramsVector); break;
    case kEntropy:    penaltyVal = RegularizationEntropy(paramsVector); break;
    case kLog:        penaltyVal = RegularizationLog(paramsVector); break;
    case kRatio:      penaltyVal = RegularizationRatio(paramsVector); break;
    case kPowerLaw:   penaltyVal = RegularizationPowerLaw(paramsVector); break;
    case kLogLog:     penaltyVal = RegularizationLogLog(paramsVector); break;
  }

  penaltyVal *= fgRegularizationWeight; //beta*PU

  // Ad
  TVectorD measGuessVector(fgMaxInput);
  measGuessVector = (*fgCorrelationMatrix) * paramsVector;

  // Ad - m
  measGuessVector -= (*fgCurrentESDVector);
  
#if 0
  // new error calcuation using error on the guess
  
  // error from guess
  TVectorD errorGuessVector(fgMaxInput);
  //errorGuessVector = (*fgCorrelationMatrixSquared) * paramsVector;
  errorGuessVector = (*fgCorrelationMatrix) * paramsVector;

  Double_t chi2FromFit = 0;
  for (Int_t i=0; i<fgMaxInput; ++i)
  {
    Float_t error = 1;
    if (errorGuessVector(i) > 0)
      error = errorGuessVector(i);
    chi2FromFit += measGuessVector(i) * measGuessVector(i) / error;
  }

#else
  // old error calcuation using the error on the measured
  TVectorD copy(measGuessVector);

  // (Ad - m) W
  // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
  // normal way is like this:
  // measGuessVector *= (*fgCorrelationCovarianceMatrix);
  // optimized way like this:
  for (Int_t i=0; i<fgMaxInput; ++i)
    measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);


  if (fgSkipBin0InChi2)
    measGuessVector[0] = 0;

  // (Ad - m) W (Ad - m)
  // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
  // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see UnfoldWithMinuit)
  Double_t chi2FromFit = measGuessVector * copy * 1e6;
#endif

  Double_t notFoundEventsConstraint = 0;
  Double_t currentNotFoundEvents = 0;
  Double_t errorNotFoundEvents = 0;
  
  if (fgNotFoundEvents > 0)
  {
    for (Int_t i=0; i<fgMaxParams; ++i)
    {
      Double_t eff = (1.0 / (*fgEfficiency)(i) - 1);
      if (i == 0)
	eff = (1.0 / params[fgMaxParams] - 1);
      currentNotFoundEvents += eff * paramsVector(i);
      errorNotFoundEvents += eff * eff * paramsVector(i); // error due to guess (paramsVector)
      if (i <= 3)
	errorNotFoundEvents += (eff * 0.03) * (eff * 0.03) * paramsVector(i) * paramsVector(i); // error on eff
      //      if ((fgCallCount % 10000) == 0)
	//Printf("%d %f %f %f", i, (*fgEfficiency)(i), paramsVector(i), currentNotFoundEvents);
    }
    errorNotFoundEvents += fgNotFoundEvents;
    // TODO add error on background, background estimate
    
    notFoundEventsConstraint = (currentNotFoundEvents - fgNotFoundEvents) * (currentNotFoundEvents - fgNotFoundEvents) / errorNotFoundEvents;
  }
  
  Float_t avg = 0;
  Float_t sum = 0;
  Float_t currentMult = 0;
  // try to match dn/deta
  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    avg += paramsVector(i) * currentMult;
    sum += paramsVector(i);
    currentMult += fgUnfoldedAxis->GetBinWidth(i);
  }
  avg /= sum;
  Float_t chi2Avg = 0; //(avg - 3.73) * (avg - 3.73) * 100;

  chi2 = chi2FromFit + penaltyVal + notFoundEventsConstraint + chi2Avg;

  if ((fgCallCount++ % 1000) == 0)
  {

    Printf("AliUnfolding::Chi2Function: Iteration %d (ev %d %d +- %f) (%f) (%f): %f %f %f %f --> %f", fgCallCount-1, (Int_t) fgNotFoundEvents, (Int_t) currentNotFoundEvents, TMath::Sqrt(errorNotFoundEvents), params[fgMaxParams-1], avg, chi2FromFit, penaltyVal, notFoundEventsConstraint, chi2Avg, chi2);

    //for (Int_t i=0; i<fgMaxInput; ++i)
    //  Printf("%d: %f", i, measGuessVector(i) * copy(i) * 1e6);
  }

  fChi2FromFit = chi2FromFit;
  fPenaltyVal = penaltyVal;
}

//____________________________________________________________________
void AliUnfolding::MakePads() {
    TPad *presult = new TPad("presult","result",0,0.4,1,1);
    presult->SetNumber(1);
    presult->Draw();
    presult->SetLogy();
    TPad *pres = new TPad("pres","residuals",0,0.2,1,0.4);
    pres->SetNumber(2);
    pres->Draw();
    TPad *ppen = new TPad("ppen","penalty",0,0.0,1,0.2);
    ppen->SetNumber(3);
    ppen->Draw();
 
}
//____________________________________________________________________
void AliUnfolding::DrawResults(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TCanvas *canv, Int_t reuseHists,TH1 *unfolded)
{
  // Draw histograms of
  //   - Result folded with response
  //   - Penalty factors
  //   - Chisquare contributions
  // (Useful for debugging/sanity checks and the interactive unfolder)
  //
  // If a canvas pointer is given (canv != 0), it will be used for all
  // plots; 3 pads are made if needed.


  Int_t blankCanvas = 0;
  TVirtualPad *presult = 0;
  TVirtualPad *pres = 0; 
  TVirtualPad *ppen = 0;
  
  if (canv) {
    canv->cd();
    presult = canv->GetPad(1);
    pres = canv->GetPad(2);
    ppen = canv->GetPad(3); 
    if (presult == 0 || pres == 0 || ppen == 0) {
      canv->Clear();
      MakePads();
      blankCanvas = 1;
      presult = canv->GetPad(1);
      pres = canv->GetPad(2);
      ppen = canv->GetPad(3); 
    } 
  }
  else {
    presult = new TCanvas;
    pres = new TCanvas;
    ppen = new TCanvas;
  }


  if (fgMaxInput == -1)
  {
    Printf("AliUnfolding::Unfold: WARNING. Number of measured bins not set with SetNbins. Using number of bins in measured distribution");
    fgMaxInput = measured->GetNbinsX();
  }
  if (fgMaxParams == -1)
  {
    Printf("AliUnfolding::Unfold: WARNING. Number of unfolded bins not set with SetNbins. Using number of bins in measured distribution");
    fgMaxParams = initialConditions->GetNbinsX();
  }

  if (fgOverflowBinLimit > 0)
    CreateOverflowBin(correlation, measured);

  // Uses Minuit methods

  SetStaticVariables(correlation, measured, efficiency);

  Double_t* params = new Double_t[fgMaxParams+1];
  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    params[i] = initialConditions->GetBinContent(i+1) * efficiency->GetBinContent(i+1);

    if (params[i] < 0)
      params[i] = -params[i];
    params[i]=TMath::Sqrt(params[i]);
  } 

  Double_t chi2;
  Int_t dummy;

  //fgPrintChi2Details = kTRUE;
  fgCallCount = 0; // To make sure that Chi2Function prints the components
  Chi2Function(dummy, 0, chi2, params, 0);

  presult->cd();
  TH1 *meas2 = (TH1*)measured->Clone("meas2");
  meas2->SetTitle("meas2");
  meas2->SetName("meas2");
  meas2->SetLineColor(1);
  meas2->SetLineWidth(2);
  meas2->SetMarkerColor(meas2->GetLineColor());
  meas2->SetMarkerStyle(20);
  Float_t scale = unfolded->GetXaxis()->GetBinWidth(1)/meas2->GetXaxis()->GetBinWidth(1);
  meas2->Scale(scale);
  if (blankCanvas)
    meas2->DrawCopy();
  else
    meas2->DrawCopy("same");
  //meas2->GetXaxis()->SetLimits(0,fgMaxInput);
  meas2->SetBit(kCannotPick);
  DrawGuess(params, presult, pres, ppen, reuseHists,unfolded);
  delete [] params;
}
//____________________________________________________________________
void AliUnfolding::RedrawInteractive() {
  // 
  // Helper function for interactive unfolding
  //
  DrawResults(fghCorrelation,fghEfficiency,fghMeasured,fghUnfolded,fgCanvas,1,fghUnfolded);
}
//____________________________________________________________________
void AliUnfolding::InteractiveUnfold(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions) 
{
  //
  // Function to do interactive unfolding
  // A canvas is drawn with the unfolding result
  // Change the histogram with your mouse and all histograms 
  // will be updated automatically

  fgCanvas = new TCanvas("UnfoldingCanvas","Interactive unfolding",500,800);
  fgCanvas->cd();

  MakePads();

  if (fghUnfolded) {
    delete fghUnfolded; fghUnfolded = 0;
  }
  
  gROOT->SetEditHistograms(kTRUE);

  fghUnfolded = new TH1F("AliUnfoldingfghUnfolded","Unfolded result (Interactive unfolder",efficiency->GetNbinsX(),efficiency->GetXaxis()->GetXmin(),efficiency->GetXaxis()->GetXmax());

  fghCorrelation = correlation;
  fghEfficiency = efficiency;
  fghMeasured = measured;

  if(fgMinuitMaxIterations>0)
    UnfoldWithMinuit(correlation,efficiency,measured,initialConditions,fghUnfolded,kFALSE);
  else
    fghUnfolded = initialConditions;

  fghUnfolded->SetLineColor(2);
  fghUnfolded->SetMarkerColor(2);
  fghUnfolded->SetLineWidth(2);


  fgCanvas->cd(1);
  fghUnfolded->Draw();
  DrawResults(correlation,efficiency,measured,fghUnfolded,fgCanvas,kFALSE,fghUnfolded);

  TExec *execRedraw = new TExec("redraw","AliUnfolding::RedrawInteractive()");
  fghUnfolded->GetListOfFunctions()->Add(execRedraw);
}
//____________________________________________________________________
void AliUnfolding::DrawGuess(Double_t *params, TVirtualPad *pfolded, TVirtualPad *pres, TVirtualPad *ppen, Int_t reuseHists,TH1* unfolded)
{
  //
  // draws residuals of solution suggested by params and effect of regularization
  //

  if (pfolded == 0)
    pfolded = new TCanvas;
  if (ppen == 0)
    ppen = new TCanvas;
  if (pres == 0)
    pres = new TCanvas;
  
  // d
  TVectorD paramsVector(fgMaxParams);
  for (Int_t i=0; i<fgMaxParams; ++i)
    paramsVector[i] = params[i] * params[i];

  // Ad
  TVectorD measGuessVector(fgMaxInput);
  measGuessVector = (*fgCorrelationMatrix) * paramsVector;

  TH1* folded = 0;
  if (reuseHists) 
    folded = dynamic_cast<TH1*>(gROOT->FindObject("__hfolded"));
  if (!reuseHists || folded == 0) {
    if (fgMeasuredAxis->GetXbins()->GetArray()) // variable bins
      folded = new TH1F("__hfolded","Folded histo from AliUnfolding",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXbins()->GetArray());
    else
      folded = new TH1F("__hfolded","Folded histo from AliUnfolding",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXmin(),fgMeasuredAxis->GetXmax());
  }

  folded->SetBit(kCannotPick);
  folded->SetLineColor(4);
  folded->SetLineWidth(2);

  for (Int_t ibin =0; ibin < fgMaxInput; ibin++) 
    folded->SetBinContent(ibin+1, measGuessVector[ibin]);

  Float_t scale = unfolded->GetXaxis()->GetBinWidth(1)/folded->GetXaxis()->GetBinWidth(1);
  folded->Scale(scale);

  pfolded->cd();

  folded->Draw("same");

  // Ad - m
  measGuessVector -= (*fgCurrentESDVector);

  TVectorD copy(measGuessVector);

  // (Ad - m) W
  // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
  // normal way is like this:
  // measGuessVector *= (*fgCorrelationCovarianceMatrix);
  // optimized way like this:
  for (Int_t i=0; i<fgMaxInput; ++i)
    measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);

  // (Ad - m) W (Ad - m)
  // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
  // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see ApplyMinuitFit)
  //Double_t chi2FromFit = measGuessVector * copy * 1e6;

  // draw residuals
  // Double_t pTarray[fgMaxParams+1];
  // for(int i=0; i<fgMaxParams; i++) {
  //   pTarray[i] = fgUnfoldedAxis->GetBinCenter(i)-0.5*fgUnfoldedAxis->GetBinWidth(i);
  // }
  //  pTarray[fgMaxParams] = fgUnfoldedAxis->GetBinCenter(fgMaxParams-1)+0.5*fgUnfoldedAxis->GetBinWidth(fgMaxParams-1);
  //  TH1* residuals = new TH1F("residuals", "residuals", fgMaxParams,pTarray);
  //  TH1* residuals = new TH1F("residuals", "residuals", fgMaxInput, -0.5, fgMaxInput - 0.5);
  // for (Int_t i=0; i<fgMaxInput; ++i)
  //   residuals->SetBinContent(i+1, measGuessVector(i) * copy(i) * 1e6);7
  TH1* residuals = GetResidualsPlot(params);
  residuals->GetXaxis()->SetTitleSize(0.06);
  residuals->GetXaxis()->SetTitleOffset(0.7);
  residuals->GetXaxis()->SetLabelSize(0.07);
  residuals->GetYaxis()->SetTitleSize(0.08);
  residuals->GetYaxis()->SetTitleOffset(0.5);
  residuals->GetYaxis()->SetLabelSize(0.06);
  pres->cd(); residuals->DrawCopy(); gPad->SetLogy();
    

  // draw penalty
  TH1* penalty = GetPenaltyPlot(params);
  penalty->GetXaxis()->SetTitleSize(0.06);
  penalty->GetXaxis()->SetTitleOffset(0.7);
  penalty->GetXaxis()->SetLabelSize(0.07);
  penalty->GetYaxis()->SetTitleSize(0.08);
  penalty->GetYaxis()->SetTitleOffset(0.5);
  penalty->GetYaxis()->SetLabelSize(0.06);
  ppen->cd(); penalty->DrawCopy(); gPad->SetLogy();
}

//____________________________________________________________________
TH1* AliUnfolding::GetResidualsPlot(TH1* corrected)
{
  //
  // MvL: THIS MUST BE INCORRECT. 
  // Need heff to calculate params from TH1 'corrected'
  //

  //
  // fill residuals histogram of solution suggested by params and effect of regularization
  //

  Double_t* params = new Double_t[fgMaxParams];
  for (Int_t i=0; i<fgMaxParams; i++)
    params[i] = 0;

  for (Int_t i=0; i<TMath::Min(fgMaxParams, corrected->GetNbinsX()); i++)
    params[i] = TMath::Sqrt(TMath::Abs(corrected->GetBinContent(i+1)*(*fgEfficiency)(i)));

  TH1 * plot = GetResidualsPlot(params);
  delete [] params;
  return plot;
}

//____________________________________________________________________
TH1* AliUnfolding::GetResidualsPlot(Double_t* params)
{
  //
  // fill residuals histogram of solution suggested by params and effect of regularization
  //

  // d
  TVectorD paramsVector(fgMaxParams);
  for (Int_t i=0; i<fgMaxParams; ++i)
    paramsVector[i] = params[i] * params[i];

  // Ad
  TVectorD measGuessVector(fgMaxInput);
  measGuessVector = (*fgCorrelationMatrix) * paramsVector;

  // Ad - m
  measGuessVector -= (*fgCurrentESDVector);

  TVectorD copy(measGuessVector);

  // (Ad - m) W
  // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
  // normal way is like this:
  // measGuessVector *= (*fgCorrelationCovarianceMatrix);
  // optimized way like this:
  for (Int_t i=0; i<fgMaxInput; ++i)
    measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);

  // (Ad - m) W (Ad - m)
  // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
  // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see ApplyMinuitFit)
  //Double_t chi2FromFit = measGuessVector * copy * 1e6;

  // draw residuals
  TH1* residuals = 0;
  if (fgMeasuredAxis->GetXbins()->GetArray()) // variable bins
    residuals = new TH1F("residuals", "residuals;unfolded pos;residual",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXbins()->GetArray());
  else
    residuals = new TH1F("residuals", "residuals;unfolded pos;residual",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXmin(), fgMeasuredAxis->GetXmax());
  //  TH1* residuals = new TH1F("residuals", "residuals", fgMaxInput, -0.5, fgMaxInput - 0.5);

  Double_t sumResiduals = 0.; 
  for (Int_t i=0; i<fgMaxInput; ++i) {
    residuals->SetBinContent(i+1, measGuessVector(i) * copy(i) * 1e6);
    sumResiduals += measGuessVector(i) * copy(i) * 1e6;
  }
  fAvgResidual = sumResiduals/(double)fgMaxInput;
 
  return residuals;
}

//____________________________________________________________________
TH1* AliUnfolding::GetPenaltyPlot(TH1* corrected)
{
  // draws the penalty factors as function of multiplicity of the current selected regularization

  Double_t* params = new Double_t[fgMaxParams];
  for (Int_t i=0; i<fgMaxParams; i++)
    params[i] = 0;

  for (Int_t i=0; i<TMath::Min(fgMaxParams, corrected->GetNbinsX()); i++)
    params[i] = (*fgEfficiency)(i)*corrected->GetBinContent(i+1);
  
  TH1* penalty = GetPenaltyPlot(params);
  
  delete[] params;
  
  return penalty;
}

//____________________________________________________________________
TH1* AliUnfolding::GetPenaltyPlot(Double_t* params)
{
  // draws the penalty factors as function of multiplicity of the current selected regularization
  
  //TH1* penalty = new TH1F("penalty", ";unfolded multiplicity;penalty factor", fgMaxParams, -0.5, fgMaxParams - 0.5);
  //  TH1* penalty = new TH1F("penalty", ";unfolded pos;penalty factor", fgMaxParams, fgUnfoldedAxis->GetBinCenter(0)-0.5*fgUnfoldedAxis->GetBinWidth(0),fgUnfoldedAxis->GetBinCenter(fgMaxParams)+0.5*fgUnfoldedAxis->GetBinWidth(fgMaxParams) );

  TH1* penalty = 0;
  if (fgUnfoldedAxis->GetXbins()->GetArray())
    penalty = new TH1F("penalty", ";unfolded pos;penalty factor", fgUnfoldedAxis->GetNbins(),fgUnfoldedAxis->GetXbins()->GetArray());
  else
    penalty = new TH1F("penalty", ";unfolded pos;penalty factor", fgUnfoldedAxis->GetNbins(),fgUnfoldedAxis->GetXmin(),fgUnfoldedAxis->GetXmax());

  for (Int_t i=1+fgSkipBinsBegin; i<fgMaxParams; ++i)
  {
    Double_t diff = 0;
    if (fgRegularizationType == kPol0)
    {
      Double_t right  = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
      Double_t left   = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);

      if (left != 0)
      {
        Double_t diffTmp = (right - left);
        diff = diffTmp * diffTmp / left / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2) / 100;
      }
    } 
    if (fgRegularizationType == kPol1 && i > 1+fgSkipBinsBegin) 
    {
      if (params[i-1] == 0)
        continue;

      Double_t right  = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
      Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
      Double_t left   = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);

      Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
      Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);

      diff = (der1 - der2) * (der1 - der2) / middle;
    }

    if (fgRegularizationType == kLog && i > 1+fgSkipBinsBegin) 
    {
      if (params[i-1] == 0)
        continue;
  
      Double_t right  = log(params[i]);
      Double_t middle = log(params[i-1]);
      Double_t left   = log(params[i-2]);
      
      Double_t der1 = (right - middle);
      Double_t der2 = (middle - left);
  
      //Double_t error = 1. / params[i] + 4. / params[i-1] + 1. / params[i-2];
      //Printf("%d %f %f", i, (der1 - der2) * (der1 - der2), error);
       
      diff = (der1 - der2) * (der1 - der2);// / error;
    }
    if (fgRegularizationType == kCurvature && i > 1+fgSkipBinsBegin)
    {
      Double_t right  = params[i];    // params are sqrt
      Double_t middle = params[i-1];
      Double_t left   = params[i-2];
      
      Double_t der1 = (right - middle)/0.5/(fgUnfoldedAxis->GetBinWidth(i-1) + fgUnfoldedAxis->GetBinWidth(i));
      Double_t der2 = (middle - left)/0.5/(fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i+1));
      
      diff = (der1 - der2)/(fgUnfoldedAxis->GetBinWidth(i-1) + fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1))*3.0;
      diff = 1e4*diff*diff;
    }
    if (fgRegularizationType == kPowerLaw && i > 1+fgSkipBinsBegin) 
    {

      if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
	continue;
      
      if (fgUnfoldedAxis->GetBinWidth(i+1) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8)
	continue;
      
      double middle = TMath::Power(params[i-1] / fgUnfoldedAxis->GetBinWidth(i),fgPowern);

      if(middle>0) {
	double right  = TMath::Power(params[i] / fgUnfoldedAxis->GetBinWidth(i+1),fgPowern)/middle;
	
	double left   = TMath::Power(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1),fgPowern)/middle;
	
	middle = 1.;
	
	Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
	Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
	
	diff = (der1 - der2) * (der1 - der2);// / error;
      }
    }

    if (fgRegularizationType == kLogLog && i > 1+fgSkipBinsBegin) 
    {

      if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
	continue;

      Double_t right  = log(params[i] / (*fgEfficiency)(i) / fgUnfoldedAxis->GetBinWidth(i+1));
      Double_t middle = log(params[i-1] / (*fgEfficiency)(i-1) / fgUnfoldedAxis->GetBinWidth(i));
      Double_t left   = log(params[i-2] / (*fgEfficiency)(i-2) / fgUnfoldedAxis->GetBinWidth(i-1));
      
      Double_t der1 = (right - middle) / ( log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i)) );
      Double_t der2 = (middle - left) /( log(fgUnfoldedAxis->GetBinCenter(i)) - log(fgUnfoldedAxis->GetBinCenter(i-1)) );
      
      double tmp = (log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i-1)))/2.;
      Double_t dder = (der1-der2) / tmp;
      
      diff = dder * dder;
    }
    
    penalty->SetBinContent(i, diff*fgRegularizationWeight);
    
    //Printf("%d %f %f %f %f", i-1, left, middle, right, diff);
  }
  
  return penalty;
}

//____________________________________________________________________
void AliUnfolding::TF1Function(Int_t& unused1, Double_t* unused2, Double_t& chi2, Double_t *params, Int_t unused3)
{
  //
  // fit function for minuit
  // uses the TF1 stored in fgFitFunction
  //

  for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
    fgFitFunction->SetParameter(i, params[i]);

  Double_t* params2 = new Double_t[fgMaxParams];

  for (Int_t i=0; i<fgMaxParams; ++i)
    params2[i] = fgFitFunction->Eval(i);

  Chi2Function(unused1, unused2, chi2, params2, unused3);
  
  delete[] params2;

  if (fgDebug)
    Printf("%f", chi2);
}

//____________________________________________________________________
Int_t AliUnfolding::UnfoldWithFunction(TH2* correlation, TH1* efficiency, TH1* measured, TH1* /* initialConditions */, TH1* aResult)
{
  //
  // correct spectrum using minuit chi2 method applying a functional fit
  //
  
  if (!fgFitFunction)
  {
    Printf("AliUnfolding::UnfoldWithFunction: ERROR fit function not set. Exiting.");
    return -1;
  }    
  
  SetChi2Regularization(kNone, 0);
  
  SetStaticVariables(correlation, measured, efficiency);
  
  // Initialize TMinuit via generic fitter interface
  TVirtualFitter *minuit = TVirtualFitter::Fitter(0, fgFitFunction->GetNpar());

  minuit->SetFCN(TF1Function);
  for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
  {
    Double_t lower, upper;
    fgFitFunction->GetParLimits(i, lower, upper);
    minuit->SetParameter(i, Form("param%d",i), fgFitFunction->GetParameter(i), fgMinuitStepSize, lower, upper);
  }

  Double_t arglist[100];
  arglist[0] = 0;
  minuit->ExecuteCommand("SET PRINT", arglist, 1);
  minuit->ExecuteCommand("SCAN", arglist, 0);
  minuit->ExecuteCommand("MIGRAD", arglist, 0);
  //minuit->ExecuteCommand("MINOS", arglist, 0);

  for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
    fgFitFunction->SetParameter(i, minuit->GetParameter(i));

  for (Int_t i=0; i<fgMaxParams; ++i)
  {
    Double_t value = fgFitFunction->Eval(i);
    if (fgDebug)
      Printf("%d : %f", i, value);
    
    if (efficiency)
    {
      if (efficiency->GetBinContent(i+1) > 0)
      {
        value /= efficiency->GetBinContent(i+1);
      }
      else
        value = 0;
    }
    aResult->SetBinContent(i+1, value);
    aResult->SetBinError(i+1, 0);
  }
  
  return 0;
}

//____________________________________________________________________
void AliUnfolding::CreateOverflowBin(TH2* correlation, TH1* measured)
{
  // Finds the first bin where the content is below fgStatLimit and combines all values for this bin and larger bins
  // The same limit is applied to the correlation  
  
  Int_t lastBin = 0;
  for (Int_t i=1; i<=measured->GetNbinsX(); ++i)
  {
    if (measured->GetBinContent(i) <= fgOverflowBinLimit)
    {
      lastBin = i;
      break;
    }
  }
  
  if (lastBin == 0)
  {
    Printf("AliUnfolding::CreateOverflowBin: WARNING: First bin is already below limit of %f", fgOverflowBinLimit);
    return;
  }
  
  Printf("AliUnfolding::CreateOverflowBin: Bin limit in measured spectrum determined to be %d", lastBin);
  
  TCanvas* canvas = 0;

  if (fgDebug)
  {
    canvas = new TCanvas("StatSolution", "StatSolution", 1000, 800);
    canvas->Divide(2, 2);

    canvas->cd(1);
    measured->SetStats(kFALSE);
    measured->DrawCopy();
    gPad->SetLogy();

    canvas->cd(2);
    correlation->SetStats(kFALSE);
    correlation->DrawCopy("COLZ");
  }

  measured->SetBinContent(lastBin, measured->Integral(lastBin, measured->GetNbinsX()));
  for (Int_t i=lastBin+1; i<=measured->GetNbinsX(); ++i)
  {
    measured->SetBinContent(i, 0);
    measured->SetBinError(i, 0);
  }
  // the error is set to sqrt(N), better solution possible?, sum of relative errors of all contributions???
  measured->SetBinError(lastBin, TMath::Sqrt(measured->GetBinContent(lastBin)));

  Printf("AliUnfolding::CreateOverflowBin: This bin has now %f +- %f entries", measured->GetBinContent(lastBin), measured->GetBinError(lastBin));

  for (Int_t i=1; i<=correlation->GetNbinsX(); ++i)
  {
    correlation->SetBinContent(i, lastBin, correlation->Integral(i, i, lastBin, correlation->GetNbinsY()));
    // the error is set to sqrt(N), better solution possible?, sum of relative errors of all contributions???
    correlation->SetBinError(i, lastBin, TMath::Sqrt(correlation->GetBinContent(i, lastBin)));

    for (Int_t j=lastBin+1; j<=correlation->GetNbinsY(); ++j)
    {
      correlation->SetBinContent(i, j, 0);
      correlation->SetBinError(i, j, 0);
    }
  }

  Printf("AliUnfolding::CreateOverflowBin: Adjusted correlation matrix!");

  if (canvas)
  {
    canvas->cd(3);
    measured->DrawCopy();
    gPad->SetLogy();

    canvas->cd(4);
    correlation->DrawCopy("COLZ");
  }
}

Int_t AliUnfolding::UnfoldGetBias(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result)
{
  // unfolds and assigns bias as errors with Eq. 19 of Cowan, "a survey of unfolding methods for particle physics"
  // b_i = sum_j dmu_i/dn_j (nu_j - n_j) with nu_j as folded guess, n_j as data
  // dmu_i/dn_j is found numerically by changing the bin content and re-unfolding
  //
  // return code: 0 (success) -1 (error: from Unfold(...) )

  if (Unfold(correlation, efficiency, measured, initialConditions, result) != 0)
    return -1;

  TMatrixD matrix(fgMaxInput, fgMaxParams);
  
  TH1* newResult[4];
  for (Int_t loop=0; loop<4; loop++)
    newResult[loop] = (TH1*) result->Clone(Form("newresult_%d", loop));
    
  // change bin-by-bin and built matrix of effects
  for (Int_t m=0; m<fgMaxInput; m++)
  {
    if (measured->GetBinContent(m+1) < 1)
      continue;
      
    for (Int_t loop=0; loop<4; loop++)
    {
      //Printf("%d %d", i, loop);
      Float_t factor = 1;
      switch (loop)
      {
        case 0: factor = 0.5; break;
        case 1: factor = -0.5; break;
        case 2: factor = 1; break;
        case 3: factor = -1; break;
        default: return -1;
      }
      
      TH1* measuredClone = (TH1*) measured->Clone("measuredClone");
  
      measuredClone->SetBinContent(m+1, measured->GetBinContent(m+1) + factor * TMath::Sqrt(measured->GetBinContent(m+1)));
      //new TCanvas; measuredClone->Draw("COLZ");
    
      newResult[loop]->Reset();
      Unfold(correlation, efficiency, measuredClone, measuredClone, newResult[loop]);
      // WARNING if we leave here, then nothing is calculated
        //return -1;
      
      delete measuredClone;
    }
    
    for (Int_t t=0; t<fgMaxParams; t++)
    {
      // one value
      //matrix(i, j) = (result->GetBinContent(j+1) - newResult->GetBinContent(j+1)) / TMath::Sqrt(changedHist->GetBinContent(1, i+1));
      
      // four values from the derivate (procedure taken from ROOT -- suggestion by Ruben)
      // = 1/2D * [ 8 (f(D/2) - f(-D/2)) - (f(D)-f(-D)) ]/3
      
      /*
      // check formula
      measured->SetBinContent(1, m+1, 100);
      newResult[0]->SetBinContent(t+1, 5 * 0.5 + 10);
      newResult[1]->SetBinContent(t+1, 5 * -0.5 + 10);
      newResult[2]->SetBinContent(t+1, 5 * 1 + 10);
      newResult[3]->SetBinContent(t+1, 5 * -1 + 10);
      */
      
      matrix(m, t) = 0.5 / TMath::Sqrt(measured->GetBinContent(m+1)) * 
        ( 8. * (newResult[0]->GetBinContent(t+1) - newResult[1]->GetBinContent(t+1)) - 
              (newResult[2]->GetBinContent(t+1) - newResult[3]->GetBinContent(t+1))
        ) / 3;
    }
  }
  
  for (Int_t loop=0; loop<4; loop++)
    delete newResult[loop];
  
  // ...calculate folded guess  
  TH1* convoluted = (TH1*) measured->Clone("convoluted");
  convoluted->Reset();
  for (Int_t m=0; m<fgMaxInput; m++)
  {
    Float_t value = 0;
    for (Int_t t = 0; t<fgMaxParams; t++)
    {
      Float_t tmp = correlation->GetBinContent(t+1, m+1) * result->GetBinContent(t+1);
      if (efficiency)
        tmp *= efficiency->GetBinContent(t+1);
      value += tmp;
    }
    convoluted->SetBinContent(m+1, value);
  }
  
  // rescale
  convoluted->Scale(measured->Integral() / convoluted->Integral());
  
  //new TCanvas; convoluted->Draw(); measured->Draw("SAME"); measured->SetLineColor(2);
  //return;
  
  // difference
  convoluted->Add(measured, -1);
  
  // ...and the bias
  for (Int_t t = 0; t<fgMaxParams; t++)
  {
    Double_t error = 0;
    for (Int_t m=0; m<fgMaxInput; m++)
      error += matrix(m, t) * convoluted->GetBinContent(m+1); 
    result->SetBinError(t+1, error);
  }
  
  //new TCanvas; result->Draw(); gPad->SetLogy();
  
  return 0;
}
 AliUnfolding.cxx:1
 AliUnfolding.cxx:2
 AliUnfolding.cxx:3
 AliUnfolding.cxx:4
 AliUnfolding.cxx:5
 AliUnfolding.cxx:6
 AliUnfolding.cxx:7
 AliUnfolding.cxx:8
 AliUnfolding.cxx:9
 AliUnfolding.cxx:10
 AliUnfolding.cxx:11
 AliUnfolding.cxx:12
 AliUnfolding.cxx:13
 AliUnfolding.cxx:14
 AliUnfolding.cxx:15
 AliUnfolding.cxx:16
 AliUnfolding.cxx:17
 AliUnfolding.cxx:18
 AliUnfolding.cxx:19
 AliUnfolding.cxx:20
 AliUnfolding.cxx:21
 AliUnfolding.cxx:22
 AliUnfolding.cxx:23
 AliUnfolding.cxx:24
 AliUnfolding.cxx:25
 AliUnfolding.cxx:26
 AliUnfolding.cxx:27
 AliUnfolding.cxx:28
 AliUnfolding.cxx:29
 AliUnfolding.cxx:30
 AliUnfolding.cxx:31
 AliUnfolding.cxx:32
 AliUnfolding.cxx:33
 AliUnfolding.cxx:34
 AliUnfolding.cxx:35
 AliUnfolding.cxx:36
 AliUnfolding.cxx:37
 AliUnfolding.cxx:38
 AliUnfolding.cxx:39
 AliUnfolding.cxx:40
 AliUnfolding.cxx:41
 AliUnfolding.cxx:42
 AliUnfolding.cxx:43
 AliUnfolding.cxx:44
 AliUnfolding.cxx:45
 AliUnfolding.cxx:46
 AliUnfolding.cxx:47
 AliUnfolding.cxx:48
 AliUnfolding.cxx:49
 AliUnfolding.cxx:50
 AliUnfolding.cxx:51
 AliUnfolding.cxx:52
 AliUnfolding.cxx:53
 AliUnfolding.cxx:54
 AliUnfolding.cxx:55
 AliUnfolding.cxx:56
 AliUnfolding.cxx:57
 AliUnfolding.cxx:58
 AliUnfolding.cxx:59
 AliUnfolding.cxx:60
 AliUnfolding.cxx:61
 AliUnfolding.cxx:62
 AliUnfolding.cxx:63
 AliUnfolding.cxx:64
 AliUnfolding.cxx:65
 AliUnfolding.cxx:66
 AliUnfolding.cxx:67
 AliUnfolding.cxx:68
 AliUnfolding.cxx:69
 AliUnfolding.cxx:70
 AliUnfolding.cxx:71
 AliUnfolding.cxx:72
 AliUnfolding.cxx:73
 AliUnfolding.cxx:74
 AliUnfolding.cxx:75
 AliUnfolding.cxx:76
 AliUnfolding.cxx:77
 AliUnfolding.cxx:78
 AliUnfolding.cxx:79
 AliUnfolding.cxx:80
 AliUnfolding.cxx:81
 AliUnfolding.cxx:82
 AliUnfolding.cxx:83
 AliUnfolding.cxx:84
 AliUnfolding.cxx:85
 AliUnfolding.cxx:86
 AliUnfolding.cxx:87
 AliUnfolding.cxx:88
 AliUnfolding.cxx:89
 AliUnfolding.cxx:90
 AliUnfolding.cxx:91
 AliUnfolding.cxx:92
 AliUnfolding.cxx:93
 AliUnfolding.cxx:94
 AliUnfolding.cxx:95
 AliUnfolding.cxx:96
 AliUnfolding.cxx:97
 AliUnfolding.cxx:98
 AliUnfolding.cxx:99
 AliUnfolding.cxx:100
 AliUnfolding.cxx:101
 AliUnfolding.cxx:102
 AliUnfolding.cxx:103
 AliUnfolding.cxx:104
 AliUnfolding.cxx:105
 AliUnfolding.cxx:106
 AliUnfolding.cxx:107
 AliUnfolding.cxx:108
 AliUnfolding.cxx:109
 AliUnfolding.cxx:110
 AliUnfolding.cxx:111
 AliUnfolding.cxx:112
 AliUnfolding.cxx:113
 AliUnfolding.cxx:114
 AliUnfolding.cxx:115
 AliUnfolding.cxx:116
 AliUnfolding.cxx:117
 AliUnfolding.cxx:118
 AliUnfolding.cxx:119
 AliUnfolding.cxx:120
 AliUnfolding.cxx:121
 AliUnfolding.cxx:122
 AliUnfolding.cxx:123
 AliUnfolding.cxx:124
 AliUnfolding.cxx:125
 AliUnfolding.cxx:126
 AliUnfolding.cxx:127
 AliUnfolding.cxx:128
 AliUnfolding.cxx:129
 AliUnfolding.cxx:130
 AliUnfolding.cxx:131
 AliUnfolding.cxx:132
 AliUnfolding.cxx:133
 AliUnfolding.cxx:134
 AliUnfolding.cxx:135
 AliUnfolding.cxx:136
 AliUnfolding.cxx:137
 AliUnfolding.cxx:138
 AliUnfolding.cxx:139
 AliUnfolding.cxx:140
 AliUnfolding.cxx:141
 AliUnfolding.cxx:142
 AliUnfolding.cxx:143
 AliUnfolding.cxx:144
 AliUnfolding.cxx:145
 AliUnfolding.cxx:146
 AliUnfolding.cxx:147
 AliUnfolding.cxx:148
 AliUnfolding.cxx:149
 AliUnfolding.cxx:150
 AliUnfolding.cxx:151
 AliUnfolding.cxx:152
 AliUnfolding.cxx:153
 AliUnfolding.cxx:154
 AliUnfolding.cxx:155
 AliUnfolding.cxx:156
 AliUnfolding.cxx:157
 AliUnfolding.cxx:158
 AliUnfolding.cxx:159
 AliUnfolding.cxx:160
 AliUnfolding.cxx:161
 AliUnfolding.cxx:162
 AliUnfolding.cxx:163
 AliUnfolding.cxx:164
 AliUnfolding.cxx:165
 AliUnfolding.cxx:166
 AliUnfolding.cxx:167
 AliUnfolding.cxx:168
 AliUnfolding.cxx:169
 AliUnfolding.cxx:170
 AliUnfolding.cxx:171
 AliUnfolding.cxx:172
 AliUnfolding.cxx:173
 AliUnfolding.cxx:174
 AliUnfolding.cxx:175
 AliUnfolding.cxx:176
 AliUnfolding.cxx:177
 AliUnfolding.cxx:178
 AliUnfolding.cxx:179
 AliUnfolding.cxx:180
 AliUnfolding.cxx:181
 AliUnfolding.cxx:182
 AliUnfolding.cxx:183
 AliUnfolding.cxx:184
 AliUnfolding.cxx:185
 AliUnfolding.cxx:186
 AliUnfolding.cxx:187
 AliUnfolding.cxx:188
 AliUnfolding.cxx:189
 AliUnfolding.cxx:190
 AliUnfolding.cxx:191
 AliUnfolding.cxx:192
 AliUnfolding.cxx:193
 AliUnfolding.cxx:194
 AliUnfolding.cxx:195
 AliUnfolding.cxx:196
 AliUnfolding.cxx:197
 AliUnfolding.cxx:198
 AliUnfolding.cxx:199
 AliUnfolding.cxx:200
 AliUnfolding.cxx:201
 AliUnfolding.cxx:202
 AliUnfolding.cxx:203
 AliUnfolding.cxx:204
 AliUnfolding.cxx:205
 AliUnfolding.cxx:206
 AliUnfolding.cxx:207
 AliUnfolding.cxx:208
 AliUnfolding.cxx:209
 AliUnfolding.cxx:210
 AliUnfolding.cxx:211
 AliUnfolding.cxx:212
 AliUnfolding.cxx:213
 AliUnfolding.cxx:214
 AliUnfolding.cxx:215
 AliUnfolding.cxx:216
 AliUnfolding.cxx:217
 AliUnfolding.cxx:218
 AliUnfolding.cxx:219
 AliUnfolding.cxx:220
 AliUnfolding.cxx:221
 AliUnfolding.cxx:222
 AliUnfolding.cxx:223
 AliUnfolding.cxx:224
 AliUnfolding.cxx:225
 AliUnfolding.cxx:226
 AliUnfolding.cxx:227
 AliUnfolding.cxx:228
 AliUnfolding.cxx:229
 AliUnfolding.cxx:230
 AliUnfolding.cxx:231
 AliUnfolding.cxx:232
 AliUnfolding.cxx:233
 AliUnfolding.cxx:234
 AliUnfolding.cxx:235
 AliUnfolding.cxx:236
 AliUnfolding.cxx:237
 AliUnfolding.cxx:238
 AliUnfolding.cxx:239
 AliUnfolding.cxx:240
 AliUnfolding.cxx:241
 AliUnfolding.cxx:242
 AliUnfolding.cxx:243
 AliUnfolding.cxx:244
 AliUnfolding.cxx:245
 AliUnfolding.cxx:246
 AliUnfolding.cxx:247
 AliUnfolding.cxx:248
 AliUnfolding.cxx:249
 AliUnfolding.cxx:250
 AliUnfolding.cxx:251
 AliUnfolding.cxx:252
 AliUnfolding.cxx:253
 AliUnfolding.cxx:254
 AliUnfolding.cxx:255
 AliUnfolding.cxx:256
 AliUnfolding.cxx:257
 AliUnfolding.cxx:258
 AliUnfolding.cxx:259
 AliUnfolding.cxx:260
 AliUnfolding.cxx:261
 AliUnfolding.cxx:262
 AliUnfolding.cxx:263
 AliUnfolding.cxx:264
 AliUnfolding.cxx:265
 AliUnfolding.cxx:266
 AliUnfolding.cxx:267
 AliUnfolding.cxx:268
 AliUnfolding.cxx:269
 AliUnfolding.cxx:270
 AliUnfolding.cxx:271
 AliUnfolding.cxx:272
 AliUnfolding.cxx:273
 AliUnfolding.cxx:274
 AliUnfolding.cxx:275
 AliUnfolding.cxx:276
 AliUnfolding.cxx:277
 AliUnfolding.cxx:278
 AliUnfolding.cxx:279
 AliUnfolding.cxx:280
 AliUnfolding.cxx:281
 AliUnfolding.cxx:282
 AliUnfolding.cxx:283
 AliUnfolding.cxx:284
 AliUnfolding.cxx:285
 AliUnfolding.cxx:286
 AliUnfolding.cxx:287
 AliUnfolding.cxx:288
 AliUnfolding.cxx:289
 AliUnfolding.cxx:290
 AliUnfolding.cxx:291
 AliUnfolding.cxx:292
 AliUnfolding.cxx:293
 AliUnfolding.cxx:294
 AliUnfolding.cxx:295
 AliUnfolding.cxx:296
 AliUnfolding.cxx:297
 AliUnfolding.cxx:298
 AliUnfolding.cxx:299
 AliUnfolding.cxx:300
 AliUnfolding.cxx:301
 AliUnfolding.cxx:302
 AliUnfolding.cxx:303
 AliUnfolding.cxx:304
 AliUnfolding.cxx:305
 AliUnfolding.cxx:306
 AliUnfolding.cxx:307
 AliUnfolding.cxx:308
 AliUnfolding.cxx:309
 AliUnfolding.cxx:310
 AliUnfolding.cxx:311
 AliUnfolding.cxx:312
 AliUnfolding.cxx:313
 AliUnfolding.cxx:314
 AliUnfolding.cxx:315
 AliUnfolding.cxx:316
 AliUnfolding.cxx:317
 AliUnfolding.cxx:318
 AliUnfolding.cxx:319
 AliUnfolding.cxx:320
 AliUnfolding.cxx:321
 AliUnfolding.cxx:322
 AliUnfolding.cxx:323
 AliUnfolding.cxx:324
 AliUnfolding.cxx:325
 AliUnfolding.cxx:326
 AliUnfolding.cxx:327
 AliUnfolding.cxx:328
 AliUnfolding.cxx:329
 AliUnfolding.cxx:330
 AliUnfolding.cxx:331
 AliUnfolding.cxx:332
 AliUnfolding.cxx:333
 AliUnfolding.cxx:334
 AliUnfolding.cxx:335
 AliUnfolding.cxx:336
 AliUnfolding.cxx:337
 AliUnfolding.cxx:338
 AliUnfolding.cxx:339
 AliUnfolding.cxx:340
 AliUnfolding.cxx:341
 AliUnfolding.cxx:342
 AliUnfolding.cxx:343
 AliUnfolding.cxx:344
 AliUnfolding.cxx:345
 AliUnfolding.cxx:346
 AliUnfolding.cxx:347
 AliUnfolding.cxx:348
 AliUnfolding.cxx:349
 AliUnfolding.cxx:350
 AliUnfolding.cxx:351
 AliUnfolding.cxx:352
 AliUnfolding.cxx:353
 AliUnfolding.cxx:354
 AliUnfolding.cxx:355
 AliUnfolding.cxx:356
 AliUnfolding.cxx:357
 AliUnfolding.cxx:358
 AliUnfolding.cxx:359
 AliUnfolding.cxx:360
 AliUnfolding.cxx:361
 AliUnfolding.cxx:362
 AliUnfolding.cxx:363
 AliUnfolding.cxx:364
 AliUnfolding.cxx:365
 AliUnfolding.cxx:366
 AliUnfolding.cxx:367
 AliUnfolding.cxx:368
 AliUnfolding.cxx:369
 AliUnfolding.cxx:370
 AliUnfolding.cxx:371
 AliUnfolding.cxx:372
 AliUnfolding.cxx:373
 AliUnfolding.cxx:374
 AliUnfolding.cxx:375
 AliUnfolding.cxx:376
 AliUnfolding.cxx:377
 AliUnfolding.cxx:378
 AliUnfolding.cxx:379
 AliUnfolding.cxx:380
 AliUnfolding.cxx:381
 AliUnfolding.cxx:382
 AliUnfolding.cxx:383
 AliUnfolding.cxx:384
 AliUnfolding.cxx:385
 AliUnfolding.cxx:386
 AliUnfolding.cxx:387
 AliUnfolding.cxx:388
 AliUnfolding.cxx:389
 AliUnfolding.cxx:390
 AliUnfolding.cxx:391
 AliUnfolding.cxx:392
 AliUnfolding.cxx:393
 AliUnfolding.cxx:394
 AliUnfolding.cxx:395
 AliUnfolding.cxx:396
 AliUnfolding.cxx:397
 AliUnfolding.cxx:398
 AliUnfolding.cxx:399
 AliUnfolding.cxx:400
 AliUnfolding.cxx:401
 AliUnfolding.cxx:402
 AliUnfolding.cxx:403
 AliUnfolding.cxx:404
 AliUnfolding.cxx:405
 AliUnfolding.cxx:406
 AliUnfolding.cxx:407
 AliUnfolding.cxx:408
 AliUnfolding.cxx:409
 AliUnfolding.cxx:410
 AliUnfolding.cxx:411
 AliUnfolding.cxx:412
 AliUnfolding.cxx:413
 AliUnfolding.cxx:414
 AliUnfolding.cxx:415
 AliUnfolding.cxx:416
 AliUnfolding.cxx:417
 AliUnfolding.cxx:418
 AliUnfolding.cxx:419
 AliUnfolding.cxx:420
 AliUnfolding.cxx:421
 AliUnfolding.cxx:422
 AliUnfolding.cxx:423
 AliUnfolding.cxx:424
 AliUnfolding.cxx:425
 AliUnfolding.cxx:426
 AliUnfolding.cxx:427
 AliUnfolding.cxx:428
 AliUnfolding.cxx:429
 AliUnfolding.cxx:430
 AliUnfolding.cxx:431
 AliUnfolding.cxx:432
 AliUnfolding.cxx:433
 AliUnfolding.cxx:434
 AliUnfolding.cxx:435
 AliUnfolding.cxx:436
 AliUnfolding.cxx:437
 AliUnfolding.cxx:438
 AliUnfolding.cxx:439
 AliUnfolding.cxx:440
 AliUnfolding.cxx:441
 AliUnfolding.cxx:442
 AliUnfolding.cxx:443
 AliUnfolding.cxx:444
 AliUnfolding.cxx:445
 AliUnfolding.cxx:446
 AliUnfolding.cxx:447
 AliUnfolding.cxx:448
 AliUnfolding.cxx:449
 AliUnfolding.cxx:450
 AliUnfolding.cxx:451
 AliUnfolding.cxx:452
 AliUnfolding.cxx:453
 AliUnfolding.cxx:454
 AliUnfolding.cxx:455
 AliUnfolding.cxx:456
 AliUnfolding.cxx:457
 AliUnfolding.cxx:458
 AliUnfolding.cxx:459
 AliUnfolding.cxx:460
 AliUnfolding.cxx:461
 AliUnfolding.cxx:462
 AliUnfolding.cxx:463
 AliUnfolding.cxx:464
 AliUnfolding.cxx:465
 AliUnfolding.cxx:466
 AliUnfolding.cxx:467
 AliUnfolding.cxx:468
 AliUnfolding.cxx:469
 AliUnfolding.cxx:470
 AliUnfolding.cxx:471
 AliUnfolding.cxx:472
 AliUnfolding.cxx:473
 AliUnfolding.cxx:474
 AliUnfolding.cxx:475
 AliUnfolding.cxx:476
 AliUnfolding.cxx:477
 AliUnfolding.cxx:478
 AliUnfolding.cxx:479
 AliUnfolding.cxx:480
 AliUnfolding.cxx:481
 AliUnfolding.cxx:482
 AliUnfolding.cxx:483
 AliUnfolding.cxx:484
 AliUnfolding.cxx:485
 AliUnfolding.cxx:486
 AliUnfolding.cxx:487
 AliUnfolding.cxx:488
 AliUnfolding.cxx:489
 AliUnfolding.cxx:490
 AliUnfolding.cxx:491
 AliUnfolding.cxx:492
 AliUnfolding.cxx:493
 AliUnfolding.cxx:494
 AliUnfolding.cxx:495
 AliUnfolding.cxx:496
 AliUnfolding.cxx:497
 AliUnfolding.cxx:498
 AliUnfolding.cxx:499
 AliUnfolding.cxx:500
 AliUnfolding.cxx:501
 AliUnfolding.cxx:502
 AliUnfolding.cxx:503
 AliUnfolding.cxx:504
 AliUnfolding.cxx:505
 AliUnfolding.cxx:506
 AliUnfolding.cxx:507
 AliUnfolding.cxx:508
 AliUnfolding.cxx:509
 AliUnfolding.cxx:510
 AliUnfolding.cxx:511
 AliUnfolding.cxx:512
 AliUnfolding.cxx:513
 AliUnfolding.cxx:514
 AliUnfolding.cxx:515
 AliUnfolding.cxx:516
 AliUnfolding.cxx:517
 AliUnfolding.cxx:518
 AliUnfolding.cxx:519
 AliUnfolding.cxx:520
 AliUnfolding.cxx:521
 AliUnfolding.cxx:522
 AliUnfolding.cxx:523
 AliUnfolding.cxx:524
 AliUnfolding.cxx:525
 AliUnfolding.cxx:526
 AliUnfolding.cxx:527
 AliUnfolding.cxx:528
 AliUnfolding.cxx:529
 AliUnfolding.cxx:530
 AliUnfolding.cxx:531
 AliUnfolding.cxx:532
 AliUnfolding.cxx:533
 AliUnfolding.cxx:534
 AliUnfolding.cxx:535
 AliUnfolding.cxx:536
 AliUnfolding.cxx:537
 AliUnfolding.cxx:538
 AliUnfolding.cxx:539
 AliUnfolding.cxx:540
 AliUnfolding.cxx:541
 AliUnfolding.cxx:542
 AliUnfolding.cxx:543
 AliUnfolding.cxx:544
 AliUnfolding.cxx:545
 AliUnfolding.cxx:546
 AliUnfolding.cxx:547
 AliUnfolding.cxx:548
 AliUnfolding.cxx:549
 AliUnfolding.cxx:550
 AliUnfolding.cxx:551
 AliUnfolding.cxx:552
 AliUnfolding.cxx:553
 AliUnfolding.cxx:554
 AliUnfolding.cxx:555
 AliUnfolding.cxx:556
 AliUnfolding.cxx:557
 AliUnfolding.cxx:558
 AliUnfolding.cxx:559
 AliUnfolding.cxx:560
 AliUnfolding.cxx:561
 AliUnfolding.cxx:562
 AliUnfolding.cxx:563
 AliUnfolding.cxx:564
 AliUnfolding.cxx:565
 AliUnfolding.cxx:566
 AliUnfolding.cxx:567
 AliUnfolding.cxx:568
 AliUnfolding.cxx:569
 AliUnfolding.cxx:570
 AliUnfolding.cxx:571
 AliUnfolding.cxx:572
 AliUnfolding.cxx:573
 AliUnfolding.cxx:574
 AliUnfolding.cxx:575
 AliUnfolding.cxx:576
 AliUnfolding.cxx:577
 AliUnfolding.cxx:578
 AliUnfolding.cxx:579
 AliUnfolding.cxx:580
 AliUnfolding.cxx:581
 AliUnfolding.cxx:582
 AliUnfolding.cxx:583
 AliUnfolding.cxx:584
 AliUnfolding.cxx:585
 AliUnfolding.cxx:586
 AliUnfolding.cxx:587
 AliUnfolding.cxx:588
 AliUnfolding.cxx:589
 AliUnfolding.cxx:590
 AliUnfolding.cxx:591
 AliUnfolding.cxx:592
 AliUnfolding.cxx:593
 AliUnfolding.cxx:594
 AliUnfolding.cxx:595
 AliUnfolding.cxx:596
 AliUnfolding.cxx:597
 AliUnfolding.cxx:598
 AliUnfolding.cxx:599
 AliUnfolding.cxx:600
 AliUnfolding.cxx:601
 AliUnfolding.cxx:602
 AliUnfolding.cxx:603
 AliUnfolding.cxx:604
 AliUnfolding.cxx:605
 AliUnfolding.cxx:606
 AliUnfolding.cxx:607
 AliUnfolding.cxx:608
 AliUnfolding.cxx:609
 AliUnfolding.cxx:610
 AliUnfolding.cxx:611
 AliUnfolding.cxx:612
 AliUnfolding.cxx:613
 AliUnfolding.cxx:614
 AliUnfolding.cxx:615
 AliUnfolding.cxx:616
 AliUnfolding.cxx:617
 AliUnfolding.cxx:618
 AliUnfolding.cxx:619
 AliUnfolding.cxx:620
 AliUnfolding.cxx:621
 AliUnfolding.cxx:622
 AliUnfolding.cxx:623
 AliUnfolding.cxx:624
 AliUnfolding.cxx:625
 AliUnfolding.cxx:626
 AliUnfolding.cxx:627
 AliUnfolding.cxx:628
 AliUnfolding.cxx:629
 AliUnfolding.cxx:630
 AliUnfolding.cxx:631
 AliUnfolding.cxx:632
 AliUnfolding.cxx:633
 AliUnfolding.cxx:634
 AliUnfolding.cxx:635
 AliUnfolding.cxx:636
 AliUnfolding.cxx:637
 AliUnfolding.cxx:638
 AliUnfolding.cxx:639
 AliUnfolding.cxx:640
 AliUnfolding.cxx:641
 AliUnfolding.cxx:642
 AliUnfolding.cxx:643
 AliUnfolding.cxx:644
 AliUnfolding.cxx:645
 AliUnfolding.cxx:646
 AliUnfolding.cxx:647
 AliUnfolding.cxx:648
 AliUnfolding.cxx:649
 AliUnfolding.cxx:650
 AliUnfolding.cxx:651
 AliUnfolding.cxx:652
 AliUnfolding.cxx:653
 AliUnfolding.cxx:654
 AliUnfolding.cxx:655
 AliUnfolding.cxx:656
 AliUnfolding.cxx:657
 AliUnfolding.cxx:658
 AliUnfolding.cxx:659
 AliUnfolding.cxx:660
 AliUnfolding.cxx:661
 AliUnfolding.cxx:662
 AliUnfolding.cxx:663
 AliUnfolding.cxx:664
 AliUnfolding.cxx:665
 AliUnfolding.cxx:666
 AliUnfolding.cxx:667
 AliUnfolding.cxx:668
 AliUnfolding.cxx:669
 AliUnfolding.cxx:670
 AliUnfolding.cxx:671
 AliUnfolding.cxx:672
 AliUnfolding.cxx:673
 AliUnfolding.cxx:674
 AliUnfolding.cxx:675
 AliUnfolding.cxx:676
 AliUnfolding.cxx:677
 AliUnfolding.cxx:678
 AliUnfolding.cxx:679
 AliUnfolding.cxx:680
 AliUnfolding.cxx:681
 AliUnfolding.cxx:682
 AliUnfolding.cxx:683
 AliUnfolding.cxx:684
 AliUnfolding.cxx:685
 AliUnfolding.cxx:686
 AliUnfolding.cxx:687
 AliUnfolding.cxx:688
 AliUnfolding.cxx:689
 AliUnfolding.cxx:690
 AliUnfolding.cxx:691
 AliUnfolding.cxx:692
 AliUnfolding.cxx:693
 AliUnfolding.cxx:694
 AliUnfolding.cxx:695
 AliUnfolding.cxx:696
 AliUnfolding.cxx:697
 AliUnfolding.cxx:698
 AliUnfolding.cxx:699
 AliUnfolding.cxx:700
 AliUnfolding.cxx:701
 AliUnfolding.cxx:702
 AliUnfolding.cxx:703
 AliUnfolding.cxx:704
 AliUnfolding.cxx:705
 AliUnfolding.cxx:706
 AliUnfolding.cxx:707
 AliUnfolding.cxx:708
 AliUnfolding.cxx:709
 AliUnfolding.cxx:710
 AliUnfolding.cxx:711
 AliUnfolding.cxx:712
 AliUnfolding.cxx:713
 AliUnfolding.cxx:714
 AliUnfolding.cxx:715
 AliUnfolding.cxx:716
 AliUnfolding.cxx:717
 AliUnfolding.cxx:718
 AliUnfolding.cxx:719
 AliUnfolding.cxx:720
 AliUnfolding.cxx:721
 AliUnfolding.cxx:722
 AliUnfolding.cxx:723
 AliUnfolding.cxx:724
 AliUnfolding.cxx:725
 AliUnfolding.cxx:726
 AliUnfolding.cxx:727
 AliUnfolding.cxx:728
 AliUnfolding.cxx:729
 AliUnfolding.cxx:730
 AliUnfolding.cxx:731
 AliUnfolding.cxx:732
 AliUnfolding.cxx:733
 AliUnfolding.cxx:734
 AliUnfolding.cxx:735
 AliUnfolding.cxx:736
 AliUnfolding.cxx:737
 AliUnfolding.cxx:738
 AliUnfolding.cxx:739
 AliUnfolding.cxx:740
 AliUnfolding.cxx:741
 AliUnfolding.cxx:742
 AliUnfolding.cxx:743
 AliUnfolding.cxx:744
 AliUnfolding.cxx:745
 AliUnfolding.cxx:746
 AliUnfolding.cxx:747
 AliUnfolding.cxx:748
 AliUnfolding.cxx:749
 AliUnfolding.cxx:750
 AliUnfolding.cxx:751
 AliUnfolding.cxx:752
 AliUnfolding.cxx:753
 AliUnfolding.cxx:754
 AliUnfolding.cxx:755
 AliUnfolding.cxx:756
 AliUnfolding.cxx:757
 AliUnfolding.cxx:758
 AliUnfolding.cxx:759
 AliUnfolding.cxx:760
 AliUnfolding.cxx:761
 AliUnfolding.cxx:762
 AliUnfolding.cxx:763
 AliUnfolding.cxx:764
 AliUnfolding.cxx:765
 AliUnfolding.cxx:766
 AliUnfolding.cxx:767
 AliUnfolding.cxx:768
 AliUnfolding.cxx:769
 AliUnfolding.cxx:770
 AliUnfolding.cxx:771
 AliUnfolding.cxx:772
 AliUnfolding.cxx:773
 AliUnfolding.cxx:774
 AliUnfolding.cxx:775
 AliUnfolding.cxx:776
 AliUnfolding.cxx:777
 AliUnfolding.cxx:778
 AliUnfolding.cxx:779
 AliUnfolding.cxx:780
 AliUnfolding.cxx:781
 AliUnfolding.cxx:782
 AliUnfolding.cxx:783
 AliUnfolding.cxx:784
 AliUnfolding.cxx:785
 AliUnfolding.cxx:786
 AliUnfolding.cxx:787
 AliUnfolding.cxx:788
 AliUnfolding.cxx:789
 AliUnfolding.cxx:790
 AliUnfolding.cxx:791
 AliUnfolding.cxx:792
 AliUnfolding.cxx:793
 AliUnfolding.cxx:794
 AliUnfolding.cxx:795
 AliUnfolding.cxx:796
 AliUnfolding.cxx:797
 AliUnfolding.cxx:798
 AliUnfolding.cxx:799
 AliUnfolding.cxx:800
 AliUnfolding.cxx:801
 AliUnfolding.cxx:802
 AliUnfolding.cxx:803
 AliUnfolding.cxx:804
 AliUnfolding.cxx:805
 AliUnfolding.cxx:806
 AliUnfolding.cxx:807
 AliUnfolding.cxx:808
 AliUnfolding.cxx:809
 AliUnfolding.cxx:810
 AliUnfolding.cxx:811
 AliUnfolding.cxx:812
 AliUnfolding.cxx:813
 AliUnfolding.cxx:814
 AliUnfolding.cxx:815
 AliUnfolding.cxx:816
 AliUnfolding.cxx:817
 AliUnfolding.cxx:818
 AliUnfolding.cxx:819
 AliUnfolding.cxx:820
 AliUnfolding.cxx:821
 AliUnfolding.cxx:822
 AliUnfolding.cxx:823
 AliUnfolding.cxx:824
 AliUnfolding.cxx:825
 AliUnfolding.cxx:826
 AliUnfolding.cxx:827
 AliUnfolding.cxx:828
 AliUnfolding.cxx:829
 AliUnfolding.cxx:830
 AliUnfolding.cxx:831
 AliUnfolding.cxx:832
 AliUnfolding.cxx:833
 AliUnfolding.cxx:834
 AliUnfolding.cxx:835
 AliUnfolding.cxx:836
 AliUnfolding.cxx:837
 AliUnfolding.cxx:838
 AliUnfolding.cxx:839
 AliUnfolding.cxx:840
 AliUnfolding.cxx:841
 AliUnfolding.cxx:842
 AliUnfolding.cxx:843
 AliUnfolding.cxx:844
 AliUnfolding.cxx:845
 AliUnfolding.cxx:846
 AliUnfolding.cxx:847
 AliUnfolding.cxx:848
 AliUnfolding.cxx:849
 AliUnfolding.cxx:850
 AliUnfolding.cxx:851
 AliUnfolding.cxx:852
 AliUnfolding.cxx:853
 AliUnfolding.cxx:854
 AliUnfolding.cxx:855
 AliUnfolding.cxx:856
 AliUnfolding.cxx:857
 AliUnfolding.cxx:858
 AliUnfolding.cxx:859
 AliUnfolding.cxx:860
 AliUnfolding.cxx:861
 AliUnfolding.cxx:862
 AliUnfolding.cxx:863
 AliUnfolding.cxx:864
 AliUnfolding.cxx:865
 AliUnfolding.cxx:866
 AliUnfolding.cxx:867
 AliUnfolding.cxx:868
 AliUnfolding.cxx:869
 AliUnfolding.cxx:870
 AliUnfolding.cxx:871
 AliUnfolding.cxx:872
 AliUnfolding.cxx:873
 AliUnfolding.cxx:874
 AliUnfolding.cxx:875
 AliUnfolding.cxx:876
 AliUnfolding.cxx:877
 AliUnfolding.cxx:878
 AliUnfolding.cxx:879
 AliUnfolding.cxx:880
 AliUnfolding.cxx:881
 AliUnfolding.cxx:882
 AliUnfolding.cxx:883
 AliUnfolding.cxx:884
 AliUnfolding.cxx:885
 AliUnfolding.cxx:886
 AliUnfolding.cxx:887
 AliUnfolding.cxx:888
 AliUnfolding.cxx:889
 AliUnfolding.cxx:890
 AliUnfolding.cxx:891
 AliUnfolding.cxx:892
 AliUnfolding.cxx:893
 AliUnfolding.cxx:894
 AliUnfolding.cxx:895
 AliUnfolding.cxx:896
 AliUnfolding.cxx:897
 AliUnfolding.cxx:898
 AliUnfolding.cxx:899
 AliUnfolding.cxx:900
 AliUnfolding.cxx:901
 AliUnfolding.cxx:902
 AliUnfolding.cxx:903
 AliUnfolding.cxx:904
 AliUnfolding.cxx:905
 AliUnfolding.cxx:906
 AliUnfolding.cxx:907
 AliUnfolding.cxx:908
 AliUnfolding.cxx:909
 AliUnfolding.cxx:910
 AliUnfolding.cxx:911
 AliUnfolding.cxx:912
 AliUnfolding.cxx:913
 AliUnfolding.cxx:914
 AliUnfolding.cxx:915
 AliUnfolding.cxx:916
 AliUnfolding.cxx:917
 AliUnfolding.cxx:918
 AliUnfolding.cxx:919
 AliUnfolding.cxx:920
 AliUnfolding.cxx:921
 AliUnfolding.cxx:922
 AliUnfolding.cxx:923
 AliUnfolding.cxx:924
 AliUnfolding.cxx:925
 AliUnfolding.cxx:926
 AliUnfolding.cxx:927
 AliUnfolding.cxx:928
 AliUnfolding.cxx:929
 AliUnfolding.cxx:930
 AliUnfolding.cxx:931
 AliUnfolding.cxx:932
 AliUnfolding.cxx:933
 AliUnfolding.cxx:934
 AliUnfolding.cxx:935
 AliUnfolding.cxx:936
 AliUnfolding.cxx:937
 AliUnfolding.cxx:938
 AliUnfolding.cxx:939
 AliUnfolding.cxx:940
 AliUnfolding.cxx:941
 AliUnfolding.cxx:942
 AliUnfolding.cxx:943
 AliUnfolding.cxx:944
 AliUnfolding.cxx:945
 AliUnfolding.cxx:946
 AliUnfolding.cxx:947
 AliUnfolding.cxx:948
 AliUnfolding.cxx:949
 AliUnfolding.cxx:950
 AliUnfolding.cxx:951
 AliUnfolding.cxx:952
 AliUnfolding.cxx:953
 AliUnfolding.cxx:954
 AliUnfolding.cxx:955
 AliUnfolding.cxx:956
 AliUnfolding.cxx:957
 AliUnfolding.cxx:958
 AliUnfolding.cxx:959
 AliUnfolding.cxx:960
 AliUnfolding.cxx:961
 AliUnfolding.cxx:962
 AliUnfolding.cxx:963
 AliUnfolding.cxx:964
 AliUnfolding.cxx:965
 AliUnfolding.cxx:966
 AliUnfolding.cxx:967
 AliUnfolding.cxx:968
 AliUnfolding.cxx:969
 AliUnfolding.cxx:970
 AliUnfolding.cxx:971
 AliUnfolding.cxx:972
 AliUnfolding.cxx:973
 AliUnfolding.cxx:974
 AliUnfolding.cxx:975
 AliUnfolding.cxx:976
 AliUnfolding.cxx:977
 AliUnfolding.cxx:978
 AliUnfolding.cxx:979
 AliUnfolding.cxx:980
 AliUnfolding.cxx:981
 AliUnfolding.cxx:982
 AliUnfolding.cxx:983
 AliUnfolding.cxx:984
 AliUnfolding.cxx:985
 AliUnfolding.cxx:986
 AliUnfolding.cxx:987
 AliUnfolding.cxx:988
 AliUnfolding.cxx:989
 AliUnfolding.cxx:990
 AliUnfolding.cxx:991
 AliUnfolding.cxx:992
 AliUnfolding.cxx:993
 AliUnfolding.cxx:994
 AliUnfolding.cxx:995
 AliUnfolding.cxx:996
 AliUnfolding.cxx:997
 AliUnfolding.cxx:998
 AliUnfolding.cxx:999
 AliUnfolding.cxx:1000
 AliUnfolding.cxx:1001
 AliUnfolding.cxx:1002
 AliUnfolding.cxx:1003
 AliUnfolding.cxx:1004
 AliUnfolding.cxx:1005
 AliUnfolding.cxx:1006
 AliUnfolding.cxx:1007
 AliUnfolding.cxx:1008
 AliUnfolding.cxx:1009
 AliUnfolding.cxx:1010
 AliUnfolding.cxx:1011
 AliUnfolding.cxx:1012
 AliUnfolding.cxx:1013
 AliUnfolding.cxx:1014
 AliUnfolding.cxx:1015
 AliUnfolding.cxx:1016
 AliUnfolding.cxx:1017
 AliUnfolding.cxx:1018
 AliUnfolding.cxx:1019
 AliUnfolding.cxx:1020
 AliUnfolding.cxx:1021
 AliUnfolding.cxx:1022
 AliUnfolding.cxx:1023
 AliUnfolding.cxx:1024
 AliUnfolding.cxx:1025
 AliUnfolding.cxx:1026
 AliUnfolding.cxx:1027
 AliUnfolding.cxx:1028
 AliUnfolding.cxx:1029
 AliUnfolding.cxx:1030
 AliUnfolding.cxx:1031
 AliUnfolding.cxx:1032
 AliUnfolding.cxx:1033
 AliUnfolding.cxx:1034
 AliUnfolding.cxx:1035
 AliUnfolding.cxx:1036
 AliUnfolding.cxx:1037
 AliUnfolding.cxx:1038
 AliUnfolding.cxx:1039
 AliUnfolding.cxx:1040
 AliUnfolding.cxx:1041
 AliUnfolding.cxx:1042
 AliUnfolding.cxx:1043
 AliUnfolding.cxx:1044
 AliUnfolding.cxx:1045
 AliUnfolding.cxx:1046
 AliUnfolding.cxx:1047
 AliUnfolding.cxx:1048
 AliUnfolding.cxx:1049
 AliUnfolding.cxx:1050
 AliUnfolding.cxx:1051
 AliUnfolding.cxx:1052
 AliUnfolding.cxx:1053
 AliUnfolding.cxx:1054
 AliUnfolding.cxx:1055
 AliUnfolding.cxx:1056
 AliUnfolding.cxx:1057
 AliUnfolding.cxx:1058
 AliUnfolding.cxx:1059
 AliUnfolding.cxx:1060
 AliUnfolding.cxx:1061
 AliUnfolding.cxx:1062
 AliUnfolding.cxx:1063
 AliUnfolding.cxx:1064
 AliUnfolding.cxx:1065
 AliUnfolding.cxx:1066
 AliUnfolding.cxx:1067
 AliUnfolding.cxx:1068
 AliUnfolding.cxx:1069
 AliUnfolding.cxx:1070
 AliUnfolding.cxx:1071
 AliUnfolding.cxx:1072
 AliUnfolding.cxx:1073
 AliUnfolding.cxx:1074
 AliUnfolding.cxx:1075
 AliUnfolding.cxx:1076
 AliUnfolding.cxx:1077
 AliUnfolding.cxx:1078
 AliUnfolding.cxx:1079
 AliUnfolding.cxx:1080
 AliUnfolding.cxx:1081
 AliUnfolding.cxx:1082
 AliUnfolding.cxx:1083
 AliUnfolding.cxx:1084
 AliUnfolding.cxx:1085
 AliUnfolding.cxx:1086
 AliUnfolding.cxx:1087
 AliUnfolding.cxx:1088
 AliUnfolding.cxx:1089
 AliUnfolding.cxx:1090
 AliUnfolding.cxx:1091
 AliUnfolding.cxx:1092
 AliUnfolding.cxx:1093
 AliUnfolding.cxx:1094
 AliUnfolding.cxx:1095
 AliUnfolding.cxx:1096
 AliUnfolding.cxx:1097
 AliUnfolding.cxx:1098
 AliUnfolding.cxx:1099
 AliUnfolding.cxx:1100
 AliUnfolding.cxx:1101
 AliUnfolding.cxx:1102
 AliUnfolding.cxx:1103
 AliUnfolding.cxx:1104
 AliUnfolding.cxx:1105
 AliUnfolding.cxx:1106
 AliUnfolding.cxx:1107
 AliUnfolding.cxx:1108
 AliUnfolding.cxx:1109
 AliUnfolding.cxx:1110
 AliUnfolding.cxx:1111
 AliUnfolding.cxx:1112
 AliUnfolding.cxx:1113
 AliUnfolding.cxx:1114
 AliUnfolding.cxx:1115
 AliUnfolding.cxx:1116
 AliUnfolding.cxx:1117
 AliUnfolding.cxx:1118
 AliUnfolding.cxx:1119
 AliUnfolding.cxx:1120
 AliUnfolding.cxx:1121
 AliUnfolding.cxx:1122
 AliUnfolding.cxx:1123
 AliUnfolding.cxx:1124
 AliUnfolding.cxx:1125
 AliUnfolding.cxx:1126
 AliUnfolding.cxx:1127
 AliUnfolding.cxx:1128
 AliUnfolding.cxx:1129
 AliUnfolding.cxx:1130
 AliUnfolding.cxx:1131
 AliUnfolding.cxx:1132
 AliUnfolding.cxx:1133
 AliUnfolding.cxx:1134
 AliUnfolding.cxx:1135
 AliUnfolding.cxx:1136
 AliUnfolding.cxx:1137
 AliUnfolding.cxx:1138
 AliUnfolding.cxx:1139
 AliUnfolding.cxx:1140
 AliUnfolding.cxx:1141
 AliUnfolding.cxx:1142
 AliUnfolding.cxx:1143
 AliUnfolding.cxx:1144
 AliUnfolding.cxx:1145
 AliUnfolding.cxx:1146
 AliUnfolding.cxx:1147
 AliUnfolding.cxx:1148
 AliUnfolding.cxx:1149
 AliUnfolding.cxx:1150
 AliUnfolding.cxx:1151
 AliUnfolding.cxx:1152
 AliUnfolding.cxx:1153
 AliUnfolding.cxx:1154
 AliUnfolding.cxx:1155
 AliUnfolding.cxx:1156
 AliUnfolding.cxx:1157
 AliUnfolding.cxx:1158
 AliUnfolding.cxx:1159
 AliUnfolding.cxx:1160
 AliUnfolding.cxx:1161
 AliUnfolding.cxx:1162
 AliUnfolding.cxx:1163
 AliUnfolding.cxx:1164
 AliUnfolding.cxx:1165
 AliUnfolding.cxx:1166
 AliUnfolding.cxx:1167
 AliUnfolding.cxx:1168
 AliUnfolding.cxx:1169
 AliUnfolding.cxx:1170
 AliUnfolding.cxx:1171
 AliUnfolding.cxx:1172
 AliUnfolding.cxx:1173
 AliUnfolding.cxx:1174
 AliUnfolding.cxx:1175
 AliUnfolding.cxx:1176
 AliUnfolding.cxx:1177
 AliUnfolding.cxx:1178
 AliUnfolding.cxx:1179
 AliUnfolding.cxx:1180
 AliUnfolding.cxx:1181
 AliUnfolding.cxx:1182
 AliUnfolding.cxx:1183
 AliUnfolding.cxx:1184
 AliUnfolding.cxx:1185
 AliUnfolding.cxx:1186
 AliUnfolding.cxx:1187
 AliUnfolding.cxx:1188
 AliUnfolding.cxx:1189
 AliUnfolding.cxx:1190
 AliUnfolding.cxx:1191
 AliUnfolding.cxx:1192
 AliUnfolding.cxx:1193
 AliUnfolding.cxx:1194
 AliUnfolding.cxx:1195
 AliUnfolding.cxx:1196
 AliUnfolding.cxx:1197
 AliUnfolding.cxx:1198
 AliUnfolding.cxx:1199
 AliUnfolding.cxx:1200
 AliUnfolding.cxx:1201
 AliUnfolding.cxx:1202
 AliUnfolding.cxx:1203
 AliUnfolding.cxx:1204
 AliUnfolding.cxx:1205
 AliUnfolding.cxx:1206
 AliUnfolding.cxx:1207
 AliUnfolding.cxx:1208
 AliUnfolding.cxx:1209
 AliUnfolding.cxx:1210
 AliUnfolding.cxx:1211
 AliUnfolding.cxx:1212
 AliUnfolding.cxx:1213
 AliUnfolding.cxx:1214
 AliUnfolding.cxx:1215
 AliUnfolding.cxx:1216
 AliUnfolding.cxx:1217
 AliUnfolding.cxx:1218
 AliUnfolding.cxx:1219
 AliUnfolding.cxx:1220
 AliUnfolding.cxx:1221
 AliUnfolding.cxx:1222
 AliUnfolding.cxx:1223
 AliUnfolding.cxx:1224
 AliUnfolding.cxx:1225
 AliUnfolding.cxx:1226
 AliUnfolding.cxx:1227
 AliUnfolding.cxx:1228
 AliUnfolding.cxx:1229
 AliUnfolding.cxx:1230
 AliUnfolding.cxx:1231
 AliUnfolding.cxx:1232
 AliUnfolding.cxx:1233
 AliUnfolding.cxx:1234
 AliUnfolding.cxx:1235
 AliUnfolding.cxx:1236
 AliUnfolding.cxx:1237
 AliUnfolding.cxx:1238
 AliUnfolding.cxx:1239
 AliUnfolding.cxx:1240
 AliUnfolding.cxx:1241
 AliUnfolding.cxx:1242
 AliUnfolding.cxx:1243
 AliUnfolding.cxx:1244
 AliUnfolding.cxx:1245
 AliUnfolding.cxx:1246
 AliUnfolding.cxx:1247
 AliUnfolding.cxx:1248
 AliUnfolding.cxx:1249
 AliUnfolding.cxx:1250
 AliUnfolding.cxx:1251
 AliUnfolding.cxx:1252
 AliUnfolding.cxx:1253
 AliUnfolding.cxx:1254
 AliUnfolding.cxx:1255
 AliUnfolding.cxx:1256
 AliUnfolding.cxx:1257
 AliUnfolding.cxx:1258
 AliUnfolding.cxx:1259
 AliUnfolding.cxx:1260
 AliUnfolding.cxx:1261
 AliUnfolding.cxx:1262
 AliUnfolding.cxx:1263
 AliUnfolding.cxx:1264
 AliUnfolding.cxx:1265
 AliUnfolding.cxx:1266
 AliUnfolding.cxx:1267
 AliUnfolding.cxx:1268
 AliUnfolding.cxx:1269
 AliUnfolding.cxx:1270
 AliUnfolding.cxx:1271
 AliUnfolding.cxx:1272
 AliUnfolding.cxx:1273
 AliUnfolding.cxx:1274
 AliUnfolding.cxx:1275
 AliUnfolding.cxx:1276
 AliUnfolding.cxx:1277
 AliUnfolding.cxx:1278
 AliUnfolding.cxx:1279
 AliUnfolding.cxx:1280
 AliUnfolding.cxx:1281
 AliUnfolding.cxx:1282
 AliUnfolding.cxx:1283
 AliUnfolding.cxx:1284
 AliUnfolding.cxx:1285
 AliUnfolding.cxx:1286
 AliUnfolding.cxx:1287
 AliUnfolding.cxx:1288
 AliUnfolding.cxx:1289
 AliUnfolding.cxx:1290
 AliUnfolding.cxx:1291
 AliUnfolding.cxx:1292
 AliUnfolding.cxx:1293
 AliUnfolding.cxx:1294
 AliUnfolding.cxx:1295
 AliUnfolding.cxx:1296
 AliUnfolding.cxx:1297
 AliUnfolding.cxx:1298
 AliUnfolding.cxx:1299
 AliUnfolding.cxx:1300
 AliUnfolding.cxx:1301
 AliUnfolding.cxx:1302
 AliUnfolding.cxx:1303
 AliUnfolding.cxx:1304
 AliUnfolding.cxx:1305
 AliUnfolding.cxx:1306
 AliUnfolding.cxx:1307
 AliUnfolding.cxx:1308
 AliUnfolding.cxx:1309
 AliUnfolding.cxx:1310
 AliUnfolding.cxx:1311
 AliUnfolding.cxx:1312
 AliUnfolding.cxx:1313
 AliUnfolding.cxx:1314
 AliUnfolding.cxx:1315
 AliUnfolding.cxx:1316
 AliUnfolding.cxx:1317
 AliUnfolding.cxx:1318
 AliUnfolding.cxx:1319
 AliUnfolding.cxx:1320
 AliUnfolding.cxx:1321
 AliUnfolding.cxx:1322
 AliUnfolding.cxx:1323
 AliUnfolding.cxx:1324
 AliUnfolding.cxx:1325
 AliUnfolding.cxx:1326
 AliUnfolding.cxx:1327
 AliUnfolding.cxx:1328
 AliUnfolding.cxx:1329
 AliUnfolding.cxx:1330
 AliUnfolding.cxx:1331
 AliUnfolding.cxx:1332
 AliUnfolding.cxx:1333
 AliUnfolding.cxx:1334
 AliUnfolding.cxx:1335
 AliUnfolding.cxx:1336
 AliUnfolding.cxx:1337
 AliUnfolding.cxx:1338
 AliUnfolding.cxx:1339
 AliUnfolding.cxx:1340
 AliUnfolding.cxx:1341
 AliUnfolding.cxx:1342
 AliUnfolding.cxx:1343
 AliUnfolding.cxx:1344
 AliUnfolding.cxx:1345
 AliUnfolding.cxx:1346
 AliUnfolding.cxx:1347
 AliUnfolding.cxx:1348
 AliUnfolding.cxx:1349
 AliUnfolding.cxx:1350
 AliUnfolding.cxx:1351
 AliUnfolding.cxx:1352
 AliUnfolding.cxx:1353
 AliUnfolding.cxx:1354
 AliUnfolding.cxx:1355
 AliUnfolding.cxx:1356
 AliUnfolding.cxx:1357
 AliUnfolding.cxx:1358
 AliUnfolding.cxx:1359
 AliUnfolding.cxx:1360
 AliUnfolding.cxx:1361
 AliUnfolding.cxx:1362
 AliUnfolding.cxx:1363
 AliUnfolding.cxx:1364
 AliUnfolding.cxx:1365
 AliUnfolding.cxx:1366
 AliUnfolding.cxx:1367
 AliUnfolding.cxx:1368
 AliUnfolding.cxx:1369
 AliUnfolding.cxx:1370
 AliUnfolding.cxx:1371
 AliUnfolding.cxx:1372
 AliUnfolding.cxx:1373
 AliUnfolding.cxx:1374
 AliUnfolding.cxx:1375
 AliUnfolding.cxx:1376
 AliUnfolding.cxx:1377
 AliUnfolding.cxx:1378
 AliUnfolding.cxx:1379
 AliUnfolding.cxx:1380
 AliUnfolding.cxx:1381
 AliUnfolding.cxx:1382
 AliUnfolding.cxx:1383
 AliUnfolding.cxx:1384
 AliUnfolding.cxx:1385
 AliUnfolding.cxx:1386
 AliUnfolding.cxx:1387
 AliUnfolding.cxx:1388
 AliUnfolding.cxx:1389
 AliUnfolding.cxx:1390
 AliUnfolding.cxx:1391
 AliUnfolding.cxx:1392
 AliUnfolding.cxx:1393
 AliUnfolding.cxx:1394
 AliUnfolding.cxx:1395
 AliUnfolding.cxx:1396
 AliUnfolding.cxx:1397
 AliUnfolding.cxx:1398
 AliUnfolding.cxx:1399
 AliUnfolding.cxx:1400
 AliUnfolding.cxx:1401
 AliUnfolding.cxx:1402
 AliUnfolding.cxx:1403
 AliUnfolding.cxx:1404
 AliUnfolding.cxx:1405
 AliUnfolding.cxx:1406
 AliUnfolding.cxx:1407
 AliUnfolding.cxx:1408
 AliUnfolding.cxx:1409
 AliUnfolding.cxx:1410
 AliUnfolding.cxx:1411
 AliUnfolding.cxx:1412
 AliUnfolding.cxx:1413
 AliUnfolding.cxx:1414
 AliUnfolding.cxx:1415
 AliUnfolding.cxx:1416
 AliUnfolding.cxx:1417
 AliUnfolding.cxx:1418
 AliUnfolding.cxx:1419
 AliUnfolding.cxx:1420
 AliUnfolding.cxx:1421
 AliUnfolding.cxx:1422
 AliUnfolding.cxx:1423
 AliUnfolding.cxx:1424
 AliUnfolding.cxx:1425
 AliUnfolding.cxx:1426
 AliUnfolding.cxx:1427
 AliUnfolding.cxx:1428
 AliUnfolding.cxx:1429
 AliUnfolding.cxx:1430
 AliUnfolding.cxx:1431
 AliUnfolding.cxx:1432
 AliUnfolding.cxx:1433
 AliUnfolding.cxx:1434
 AliUnfolding.cxx:1435
 AliUnfolding.cxx:1436
 AliUnfolding.cxx:1437
 AliUnfolding.cxx:1438
 AliUnfolding.cxx:1439
 AliUnfolding.cxx:1440
 AliUnfolding.cxx:1441
 AliUnfolding.cxx:1442
 AliUnfolding.cxx:1443
 AliUnfolding.cxx:1444
 AliUnfolding.cxx:1445
 AliUnfolding.cxx:1446
 AliUnfolding.cxx:1447
 AliUnfolding.cxx:1448
 AliUnfolding.cxx:1449
 AliUnfolding.cxx:1450
 AliUnfolding.cxx:1451
 AliUnfolding.cxx:1452
 AliUnfolding.cxx:1453
 AliUnfolding.cxx:1454
 AliUnfolding.cxx:1455
 AliUnfolding.cxx:1456
 AliUnfolding.cxx:1457
 AliUnfolding.cxx:1458
 AliUnfolding.cxx:1459
 AliUnfolding.cxx:1460
 AliUnfolding.cxx:1461
 AliUnfolding.cxx:1462
 AliUnfolding.cxx:1463
 AliUnfolding.cxx:1464
 AliUnfolding.cxx:1465
 AliUnfolding.cxx:1466
 AliUnfolding.cxx:1467
 AliUnfolding.cxx:1468
 AliUnfolding.cxx:1469
 AliUnfolding.cxx:1470
 AliUnfolding.cxx:1471
 AliUnfolding.cxx:1472
 AliUnfolding.cxx:1473
 AliUnfolding.cxx:1474
 AliUnfolding.cxx:1475
 AliUnfolding.cxx:1476
 AliUnfolding.cxx:1477
 AliUnfolding.cxx:1478
 AliUnfolding.cxx:1479
 AliUnfolding.cxx:1480
 AliUnfolding.cxx:1481
 AliUnfolding.cxx:1482
 AliUnfolding.cxx:1483
 AliUnfolding.cxx:1484
 AliUnfolding.cxx:1485
 AliUnfolding.cxx:1486
 AliUnfolding.cxx:1487
 AliUnfolding.cxx:1488
 AliUnfolding.cxx:1489
 AliUnfolding.cxx:1490
 AliUnfolding.cxx:1491
 AliUnfolding.cxx:1492
 AliUnfolding.cxx:1493
 AliUnfolding.cxx:1494
 AliUnfolding.cxx:1495
 AliUnfolding.cxx:1496
 AliUnfolding.cxx:1497
 AliUnfolding.cxx:1498
 AliUnfolding.cxx:1499
 AliUnfolding.cxx:1500
 AliUnfolding.cxx:1501
 AliUnfolding.cxx:1502
 AliUnfolding.cxx:1503
 AliUnfolding.cxx:1504
 AliUnfolding.cxx:1505
 AliUnfolding.cxx:1506
 AliUnfolding.cxx:1507
 AliUnfolding.cxx:1508
 AliUnfolding.cxx:1509
 AliUnfolding.cxx:1510
 AliUnfolding.cxx:1511
 AliUnfolding.cxx:1512
 AliUnfolding.cxx:1513
 AliUnfolding.cxx:1514
 AliUnfolding.cxx:1515
 AliUnfolding.cxx:1516
 AliUnfolding.cxx:1517
 AliUnfolding.cxx:1518
 AliUnfolding.cxx:1519
 AliUnfolding.cxx:1520
 AliUnfolding.cxx:1521
 AliUnfolding.cxx:1522
 AliUnfolding.cxx:1523
 AliUnfolding.cxx:1524
 AliUnfolding.cxx:1525
 AliUnfolding.cxx:1526
 AliUnfolding.cxx:1527
 AliUnfolding.cxx:1528
 AliUnfolding.cxx:1529
 AliUnfolding.cxx:1530
 AliUnfolding.cxx:1531
 AliUnfolding.cxx:1532
 AliUnfolding.cxx:1533
 AliUnfolding.cxx:1534
 AliUnfolding.cxx:1535
 AliUnfolding.cxx:1536
 AliUnfolding.cxx:1537
 AliUnfolding.cxx:1538
 AliUnfolding.cxx:1539
 AliUnfolding.cxx:1540
 AliUnfolding.cxx:1541
 AliUnfolding.cxx:1542
 AliUnfolding.cxx:1543
 AliUnfolding.cxx:1544
 AliUnfolding.cxx:1545
 AliUnfolding.cxx:1546
 AliUnfolding.cxx:1547
 AliUnfolding.cxx:1548
 AliUnfolding.cxx:1549
 AliUnfolding.cxx:1550
 AliUnfolding.cxx:1551
 AliUnfolding.cxx:1552
 AliUnfolding.cxx:1553
 AliUnfolding.cxx:1554
 AliUnfolding.cxx:1555
 AliUnfolding.cxx:1556
 AliUnfolding.cxx:1557
 AliUnfolding.cxx:1558
 AliUnfolding.cxx:1559
 AliUnfolding.cxx:1560
 AliUnfolding.cxx:1561
 AliUnfolding.cxx:1562
 AliUnfolding.cxx:1563
 AliUnfolding.cxx:1564
 AliUnfolding.cxx:1565
 AliUnfolding.cxx:1566
 AliUnfolding.cxx:1567
 AliUnfolding.cxx:1568
 AliUnfolding.cxx:1569
 AliUnfolding.cxx:1570
 AliUnfolding.cxx:1571
 AliUnfolding.cxx:1572
 AliUnfolding.cxx:1573
 AliUnfolding.cxx:1574
 AliUnfolding.cxx:1575
 AliUnfolding.cxx:1576
 AliUnfolding.cxx:1577
 AliUnfolding.cxx:1578
 AliUnfolding.cxx:1579
 AliUnfolding.cxx:1580
 AliUnfolding.cxx:1581
 AliUnfolding.cxx:1582
 AliUnfolding.cxx:1583
 AliUnfolding.cxx:1584
 AliUnfolding.cxx:1585
 AliUnfolding.cxx:1586
 AliUnfolding.cxx:1587
 AliUnfolding.cxx:1588
 AliUnfolding.cxx:1589
 AliUnfolding.cxx:1590
 AliUnfolding.cxx:1591
 AliUnfolding.cxx:1592
 AliUnfolding.cxx:1593
 AliUnfolding.cxx:1594
 AliUnfolding.cxx:1595
 AliUnfolding.cxx:1596
 AliUnfolding.cxx:1597
 AliUnfolding.cxx:1598
 AliUnfolding.cxx:1599
 AliUnfolding.cxx:1600
 AliUnfolding.cxx:1601
 AliUnfolding.cxx:1602
 AliUnfolding.cxx:1603
 AliUnfolding.cxx:1604
 AliUnfolding.cxx:1605
 AliUnfolding.cxx:1606
 AliUnfolding.cxx:1607
 AliUnfolding.cxx:1608
 AliUnfolding.cxx:1609
 AliUnfolding.cxx:1610
 AliUnfolding.cxx:1611
 AliUnfolding.cxx:1612
 AliUnfolding.cxx:1613
 AliUnfolding.cxx:1614
 AliUnfolding.cxx:1615
 AliUnfolding.cxx:1616
 AliUnfolding.cxx:1617
 AliUnfolding.cxx:1618
 AliUnfolding.cxx:1619
 AliUnfolding.cxx:1620
 AliUnfolding.cxx:1621
 AliUnfolding.cxx:1622
 AliUnfolding.cxx:1623
 AliUnfolding.cxx:1624
 AliUnfolding.cxx:1625
 AliUnfolding.cxx:1626
 AliUnfolding.cxx:1627
 AliUnfolding.cxx:1628
 AliUnfolding.cxx:1629
 AliUnfolding.cxx:1630
 AliUnfolding.cxx:1631
 AliUnfolding.cxx:1632
 AliUnfolding.cxx:1633
 AliUnfolding.cxx:1634
 AliUnfolding.cxx:1635
 AliUnfolding.cxx:1636
 AliUnfolding.cxx:1637
 AliUnfolding.cxx:1638
 AliUnfolding.cxx:1639
 AliUnfolding.cxx:1640
 AliUnfolding.cxx:1641
 AliUnfolding.cxx:1642
 AliUnfolding.cxx:1643
 AliUnfolding.cxx:1644
 AliUnfolding.cxx:1645
 AliUnfolding.cxx:1646
 AliUnfolding.cxx:1647
 AliUnfolding.cxx:1648
 AliUnfolding.cxx:1649
 AliUnfolding.cxx:1650
 AliUnfolding.cxx:1651
 AliUnfolding.cxx:1652
 AliUnfolding.cxx:1653
 AliUnfolding.cxx:1654
 AliUnfolding.cxx:1655
 AliUnfolding.cxx:1656
 AliUnfolding.cxx:1657
 AliUnfolding.cxx:1658
 AliUnfolding.cxx:1659
 AliUnfolding.cxx:1660
 AliUnfolding.cxx:1661
 AliUnfolding.cxx:1662
 AliUnfolding.cxx:1663
 AliUnfolding.cxx:1664
 AliUnfolding.cxx:1665
 AliUnfolding.cxx:1666
 AliUnfolding.cxx:1667
 AliUnfolding.cxx:1668
 AliUnfolding.cxx:1669
 AliUnfolding.cxx:1670
 AliUnfolding.cxx:1671
 AliUnfolding.cxx:1672
 AliUnfolding.cxx:1673
 AliUnfolding.cxx:1674
 AliUnfolding.cxx:1675
 AliUnfolding.cxx:1676
 AliUnfolding.cxx:1677
 AliUnfolding.cxx:1678
 AliUnfolding.cxx:1679
 AliUnfolding.cxx:1680
 AliUnfolding.cxx:1681
 AliUnfolding.cxx:1682
 AliUnfolding.cxx:1683
 AliUnfolding.cxx:1684
 AliUnfolding.cxx:1685
 AliUnfolding.cxx:1686
 AliUnfolding.cxx:1687
 AliUnfolding.cxx:1688
 AliUnfolding.cxx:1689
 AliUnfolding.cxx:1690
 AliUnfolding.cxx:1691
 AliUnfolding.cxx:1692
 AliUnfolding.cxx:1693
 AliUnfolding.cxx:1694
 AliUnfolding.cxx:1695
 AliUnfolding.cxx:1696
 AliUnfolding.cxx:1697
 AliUnfolding.cxx:1698
 AliUnfolding.cxx:1699
 AliUnfolding.cxx:1700
 AliUnfolding.cxx:1701
 AliUnfolding.cxx:1702
 AliUnfolding.cxx:1703
 AliUnfolding.cxx:1704
 AliUnfolding.cxx:1705
 AliUnfolding.cxx:1706
 AliUnfolding.cxx:1707
 AliUnfolding.cxx:1708
 AliUnfolding.cxx:1709
 AliUnfolding.cxx:1710
 AliUnfolding.cxx:1711
 AliUnfolding.cxx:1712
 AliUnfolding.cxx:1713
 AliUnfolding.cxx:1714
 AliUnfolding.cxx:1715
 AliUnfolding.cxx:1716
 AliUnfolding.cxx:1717
 AliUnfolding.cxx:1718
 AliUnfolding.cxx:1719
 AliUnfolding.cxx:1720
 AliUnfolding.cxx:1721
 AliUnfolding.cxx:1722
 AliUnfolding.cxx:1723
 AliUnfolding.cxx:1724
 AliUnfolding.cxx:1725
 AliUnfolding.cxx:1726
 AliUnfolding.cxx:1727
 AliUnfolding.cxx:1728
 AliUnfolding.cxx:1729
 AliUnfolding.cxx:1730
 AliUnfolding.cxx:1731
 AliUnfolding.cxx:1732
 AliUnfolding.cxx:1733
 AliUnfolding.cxx:1734
 AliUnfolding.cxx:1735
 AliUnfolding.cxx:1736
 AliUnfolding.cxx:1737
 AliUnfolding.cxx:1738
 AliUnfolding.cxx:1739
 AliUnfolding.cxx:1740
 AliUnfolding.cxx:1741
 AliUnfolding.cxx:1742
 AliUnfolding.cxx:1743
 AliUnfolding.cxx:1744
 AliUnfolding.cxx:1745
 AliUnfolding.cxx:1746
 AliUnfolding.cxx:1747
 AliUnfolding.cxx:1748
 AliUnfolding.cxx:1749
 AliUnfolding.cxx:1750
 AliUnfolding.cxx:1751
 AliUnfolding.cxx:1752
 AliUnfolding.cxx:1753
 AliUnfolding.cxx:1754
 AliUnfolding.cxx:1755
 AliUnfolding.cxx:1756
 AliUnfolding.cxx:1757
 AliUnfolding.cxx:1758
 AliUnfolding.cxx:1759
 AliUnfolding.cxx:1760
 AliUnfolding.cxx:1761
 AliUnfolding.cxx:1762
 AliUnfolding.cxx:1763
 AliUnfolding.cxx:1764
 AliUnfolding.cxx:1765
 AliUnfolding.cxx:1766
 AliUnfolding.cxx:1767
 AliUnfolding.cxx:1768
 AliUnfolding.cxx:1769
 AliUnfolding.cxx:1770
 AliUnfolding.cxx:1771
 AliUnfolding.cxx:1772
 AliUnfolding.cxx:1773
 AliUnfolding.cxx:1774
 AliUnfolding.cxx:1775
 AliUnfolding.cxx:1776
 AliUnfolding.cxx:1777
 AliUnfolding.cxx:1778
 AliUnfolding.cxx:1779
 AliUnfolding.cxx:1780
 AliUnfolding.cxx:1781
 AliUnfolding.cxx:1782
 AliUnfolding.cxx:1783
 AliUnfolding.cxx:1784
 AliUnfolding.cxx:1785
 AliUnfolding.cxx:1786
 AliUnfolding.cxx:1787
 AliUnfolding.cxx:1788
 AliUnfolding.cxx:1789
 AliUnfolding.cxx:1790
 AliUnfolding.cxx:1791
 AliUnfolding.cxx:1792
 AliUnfolding.cxx:1793
 AliUnfolding.cxx:1794
 AliUnfolding.cxx:1795
 AliUnfolding.cxx:1796
 AliUnfolding.cxx:1797
 AliUnfolding.cxx:1798
 AliUnfolding.cxx:1799
 AliUnfolding.cxx:1800
 AliUnfolding.cxx:1801
 AliUnfolding.cxx:1802
 AliUnfolding.cxx:1803
 AliUnfolding.cxx:1804
 AliUnfolding.cxx:1805
 AliUnfolding.cxx:1806
 AliUnfolding.cxx:1807
 AliUnfolding.cxx:1808
 AliUnfolding.cxx:1809
 AliUnfolding.cxx:1810
 AliUnfolding.cxx:1811
 AliUnfolding.cxx:1812
 AliUnfolding.cxx:1813
 AliUnfolding.cxx:1814
 AliUnfolding.cxx:1815
 AliUnfolding.cxx:1816
 AliUnfolding.cxx:1817
 AliUnfolding.cxx:1818
 AliUnfolding.cxx:1819
 AliUnfolding.cxx:1820
 AliUnfolding.cxx:1821
 AliUnfolding.cxx:1822
 AliUnfolding.cxx:1823
 AliUnfolding.cxx:1824
 AliUnfolding.cxx:1825
 AliUnfolding.cxx:1826
 AliUnfolding.cxx:1827
 AliUnfolding.cxx:1828
 AliUnfolding.cxx:1829
 AliUnfolding.cxx:1830
 AliUnfolding.cxx:1831
 AliUnfolding.cxx:1832
 AliUnfolding.cxx:1833
 AliUnfolding.cxx:1834
 AliUnfolding.cxx:1835
 AliUnfolding.cxx:1836
 AliUnfolding.cxx:1837
 AliUnfolding.cxx:1838
 AliUnfolding.cxx:1839
 AliUnfolding.cxx:1840
 AliUnfolding.cxx:1841
 AliUnfolding.cxx:1842
 AliUnfolding.cxx:1843
 AliUnfolding.cxx:1844
 AliUnfolding.cxx:1845
 AliUnfolding.cxx:1846
 AliUnfolding.cxx:1847
 AliUnfolding.cxx:1848
 AliUnfolding.cxx:1849
 AliUnfolding.cxx:1850
 AliUnfolding.cxx:1851
 AliUnfolding.cxx:1852
 AliUnfolding.cxx:1853
 AliUnfolding.cxx:1854
 AliUnfolding.cxx:1855
 AliUnfolding.cxx:1856
 AliUnfolding.cxx:1857
 AliUnfolding.cxx:1858
 AliUnfolding.cxx:1859
 AliUnfolding.cxx:1860
 AliUnfolding.cxx:1861
 AliUnfolding.cxx:1862
 AliUnfolding.cxx:1863
 AliUnfolding.cxx:1864
 AliUnfolding.cxx:1865
 AliUnfolding.cxx:1866
 AliUnfolding.cxx:1867
 AliUnfolding.cxx:1868
 AliUnfolding.cxx:1869
 AliUnfolding.cxx:1870
 AliUnfolding.cxx:1871
 AliUnfolding.cxx:1872
 AliUnfolding.cxx:1873
 AliUnfolding.cxx:1874
 AliUnfolding.cxx:1875
 AliUnfolding.cxx:1876
 AliUnfolding.cxx:1877
 AliUnfolding.cxx:1878
 AliUnfolding.cxx:1879
 AliUnfolding.cxx:1880
 AliUnfolding.cxx:1881
 AliUnfolding.cxx:1882
 AliUnfolding.cxx:1883
 AliUnfolding.cxx:1884
 AliUnfolding.cxx:1885
 AliUnfolding.cxx:1886
 AliUnfolding.cxx:1887
 AliUnfolding.cxx:1888
 AliUnfolding.cxx:1889
 AliUnfolding.cxx:1890
 AliUnfolding.cxx:1891
 AliUnfolding.cxx:1892
 AliUnfolding.cxx:1893
 AliUnfolding.cxx:1894
 AliUnfolding.cxx:1895
 AliUnfolding.cxx:1896
 AliUnfolding.cxx:1897
 AliUnfolding.cxx:1898
 AliUnfolding.cxx:1899
 AliUnfolding.cxx:1900
 AliUnfolding.cxx:1901
 AliUnfolding.cxx:1902
 AliUnfolding.cxx:1903
 AliUnfolding.cxx:1904
 AliUnfolding.cxx:1905
 AliUnfolding.cxx:1906
 AliUnfolding.cxx:1907
 AliUnfolding.cxx:1908
 AliUnfolding.cxx:1909
 AliUnfolding.cxx:1910
 AliUnfolding.cxx:1911
 AliUnfolding.cxx:1912
 AliUnfolding.cxx:1913
 AliUnfolding.cxx:1914
 AliUnfolding.cxx:1915
 AliUnfolding.cxx:1916
 AliUnfolding.cxx:1917
 AliUnfolding.cxx:1918
 AliUnfolding.cxx:1919
 AliUnfolding.cxx:1920
 AliUnfolding.cxx:1921
 AliUnfolding.cxx:1922
 AliUnfolding.cxx:1923
 AliUnfolding.cxx:1924
 AliUnfolding.cxx:1925
 AliUnfolding.cxx:1926
 AliUnfolding.cxx:1927
 AliUnfolding.cxx:1928
 AliUnfolding.cxx:1929
 AliUnfolding.cxx:1930
 AliUnfolding.cxx:1931
 AliUnfolding.cxx:1932
 AliUnfolding.cxx:1933
 AliUnfolding.cxx:1934
 AliUnfolding.cxx:1935
 AliUnfolding.cxx:1936
 AliUnfolding.cxx:1937
 AliUnfolding.cxx:1938
 AliUnfolding.cxx:1939
 AliUnfolding.cxx:1940
 AliUnfolding.cxx:1941
 AliUnfolding.cxx:1942
 AliUnfolding.cxx:1943
 AliUnfolding.cxx:1944
 AliUnfolding.cxx:1945
 AliUnfolding.cxx:1946
 AliUnfolding.cxx:1947
 AliUnfolding.cxx:1948
 AliUnfolding.cxx:1949
 AliUnfolding.cxx:1950
 AliUnfolding.cxx:1951
 AliUnfolding.cxx:1952
 AliUnfolding.cxx:1953
 AliUnfolding.cxx:1954
 AliUnfolding.cxx:1955
 AliUnfolding.cxx:1956
 AliUnfolding.cxx:1957
 AliUnfolding.cxx:1958
 AliUnfolding.cxx:1959
 AliUnfolding.cxx:1960
 AliUnfolding.cxx:1961
 AliUnfolding.cxx:1962
 AliUnfolding.cxx:1963
 AliUnfolding.cxx:1964
 AliUnfolding.cxx:1965
 AliUnfolding.cxx:1966
 AliUnfolding.cxx:1967
 AliUnfolding.cxx:1968
 AliUnfolding.cxx:1969
 AliUnfolding.cxx:1970
 AliUnfolding.cxx:1971
 AliUnfolding.cxx:1972
 AliUnfolding.cxx:1973
 AliUnfolding.cxx:1974
 AliUnfolding.cxx:1975
 AliUnfolding.cxx:1976
 AliUnfolding.cxx:1977
 AliUnfolding.cxx:1978
 AliUnfolding.cxx:1979
 AliUnfolding.cxx:1980
 AliUnfolding.cxx:1981
 AliUnfolding.cxx:1982
 AliUnfolding.cxx:1983
 AliUnfolding.cxx:1984
 AliUnfolding.cxx:1985
 AliUnfolding.cxx:1986
 AliUnfolding.cxx:1987
 AliUnfolding.cxx:1988
 AliUnfolding.cxx:1989
 AliUnfolding.cxx:1990
 AliUnfolding.cxx:1991
 AliUnfolding.cxx:1992
 AliUnfolding.cxx:1993
 AliUnfolding.cxx:1994
 AliUnfolding.cxx:1995
 AliUnfolding.cxx:1996
 AliUnfolding.cxx:1997
 AliUnfolding.cxx:1998
 AliUnfolding.cxx:1999
 AliUnfolding.cxx:2000
 AliUnfolding.cxx:2001
 AliUnfolding.cxx:2002
 AliUnfolding.cxx:2003
 AliUnfolding.cxx:2004
 AliUnfolding.cxx:2005
 AliUnfolding.cxx:2006
 AliUnfolding.cxx:2007
 AliUnfolding.cxx:2008
 AliUnfolding.cxx:2009
 AliUnfolding.cxx:2010
 AliUnfolding.cxx:2011
 AliUnfolding.cxx:2012
 AliUnfolding.cxx:2013