#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;
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;
Int_t AliUnfolding::fgMaxParams = -1;
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;
Float_t AliUnfolding::fgMinuitPrecision = 1e-6;
Int_t AliUnfolding::fgMinuitMaxIterations = 1000000;
Double_t AliUnfolding::fgMinuitStrategy = 1.;
Bool_t AliUnfolding::fgMinimumInitialValue = kFALSE;
Float_t AliUnfolding::fgMinimumInitialValueFix = -1;
Bool_t AliUnfolding::fgNormalizeInput = kFALSE;
Float_t AliUnfolding::fgNotFoundEvents = 0;
Bool_t AliUnfolding::fgSkipBin0InChi2 = kFALSE;
Float_t AliUnfolding::fgBayesianSmoothing = 1;
Int_t AliUnfolding::fgBayesianIterations = 10;
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)
{
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)
{
fgOverflowBinLimit = overflowBinLimit;
Printf("AliUnfolding::SetCreateOverflowBin: overflow bin limit set to %f", overflowBinLimit);
}
void AliUnfolding::SetSkipBinsBegin(Int_t nBins)
{
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)
{
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)
{
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)
{
fgBayesianSmoothing = smoothing;
fgBayesianIterations = nIterations;
Printf("AliUnfolding::SetBayesianParameters --> Parameters set to %d iterations with smoothing %f", fgBayesianIterations, fgBayesianSmoothing);
}
void AliUnfolding::SetFunction(TF1* 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)
{
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)
{
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();
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;
for (Int_t j=1; j<=correlation->GetNbinsY(); ++j)
{
if (maxValue < correlation->GetBinContent(i, j))
{
maxValue = correlation->GetBinContent(i, j);
}
correlation->SetBinContent(i, j, correlation->GetBinContent(i, j) / sum);
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);
}
}
}
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
(*fgCorrelationCovarianceMatrix)(i, i) = (Double_t) 1e-6 / smallestError / smallestError;
if ((*fgCorrelationCovarianceMatrix)(i, i) > 1e7)
(*fgCorrelationCovarianceMatrix)(i, i) = 0;
}
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)
{
SetStaticVariables(correlation, measured, efficiency);
Int_t params = fgMaxParams;
if (fgNotFoundEvents > 0)
params++;
TVirtualFitter *minuit = TVirtualFitter::Fitter(0, params);
Double_t arglist[100];
arglist[0] = 0;
minuit->ExecuteCommand("SET PRINT", arglist, 1);
minuit->SetFCN(Chi2Function);
minuit->SetPrecision(fgMinuitPrecision);
minuit->SetMaxIterations(fgMinuitMaxIterations);
minuit->ExecuteCommand("SET STRATEGY",&fgMinuitStrategy,1);
for (Int_t i=0; i<fgMaxParams; i++)
(*fgEntropyAPriori)[i] = 1;
if (!initialConditions) {
initialConditions = measured;
} else {
Printf("AliUnfolding::UnfoldWithMinuit: Using different initial conditions...");
if (fgNormalizeInput)
initialConditions->Scale(1.0 / initialConditions->Integral());
}
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;
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;
}
arglist[0] = (float)fgMinuitMaxIterations;
Int_t status = minuit->ExecuteCommand("MIGRAD", arglist, 1);
Printf("AliUnfolding::UnfoldWithMinuit: MINUIT status is %d", status);
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];
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)
{
if (efficiency->GetBinContent(i+1) > 0)
{
value /= efficiency->GetBinContent(i+1);
error /= efficiency->GetBinContent(i+1);
}
else
{
value = 0;
error = 0;
}
}
result->SetBinContent(i+1, value);
result->SetBinError(i+1, error);
}
Int_t tmpCallCount = fgCallCount;
fgCallCount = 0;
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)
{
if (measured->Integral() <= 0)
{
Printf("AliUnfolding::UnfoldWithBayesian: ERROR: The measured spectrum is empty");
return -1;
}
const Int_t kStartBin = 0;
Int_t kMaxM = fgMaxInput;
Int_t kMaxT = fgMaxParams;
const Double_t kConvergenceLimit = kMaxT * 1e-6;
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;
}
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++) {
if (normResponse[t] != 0)
response[t][m] /= normResponse[t];
else
Printf("AliUnfolding::UnfoldWithBayesian: Empty row,column in response matrix, for truth bin %d",t);
}
}
if (initialConditions)
{
printf("Using different starting conditions...\n");
Float_t inputDistIntegral = initialConditions->Integral();
for (Int_t i=0; i<kMaxT; i++)
prior[i] = initialConditions->GetBinContent(i+1) / inputDistIntegral;
}
for (Int_t i=0; i<fgBayesianIterations || fgBayesianIterations < 0; i++)
{
if (fgDebug)
Printf("AliUnfolding::UnfoldWithBayesian: iteration %i", i);
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];
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;
}
}
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;
}
Double_t chi2LastIter = 0;
for (Int_t t=kStartBin; t<kMaxT; t++)
{
Float_t newValue = 0;
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];
newValue = (1 - fgBayesianSmoothing) * result[t] + fgBayesianSmoothing * average;
}
else
newValue = result[t];
if (prior[t] > 1e-5)
{
Double_t diff = (prior[t] - newValue) / prior[t];
chi2LastIter += diff * diff;
}
prior[t] = newValue;
}
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;
}
}
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;
}
Double_t AliUnfolding::RegularizationPol0(TVectorD& params)
{
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)
{
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);
chi2 += (der1 - der2) * (der1 - der2) / middle * fgUnfoldedAxis->GetBinWidth(i);
}
return chi2;
}
Double_t AliUnfolding::RegularizationLog(TVectorD& params)
{
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);
chi2 += (der1 - der2) * (der1 - der2);
}
return chi2;
}
Double_t AliUnfolding::RegularizationTotalCurvature(TVectorD& params)
{
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)
{
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;
if (tmp > 0 && (*fgEntropyAPriori)[i] > 0)
{
chi2 += tmp * TMath::Log(tmp / (*fgEntropyAPriori)[i]);
}
else
chi2 += 100;
}
return -chi2;
}
Double_t AliUnfolding::RegularizationRatio(TVectorD& params)
{
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 = 4 * left2 / left - middle / right - left / middle - left3 / left2 - left4 / left3;
chi2 += diff * diff;
}
return chi2;
}
Double_t AliUnfolding::RegularizationPowerLaw(TVectorD& params)
{
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.);
}
}
return chi2;
}
Double_t AliUnfolding::RegularizationLogLog(TVectorD& params)
{
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)
{
TVectorD paramsVector(fgMaxParams);
for (Int_t i=0; i<fgMaxParams; ++i)
paramsVector[i] = params[i] * params[i];
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;
TVectorD measGuessVector(fgMaxInput);
measGuessVector = (*fgCorrelationMatrix) * paramsVector;
measGuessVector -= (*fgCurrentESDVector);
#if 0
TVectorD errorGuessVector(fgMaxInput);
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
TVectorD copy(measGuessVector);
for (Int_t i=0; i<fgMaxInput; ++i)
measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
if (fgSkipBin0InChi2)
measGuessVector[0] = 0;
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);
if (i <= 3)
errorNotFoundEvents += (eff * 0.03) * (eff * 0.03) * paramsVector(i) * paramsVector(i);
}
errorNotFoundEvents += fgNotFoundEvents;
notFoundEventsConstraint = (currentNotFoundEvents - fgNotFoundEvents) * (currentNotFoundEvents - fgNotFoundEvents) / errorNotFoundEvents;
}
Float_t avg = 0;
Float_t sum = 0;
Float_t currentMult = 0;
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;
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);
}
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)
{
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);
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;
fgCallCount = 0;
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->SetBit(kCannotPick);
DrawGuess(params, presult, pres, ppen, reuseHists,unfolded);
delete [] params;
}
void AliUnfolding::RedrawInteractive() {
DrawResults(fghCorrelation,fghEfficiency,fghMeasured,fghUnfolded,fgCanvas,1,fghUnfolded);
}
void AliUnfolding::InteractiveUnfold(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions)
{
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)
{
if (pfolded == 0)
pfolded = new TCanvas;
if (ppen == 0)
ppen = new TCanvas;
if (pres == 0)
pres = new TCanvas;
TVectorD paramsVector(fgMaxParams);
for (Int_t i=0; i<fgMaxParams; ++i)
paramsVector[i] = params[i] * params[i];
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())
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");
measGuessVector -= (*fgCurrentESDVector);
TVectorD copy(measGuessVector);
for (Int_t i=0; i<fgMaxInput; ++i)
measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
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();
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)
{
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)
{
TVectorD paramsVector(fgMaxParams);
for (Int_t i=0; i<fgMaxParams; ++i)
paramsVector[i] = params[i] * params[i];
TVectorD measGuessVector(fgMaxInput);
measGuessVector = (*fgCorrelationMatrix) * paramsVector;
measGuessVector -= (*fgCurrentESDVector);
TVectorD copy(measGuessVector);
for (Int_t i=0; i<fgMaxInput; ++i)
measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
TH1* residuals = 0;
if (fgMeasuredAxis->GetXbins()->GetArray())
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());
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)
{
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)
{
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);
diff = (der1 - der2) * (der1 - der2);
}
if (fgRegularizationType == kCurvature && i > 1+fgSkipBinsBegin)
{
Double_t right = params[i];
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);
}
}
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);
}
return penalty;
}
void AliUnfolding::TF1Function(Int_t& unused1, Double_t* unused2, Double_t& chi2, Double_t *params, Int_t unused3)
{
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* , TH1* aResult)
{
if (!fgFitFunction)
{
Printf("AliUnfolding::UnfoldWithFunction: ERROR fit function not set. Exiting.");
return -1;
}
SetChi2Regularization(kNone, 0);
SetStaticVariables(correlation, measured, efficiency);
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);
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)
{
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);
}
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()));
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)
{
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));
for (Int_t m=0; m<fgMaxInput; m++)
{
if (measured->GetBinContent(m+1) < 1)
continue;
for (Int_t loop=0; loop<4; 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)));
newResult[loop]->Reset();
Unfold(correlation, efficiency, measuredClone, measuredClone, newResult[loop]);
delete measuredClone;
}
for (Int_t t=0; t<fgMaxParams; t++)
{
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];
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);
}
convoluted->Scale(measured->Integral() / convoluted->Integral());
convoluted->Add(measured, -1);
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);
}
return 0;
}