1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /* $Id: AliUnfolding.cxx 31168 2009-02-23 15:18:45Z jgrosseo $ */
18 // This class allows 1-dimensional unfolding.
19 // Methods that are implemented are chi2 minimization and bayesian unfolding.
21 // Author: Jan.Fiete.Grosse-Oetringhaus@cern.ch
23 #include "AliUnfolding.h"
26 #include <TVirtualFitter.h>
31 #include "Riostream.h"
34 using namespace std; //required for resolving the 'cout' symbol
36 TMatrixD* AliUnfolding::fgCorrelationMatrix = 0;
37 TMatrixD* AliUnfolding::fgCorrelationMatrixSquared = 0;
38 TMatrixD* AliUnfolding::fgCorrelationCovarianceMatrix = 0;
39 TVectorD* AliUnfolding::fgCurrentESDVector = 0;
40 TVectorD* AliUnfolding::fgEntropyAPriori = 0;
41 TVectorD* AliUnfolding::fgEfficiency = 0;
43 TAxis* AliUnfolding::fgUnfoldedAxis = 0;
44 TAxis* AliUnfolding::fgMeasuredAxis = 0;
46 TF1* AliUnfolding::fgFitFunction = 0;
48 AliUnfolding::MethodType AliUnfolding::fgMethodType = AliUnfolding::kInvalid;
49 Int_t AliUnfolding::fgMaxInput = -1; // bins in measured histogram
50 Int_t AliUnfolding::fgMaxParams = -1; // bins in unfolded histogram = number of fit params
51 Float_t AliUnfolding::fgOverflowBinLimit = -1;
53 AliUnfolding::RegularizationType AliUnfolding::fgRegularizationType = AliUnfolding::kPol1;
54 Float_t AliUnfolding::fgRegularizationWeight = 10000;
55 Int_t AliUnfolding::fgSkipBinsBegin = 0;
56 Float_t AliUnfolding::fgMinuitStepSize = 0.1; // (usually not needed to be changed) step size in minimization
57 Float_t AliUnfolding::fgMinuitPrecision = 1e-6; // minuit precision
58 Int_t AliUnfolding::fgMinuitMaxIterations = 1000000; // minuit maximum number of iterations
59 Double_t AliUnfolding::fgMinuitStrategy = 1.; // minuit strategy
60 Bool_t AliUnfolding::fgMinimumInitialValue = kFALSE; // set all initial values at least to the smallest value among the initial values
61 Float_t AliUnfolding::fgMinimumInitialValueFix = -1;
62 Bool_t AliUnfolding::fgNormalizeInput = kFALSE; // normalize input spectrum
63 Float_t AliUnfolding::fgNotFoundEvents = 0;
64 Bool_t AliUnfolding::fgSkipBin0InChi2 = kFALSE;
66 Float_t AliUnfolding::fgBayesianSmoothing = 1; // smoothing parameter (0 = no smoothing)
67 Int_t AliUnfolding::fgBayesianIterations = 10; // number of iterations in Bayesian method
69 Bool_t AliUnfolding::fgDebug = kFALSE;
71 Int_t AliUnfolding::fgCallCount = 0;
73 Int_t AliUnfolding::fgPowern = 5;
75 Double_t AliUnfolding::fChi2FromFit = 0.;
76 Double_t AliUnfolding::fPenaltyVal = 0.;
77 Double_t AliUnfolding::fAvgResidual = 0.;
79 Int_t AliUnfolding::fgPrintChi2Details = 0;
81 TCanvas *AliUnfolding::fgCanvas = 0;
82 TH1 *AliUnfolding::fghUnfolded = 0;
83 TH2 *AliUnfolding::fghCorrelation = 0;
84 TH1 *AliUnfolding::fghEfficiency = 0;
85 TH1 *AliUnfolding::fghMeasured = 0;
87 ClassImp(AliUnfolding)
89 //____________________________________________________________________
90 void AliUnfolding::SetUnfoldingMethod(MethodType methodType)
92 // set unfolding method
93 fgMethodType = methodType;
98 case kInvalid: name = "INVALID"; break;
99 case kChi2Minimization: name = "Chi2 Minimization"; break;
100 case kBayesian: name = "Bayesian unfolding"; break;
101 case kFunction: name = "Functional fit"; break;
103 Printf("AliUnfolding::SetUnfoldingMethod: %s enabled.", name);
106 //____________________________________________________________________
107 void AliUnfolding::SetCreateOverflowBin(Float_t overflowBinLimit)
109 // enable the creation of a overflow bin that includes all statistics below the given limit
111 fgOverflowBinLimit = overflowBinLimit;
113 Printf("AliUnfolding::SetCreateOverflowBin: overflow bin limit set to %f", overflowBinLimit);
116 //____________________________________________________________________
117 void AliUnfolding::SetSkipBinsBegin(Int_t nBins)
119 // set number of skipped bins in regularization
121 fgSkipBinsBegin = nBins;
123 Printf("AliUnfolding::SetSkipBinsBegin: skipping %d bins at the beginning of the spectrum in the regularization.", fgSkipBinsBegin);
126 //____________________________________________________________________
127 void AliUnfolding::SetNbins(Int_t nMeasured, Int_t nUnfolded)
129 // set number of bins in the input (measured) distribution and in the unfolded distribution
130 fgMaxInput = nMeasured;
131 fgMaxParams = nUnfolded;
133 if (fgCorrelationMatrix)
135 delete fgCorrelationMatrix;
136 fgCorrelationMatrix = 0;
138 if (fgCorrelationMatrixSquared)
140 fgCorrelationMatrixSquared = 0;
141 delete fgCorrelationMatrixSquared;
143 if (fgCorrelationCovarianceMatrix)
145 delete fgCorrelationCovarianceMatrix;
146 fgCorrelationCovarianceMatrix = 0;
148 if (fgCurrentESDVector)
150 delete fgCurrentESDVector;
151 fgCurrentESDVector = 0;
153 if (fgEntropyAPriori)
155 delete fgEntropyAPriori;
156 fgEntropyAPriori = 0;
165 delete fgUnfoldedAxis;
170 delete fgMeasuredAxis;
174 Printf("AliUnfolding::SetNbins: Set %d measured bins and %d unfolded bins", nMeasured, nUnfolded);
177 //____________________________________________________________________
178 void AliUnfolding::SetChi2Regularization(RegularizationType type, Float_t weight)
181 // sets the parameters for chi2 minimization
184 fgRegularizationType = type;
185 fgRegularizationWeight = weight;
187 Printf("AliUnfolding::SetChi2Regularization --> Regularization set to %d with weight %f", (Int_t) type, weight);
190 //____________________________________________________________________
191 void AliUnfolding::SetBayesianParameters(Float_t smoothing, Int_t nIterations)
194 // sets the parameters for Bayesian unfolding
197 fgBayesianSmoothing = smoothing;
198 fgBayesianIterations = nIterations;
200 Printf("AliUnfolding::SetBayesianParameters --> Parameters set to %d iterations with smoothing %f", fgBayesianIterations, fgBayesianSmoothing);
203 //____________________________________________________________________
204 void AliUnfolding::SetFunction(TF1* function)
206 // set function for unfolding with a fit function
208 fgFitFunction = function;
210 Printf("AliUnfolding::SetFunction: Set fit function with %d parameters.", function->GetNpar());
213 //____________________________________________________________________
214 Int_t AliUnfolding::Unfold(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result, Bool_t check)
216 // unfolds with unfolding method fgMethodType
219 // correlation: response matrix as measured vs. generated
220 // 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.
221 // measured: the measured spectrum
222 // initialConditions: (optional) initial conditions for the unfolding. if 0 the measured spectrum is used as initial conditions.
223 // result: target for the unfolded result
224 // check: depends on the unfolding method, see comments in specific functions
226 // return code: see UnfoldWithMinuit/UnfoldWithBayesian/UnfoldWithFunction
228 if (fgMaxInput == -1)
230 Printf("AliUnfolding::Unfold: WARNING. Number of measured bins not set with SetNbins. Using number of bins in measured distribution");
231 fgMaxInput = measured->GetNbinsX();
233 if (fgMaxParams == -1)
235 Printf("AliUnfolding::Unfold: WARNING. Number of unfolded bins not set with SetNbins. Using number of bins in measured distribution");
236 fgMaxParams = measured->GetNbinsX();
239 if (fgOverflowBinLimit > 0)
240 CreateOverflowBin(correlation, measured);
242 switch (fgMethodType)
246 Printf("AliUnfolding::Unfold: ERROR: Unfolding method not set. Use SetUnfoldingMethod. Exiting...");
249 case kChi2Minimization:
250 return UnfoldWithMinuit(correlation, efficiency, measured, initialConditions, result, check);
252 return UnfoldWithBayesian(correlation, efficiency, measured, initialConditions, result);
254 return UnfoldWithFunction(correlation, efficiency, measured, initialConditions, result);
262 //____________________________________________________________________
263 void AliUnfolding::SetStaticVariables(TH2* correlation, TH1* measured, TH1* efficiency)
265 // fill static variables needed for minuit fit
267 if (!fgCorrelationMatrix)
268 fgCorrelationMatrix = new TMatrixD(fgMaxInput, fgMaxParams);
269 if (!fgCorrelationMatrixSquared)
270 fgCorrelationMatrixSquared = new TMatrixD(fgMaxInput, fgMaxParams);
271 if (!fgCorrelationCovarianceMatrix)
272 fgCorrelationCovarianceMatrix = new TMatrixD(fgMaxInput, fgMaxInput);
273 if (!fgCurrentESDVector)
274 fgCurrentESDVector = new TVectorD(fgMaxInput);
275 if (!fgEntropyAPriori)
276 fgEntropyAPriori = new TVectorD(fgMaxParams);
278 fgEfficiency = new TVectorD(fgMaxParams);
281 delete fgUnfoldedAxis;
285 fgUnfoldedAxis = new TAxis(*(correlation->GetXaxis()));
288 delete fgMeasuredAxis;
292 fgMeasuredAxis = new TAxis(*(correlation->GetYaxis()));
294 fgCorrelationMatrix->Zero();
295 fgCorrelationCovarianceMatrix->Zero();
296 fgCurrentESDVector->Zero();
297 fgEntropyAPriori->Zero();
299 // normalize correction for given nPart
300 for (Int_t i=1; i<=correlation->GetNbinsX(); ++i)
302 Double_t sum = correlation->Integral(i, i, 1, correlation->GetNbinsY());
305 Float_t maxValue = 0;
307 for (Int_t j=1; j<=correlation->GetNbinsY(); ++j)
309 // find most probably value
310 if (maxValue < correlation->GetBinContent(i, j))
312 maxValue = correlation->GetBinContent(i, j);
317 correlation->SetBinContent(i, j, correlation->GetBinContent(i, j) / sum);// * correlation->GetXaxis()->GetBinWidth(i));
318 correlation->SetBinError(i, j, correlation->GetBinError(i, j) / sum);
320 if (i <= fgMaxParams && j <= fgMaxInput)
322 (*fgCorrelationMatrix)(j-1, i-1) = correlation->GetBinContent(i, j);
323 (*fgCorrelationMatrixSquared)(j-1, i-1) = correlation->GetBinContent(i, j) * correlation->GetBinContent(i, j);
327 //printf("MPV for Ntrue = %f is %f\n", fCurrentCorrelation->GetXaxis()->GetBinCenter(i), fCurrentCorrelation->GetYaxis()->GetBinCenter(maxBin));
331 Float_t smallestError = 1;
332 if (fgNormalizeInput)
334 Float_t sumMeasured = measured->Integral();
335 measured->Scale(1.0 / sumMeasured);
336 smallestError /= sumMeasured;
339 for (Int_t i=0; i<fgMaxInput; ++i)
341 (*fgCurrentESDVector)[i] = measured->GetBinContent(i+1);
342 if (measured->GetBinError(i+1) > 0)
344 (*fgCorrelationCovarianceMatrix)(i, i) = (Double_t) 1e-6 / measured->GetBinError(i+1) / measured->GetBinError(i+1);
346 else // in this case put error of 1, otherwise 0 bins are not added to the chi2...
347 (*fgCorrelationCovarianceMatrix)(i, i) = (Double_t) 1e-6 / smallestError / smallestError;
349 if ((*fgCorrelationCovarianceMatrix)(i, i) > 1e7)
350 (*fgCorrelationCovarianceMatrix)(i, i) = 0;
351 //Printf("%d, %e", i, (*fgCorrelationCovarianceMatrix)(i, i));
354 // efficiency is expected to match bin width of result
355 for (Int_t i=0; i<fgMaxParams; ++i)
357 (*fgEfficiency)(i) = efficiency->GetBinContent(i+1);
360 if (correlation->GetNbinsX() != fgMaxParams || correlation->GetNbinsY() != fgMaxInput)
361 cout << "Response histo has incorrect dimensions; expect (" << fgMaxParams << ", " << fgMaxInput << "), got (" << correlation->GetNbinsX() << ", " << correlation->GetNbinsY() << ")" << endl;
365 //____________________________________________________________________
366 Int_t AliUnfolding::UnfoldWithMinuit(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result, Bool_t check)
369 // implementation of unfolding (internal function)
371 // unfolds <measured> using response from <correlation> and effiency <efficiency>
372 // output is in <result>
373 // <initialConditions> set the initial values for the minimization, if 0 <measured> is used
374 // negative values in initialConditions mean that the given parameter is fixed to the absolute of the value
375 // if <check> is true no unfolding is made, instead only the chi2 without unfolding is printed
377 // returns minuit status (0 = success), (-1 when check was set)
380 SetStaticVariables(correlation, measured, efficiency);
382 // Initialize TMinuit via generic fitter interface
383 Int_t params = fgMaxParams;
384 if (fgNotFoundEvents > 0)
387 TVirtualFitter *minuit = TVirtualFitter::Fitter(0, params);
388 Double_t arglist[100];
389 // minuit->SetDefaultFitter("Minuit2");
391 // disable any output (-1), unfortuantly we do not see warnings anymore then. Have to find another way...
393 minuit->ExecuteCommand("SET PRINT", arglist, 1);
395 // however, enable warnings
396 //minuit->ExecuteCommand("SET WAR", arglist, 0);
398 // set minimization function
399 minuit->SetFCN(Chi2Function);
402 minuit->SetPrecision(fgMinuitPrecision);
404 minuit->SetMaxIterations(fgMinuitMaxIterations);
406 minuit->ExecuteCommand("SET STRATEGY",&fgMinuitStrategy,1);
408 for (Int_t i=0; i<fgMaxParams; i++)
409 (*fgEntropyAPriori)[i] = 1;
411 // set initial conditions as a-priori distribution for MRX regularization
413 for (Int_t i=0; i<fgMaxParams; i++)
414 if (initialConditions && initialConditions->GetBinContent(i+1) > 0)
415 (*fgEntropyAPriori)[i] = initialConditions->GetBinContent(i+1);
418 if (!initialConditions) {
419 initialConditions = measured;
421 Printf("AliUnfolding::UnfoldWithMinuit: Using different initial conditions...");
422 //new TCanvas; initialConditions->DrawCopy();
423 if (fgNormalizeInput)
424 initialConditions->Scale(1.0 / initialConditions->Integral());
427 // extract minimum value from initial conditions (if we set a value to 0 it will stay 0)
428 Float_t minValue = 1e35;
429 if (fgMinimumInitialValueFix < 0)
431 for (Int_t i=0; i<fgMaxParams; ++i)
433 Int_t bin = initialConditions->GetXaxis()->FindBin(result->GetXaxis()->GetBinCenter(i+1));
434 if (initialConditions->GetBinContent(bin) > 0)
435 minValue = TMath::Min(minValue, (Float_t) initialConditions->GetBinContent(bin));
439 minValue = fgMinimumInitialValueFix;
441 Double_t* results = new Double_t[fgMaxParams+1];
442 for (Int_t i=0; i<fgMaxParams; ++i)
444 Int_t bin = initialConditions->GetXaxis()->FindBin(result->GetXaxis()->GetBinCenter(i+1));
445 results[i] = initialConditions->GetBinContent(bin);
451 results[i] = -results[i];
454 if (!fix && fgMinimumInitialValue && results[i] < minValue)
455 results[i] = minValue;
457 // minuit sees squared values to prevent it from going negative...
458 results[i] = TMath::Sqrt(results[i]);
460 minuit->SetParameter(i, Form("param%d", i), results[i], (fix) ? 0 : fgMinuitStepSize, 0, 0);
462 if (fgNotFoundEvents > 0)
464 results[fgMaxParams] = efficiency->GetBinContent(1);
465 minuit->SetParameter(fgMaxParams, "vtx0", results[fgMaxParams], fgMinuitStepSize / 100, 0.01, 0.80);
470 Chi2Function(dummy, 0, chi2, results, 0);
471 printf("AliUnfolding::UnfoldWithMinuit: Chi2 of initial parameters is = %f\n", chi2);
480 // first param is number of iterations, second is precision....
481 arglist[0] = (float)fgMinuitMaxIterations;
482 // arglist[1] = 1e-5;
483 // minuit->ExecuteCommand("SET PRINT", arglist, 3);
484 // minuit->ExecuteCommand("SCAN", arglist, 0);
485 Int_t status = minuit->ExecuteCommand("MIGRAD", arglist, 1);
486 Printf("AliUnfolding::UnfoldWithMinuit: MINUIT status is %d", status);
487 //printf("!!!!!!!!!!!!!! MIGRAD finished: Starting MINOS !!!!!!!!!!!!!!");
488 //minuit->ExecuteCommand("MINOS", arglist, 0);
490 if (fgNotFoundEvents > 0)
492 results[fgMaxParams] = minuit->GetParameter(fgMaxParams);
493 Printf("Efficiency for bin 0 changed from %f to %f", efficiency->GetBinContent(1), results[fgMaxParams]);
494 efficiency->SetBinContent(1, results[fgMaxParams]);
497 for (Int_t i=0; i<fgMaxParams; ++i)
499 results[i] = minuit->GetParameter(i);
500 Double_t value = results[i] * results[i];
501 // error is : 2 * (relError on results[i]) * (value) = 2 * (minuit->GetParError(i) / minuit->GetParameter(i)) * (minuit->GetParameter(i) * minuit->GetParameter(i))
503 if (TMath::IsNaN(minuit->GetParError(i)))
504 Printf("WARNING: Parameter %d error is nan", i);
506 error = 2 * minuit->GetParError(i) * results[i];
510 //printf("value before efficiency correction: %f\n",value);
511 if (efficiency->GetBinContent(i+1) > 0)
513 value /= efficiency->GetBinContent(i+1);
514 error /= efficiency->GetBinContent(i+1);
522 //printf("value after efficiency correction: %f +/- %f\n",value,error);
523 result->SetBinContent(i+1, value);
524 result->SetBinError(i+1, error);
527 Int_t tmpCallCount = fgCallCount;
528 fgCallCount = 0; // needs to be 0 so that the Chi2Function prints its output
529 Chi2Function(dummy, 0, chi2, results, 0);
531 Printf("AliUnfolding::UnfoldWithMinuit: iterations %d. Chi2 of final parameters is = %f", tmpCallCount, chi2);
538 //____________________________________________________________________
539 Int_t AliUnfolding::UnfoldWithBayesian(TH2* correlation, TH1* aEfficiency, TH1* measured, TH1* initialConditions, TH1* aResult)
542 // unfolds a spectrum using the Bayesian method
545 if (measured->Integral() <= 0)
547 Printf("AliUnfolding::UnfoldWithBayesian: ERROR: The measured spectrum is empty");
551 const Int_t kStartBin = 0;
553 Int_t kMaxM = fgMaxInput; //<= fCurrentCorrelation->GetNbinsY(); // max measured axis
554 Int_t kMaxT = fgMaxParams; //<= fCurrentCorrelation->GetNbinsX(); // max true axis
556 // convergence limit: kMaxT * 0.001^2 = kMaxT * 1e-6 (e.g. 250 bins --> 2.5 e-4)
557 const Double_t kConvergenceLimit = kMaxT * 1e-6;
559 // store information in arrays, to increase processing speed (~ factor 5)
560 Double_t* measuredCopy = new Double_t[kMaxM];
561 Double_t* measuredError = new Double_t[kMaxM];
562 Double_t* prior = new Double_t[kMaxT];
563 Double_t* result = new Double_t[kMaxT];
564 Double_t* efficiency = new Double_t[kMaxT];
565 Double_t* binWidths = new Double_t[kMaxT];
566 Double_t* normResponse = new Double_t[kMaxT];
568 Double_t** response = new Double_t*[kMaxT];
569 Double_t** inverseResponse = new Double_t*[kMaxT];
570 for (Int_t i=0; i<kMaxT; i++)
572 response[i] = new Double_t[kMaxM];
573 inverseResponse[i] = new Double_t[kMaxM];
578 Float_t measuredIntegral = measured->Integral();
579 for (Int_t m=0; m<kMaxM; m++)
581 measuredCopy[m] = measured->GetBinContent(m+1) / measuredIntegral;
582 measuredError[m] = measured->GetBinError(m+1) / measuredIntegral;
584 for (Int_t t=0; t<kMaxT; t++)
586 response[t][m] = correlation->GetBinContent(t+1, m+1);
587 inverseResponse[t][m] = 0;
588 normResponse[t] += correlation->GetBinContent(t+1, m+1);
592 for (Int_t t=0; t<kMaxT; t++)
596 efficiency[t] = aEfficiency->GetBinContent(t+1);
601 prior[t] = measuredCopy[t];
603 binWidths[t] = aResult->GetXaxis()->GetBinWidth(t+1);
605 for (Int_t m=0; m<kMaxM; m++) { // Normalise response matrix
606 if (normResponse[t] != 0)
607 response[t][m] /= normResponse[t];
609 Printf("AliUnfolding::UnfoldWithBayesian: Empty row,column in response matrix, for truth bin %d",t);
613 // pick prior distribution
614 if (initialConditions)
616 printf("Using different starting conditions...\n");
618 Float_t inputDistIntegral = initialConditions->Integral();
619 for (Int_t i=0; i<kMaxT; i++)
620 prior[i] = initialConditions->GetBinContent(i+1) / inputDistIntegral;
623 //TH1F* convergence = new TH1F("convergence", "convergence", 200, 0.5, 200.5);
627 for (Int_t i=0; i<fgBayesianIterations || fgBayesianIterations < 0; i++)
630 Printf("AliUnfolding::UnfoldWithBayesian: iteration %i", i);
632 // calculate IR from Bayes theorem
633 // IR_ji = R_ij * prior_i / sum_k(R_kj * prior_k)
635 Double_t chi2Measured = 0;
636 for (Int_t m=0; m<kMaxM; m++)
639 for (Int_t t = kStartBin; t<kMaxT; t++)
640 norm += response[t][m] * prior[t] * efficiency[t];
642 // calc. chi2: (measured - response * prior) / error
643 if (measuredError[m] > 0)
645 Double_t value = (measuredCopy[m] - norm) / measuredError[m];
646 chi2Measured += value * value;
651 for (Int_t t = kStartBin; t<kMaxT; t++)
652 inverseResponse[t][m] = response[t][m] * prior[t] * efficiency[t] / norm;
656 for (Int_t t = kStartBin; t<kMaxT; t++)
657 inverseResponse[t][m] = 0;
660 //Printf("chi2Measured of the last prior is %e", chi2Measured);
662 for (Int_t t = kStartBin; t<kMaxT; t++)
665 for (Int_t m=0; m<kMaxM; m++)
666 value += inverseResponse[t][m] * measuredCopy[m];
668 if (efficiency[t] > 0)
669 result[t] = value / efficiency[t];
675 // draw intermediate result
676 for (Int_t t=0; t<kMaxT; t++)
678 aResult->SetBinContent(t+1, result[t]);
680 aResult->SetMarkerStyle(24+i);
681 aResult->SetMarkerColor(2);
682 aResult->DrawCopy((i == 0) ? "P" : "PSAME");
685 Double_t chi2LastIter = 0;
686 // regularization (simple smoothing)
687 for (Int_t t=kStartBin; t<kMaxT; t++)
689 Float_t newValue = 0;
691 // 0 bin excluded from smoothing
692 if (t > kStartBin+2 && t<kMaxT-1)
694 Float_t average = (result[t-1] / binWidths[t-1] + result[t] / binWidths[t] + result[t+1] / binWidths[t+1]) / 3 * binWidths[t];
696 // weight the average with the regularization parameter
697 newValue = (1 - fgBayesianSmoothing) * result[t] + fgBayesianSmoothing * average;
700 newValue = result[t];
702 // calculate chi2 (change from last iteration)
705 Double_t diff = (prior[t] - newValue) / prior[t];
706 chi2LastIter += diff * diff;
711 //printf("Chi2 of %d iteration = %e\n", i, chi2LastIter);
712 //convergence->Fill(i+1, chi2LastIter);
714 if (fgBayesianIterations < 0 && chi2LastIter < kConvergenceLimit)
716 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);
719 } // end of iterations
721 //new TCanvas; convergence->DrawCopy(); gPad->SetLogy();
722 //delete convergence;
725 if (!fgNormalizeInput)
726 factor = measuredIntegral;
727 for (Int_t t=0; t<kMaxT; t++)
728 aResult->SetBinContent(t+1, result[t] * factor);
730 delete[] measuredCopy;
731 delete[] measuredError;
736 delete[] normResponse;
738 for (Int_t i=0; i<kMaxT; i++)
740 delete[] response[i];
741 delete[] inverseResponse[i];
744 delete[] inverseResponse;
749 // Calculate the covariance matrix, all arguments are taken from NIM,A362,487-498,1995
751 /*printf("Calculating covariance matrix. This may take some time...\n");
753 // check if this is the right one...
754 TH1* sumHist = GetMultiplicityMC(inputRange, eventType)->ProjectionY("sumHist", 1, GetMultiplicityMC(inputRange, eventType)->GetNbinsX());
756 Int_t xBins = hInverseResponseBayes->GetNbinsX();
757 Int_t yBins = hInverseResponseBayes->GetNbinsY();
759 // calculate "unfolding matrix" Mij
760 Float_t matrixM[251][251];
761 for (Int_t i=1; i<=xBins; i++)
763 for (Int_t j=1; j<=yBins; j++)
765 if (fCurrentEfficiency->GetBinContent(i) > 0)
766 matrixM[i-1][j-1] = hInverseResponseBayes->GetBinContent(i, j) / fCurrentEfficiency->GetBinContent(i);
768 matrixM[i-1][j-1] = 0;
772 Float_t* vectorn = new Float_t[yBins];
773 for (Int_t j=1; j<=yBins; j++)
774 vectorn[j-1] = fCurrentESD->GetBinContent(j);
776 // first part of covariance matrix, depends on input distribution n(E)
777 Float_t cov1[251][251];
779 Float_t nEvents = fCurrentESD->Integral(); // N
784 for (Int_t k=0; k<xBins; k++)
786 printf("In Cov1: %d\n", k);
787 for (Int_t l=0; l<yBins; l++)
791 // sum_j Mkj Mlj n(Ej) * (1 - n(Ej) / N)
792 for (Int_t j=0; j<yBins; j++)
793 cov1[k][l] += matrixM[k][j] * matrixM[l][j] * vectorn[j]
794 * (1.0 - vectorn[j] / nEvents);
796 // - sum_i,j (i != j) Mki Mlj n(Ei) n(Ej) / N
797 for (Int_t i=0; i<yBins; i++)
798 for (Int_t j=0; j<yBins; j++)
802 cov1[k][l] -= matrixM[k][i] * matrixM[l][j] * vectorn[i]
803 * vectorn[j] / nEvents;
808 printf("Cov1 finished\n");
810 TH2F* cov = (TH2F*) hInverseResponseBayes->Clone("cov");
813 for (Int_t i=1; i<=xBins; i++)
814 for (Int_t j=1; j<=yBins; j++)
815 cov->SetBinContent(i, j, cov1[i-1][j-1]);
820 // second part of covariance matrix, depends on response matrix
821 Float_t cov2[251][251];
823 // Cov[P(Er|Cu), P(Es|Cu)] term
824 Float_t covTerm[100][100][100];
825 for (Int_t r=0; r<yBins; r++)
826 for (Int_t u=0; u<xBins; u++)
827 for (Int_t s=0; s<yBins; s++)
830 covTerm[r][u][s] = 1.0 / sumHist->GetBinContent(u+1) * hResponse->GetBinContent(u+1, r+1)
831 * (1.0 - hResponse->GetBinContent(u+1, r+1));
833 covTerm[r][u][s] = - 1.0 / sumHist->GetBinContent(u+1) * hResponse->GetBinContent(u+1, r+1)
834 * hResponse->GetBinContent(u+1, s+1);
837 for (Int_t k=0; k<xBins; k++)
838 for (Int_t l=0; l<yBins; l++)
841 printf("In Cov2: %d %d\n", k, l);
842 for (Int_t i=0; i<yBins; i++)
843 for (Int_t j=0; j<yBins; j++)
845 //printf("In Cov2: %d %d %d %d\n", k, l, i, j);
846 // calculate Cov(Mki, Mlj) = sum{ru},{su} ...
848 for (Int_t r=0; r<yBins; r++)
849 for (Int_t u=0; u<xBins; u++)
850 for (Int_t s=0; s<yBins; s++)
852 if (hResponse->GetBinContent(u+1, r+1) == 0 || hResponse->GetBinContent(u+1, s+1) == 0
853 || hResponse->GetBinContent(u+1, i+1) == 0)
856 tmpCov += BayesCovarianceDerivate(matrixM, hResponse, fCurrentEfficiency, k, i, r, u)
857 * BayesCovarianceDerivate(matrixM, hResponse, fCurrentEfficiency, l, j, s, u)
861 cov2[k][l] += fCurrentESD->GetBinContent(i+1) * fCurrentESD->GetBinContent(j+1) * tmpCov;
865 printf("Cov2 finished\n");
867 for (Int_t i=1; i<=xBins; i++)
868 for (Int_t j=1; j<=yBins; j++)
869 cov->SetBinContent(i, j, cov1[i-1][j-1] + cov2[i-1][j-1]);
875 //____________________________________________________________________
876 Double_t AliUnfolding::RegularizationPol0(TVectorD& params)
878 // homogenity term for minuit fitting
879 // pure function of the parameters
880 // prefers constant function (pol0)
882 // Does not take into account efficiency
885 for (Int_t i=1+fgSkipBinsBegin; i<fgMaxParams; ++i)
887 Double_t right = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
888 Double_t left = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
892 Double_t diff = (right - left);
893 chi2 += diff * diff / left / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
900 //____________________________________________________________________
901 Double_t AliUnfolding::RegularizationPol1(TVectorD& params)
903 // homogenity term for minuit fitting
904 // pure function of the parameters
905 // prefers linear function (pol1)
907 // Does not take into account efficiency
910 for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
912 if (params[i-1] == 0)
915 Double_t right = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
916 Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
917 Double_t left = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);
919 Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
920 Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
922 //Double_t diff = (der1 - der2) / middle;
923 //chi2 += diff * diff;
924 chi2 += (der1 - der2) * (der1 - der2) / middle * fgUnfoldedAxis->GetBinWidth(i);
930 //____________________________________________________________________
931 Double_t AliUnfolding::RegularizationLog(TVectorD& params)
933 // homogenity term for minuit fitting
934 // pure function of the parameters
935 // prefers logarithmic function (log)
937 // Does not take into account efficiency
941 for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
943 if (params[i-1] == 0 || params[i] == 0 || params[i-2] == 0)
946 Double_t right = log(params[i] / fgUnfoldedAxis->GetBinWidth(i+1));
947 Double_t middle = log(params[i-1] / fgUnfoldedAxis->GetBinWidth(i));
948 Double_t left = log(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1));
950 Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
951 Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
953 //Double_t error = 1. / params[i] + 4. / params[i-1] + 1. / params[i-2];
955 //if (fgCallCount == 0)
956 // Printf("%d %f %f", i, (der1 - der2) * (der1 - der2), error);
957 chi2 += (der1 - der2) * (der1 - der2);// / error;
963 //____________________________________________________________________
964 Double_t AliUnfolding::RegularizationTotalCurvature(TVectorD& params)
966 // homogenity term for minuit fitting
967 // pure function of the parameters
968 // minimizes the total curvature (from Unfolding Methods In High-Energy Physics Experiments,
969 // V. Blobel (Hamburg U.) . DESY 84/118, Dec 1984. 40pp.
971 // Does not take into account efficiency
975 for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
977 Double_t right = params[i];
978 Double_t middle = params[i-1];
979 Double_t left = params[i-2];
981 Double_t der1 = (right - middle);
982 Double_t der2 = (middle - left);
984 Double_t diff = (der1 - der2);
992 //____________________________________________________________________
993 Double_t AliUnfolding::RegularizationEntropy(TVectorD& params)
995 // homogenity term for minuit fitting
996 // pure function of the parameters
997 // calculates entropy, from
998 // The method of reduced cross-entropy (M. Schmelling 1993)
1000 // Does not take into account efficiency
1002 Double_t paramSum = 0;
1004 for (Int_t i=fgSkipBinsBegin; i<fgMaxParams; ++i)
1005 paramSum += params[i];
1008 for (Int_t i=fgSkipBinsBegin; i<fgMaxParams; ++i)
1010 Double_t tmp = params[i] / paramSum;
1011 //Double_t tmp = params[i];
1012 if (tmp > 0 && (*fgEntropyAPriori)[i] > 0)
1014 chi2 += tmp * TMath::Log(tmp / (*fgEntropyAPriori)[i]);
1023 //____________________________________________________________________
1024 Double_t AliUnfolding::RegularizationRatio(TVectorD& params)
1026 // homogenity term for minuit fitting
1027 // pure function of the parameters
1029 // Does not take into account efficiency
1033 for (Int_t i=5+fgSkipBinsBegin; i<fgMaxParams; ++i)
1035 if (params[i-1] == 0 || params[i] == 0)
1038 Double_t right = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
1039 Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
1040 Double_t left = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);
1041 Double_t left2 = params[i-3] / fgUnfoldedAxis->GetBinWidth(i-2);
1042 Double_t left3 = params[i-4] / fgUnfoldedAxis->GetBinWidth(i-3);
1043 Double_t left4 = params[i-5] / fgUnfoldedAxis->GetBinWidth(i-4);
1045 //Double_t diff = left / middle - middle / right;
1046 //Double_t diff = 2 * left / middle - middle / right - left2 / left;
1047 Double_t diff = 4 * left2 / left - middle / right - left / middle - left3 / left2 - left4 / left3;
1049 chi2 += diff * diff;// / middle;
1055 //____________________________________________________________________
1056 Double_t AliUnfolding::RegularizationPowerLaw(TVectorD& params)
1058 // homogenity term for minuit fitting
1059 // pure function of the parameters
1060 // prefers power law with n = -5
1062 // Does not take into account efficiency
1066 Double_t right = 0.;
1067 Double_t middle = 0.;
1070 for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
1072 if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
1075 if (fgUnfoldedAxis->GetBinWidth(i+1) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i-1) < 1e-8)
1078 middle = TMath::Power(params[i-1] / fgUnfoldedAxis->GetBinWidth(i),fgPowern);
1081 right = TMath::Power(params[i] / fgUnfoldedAxis->GetBinWidth(i+1),fgPowern)/middle;
1083 left = TMath::Power(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1),fgPowern)/middle;
1087 Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
1088 Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
1090 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;
1091 // printf("i: %d chi2 = %f\n",i,chi2);
1099 //____________________________________________________________________
1100 Double_t AliUnfolding::RegularizationLogLog(TVectorD& params)
1102 // homogenity term for minuit fitting
1103 // pure function of the parameters
1104 // prefers a powerlaw (linear on a log-log scale)
1106 // The calculation takes into account the efficiencies
1110 for (Int_t i=2+fgSkipBinsBegin; i<fgMaxParams; ++i)
1112 if (params[i-1] == 0 || params[i] == 0 || params[i-2] == 0)
1114 if ((*fgEfficiency)(i-1) == 0 || (*fgEfficiency)(i) == 0 || (*fgEfficiency)(i-2) == 0)
1117 Double_t right = log(params[i] / (*fgEfficiency)(i) / fgUnfoldedAxis->GetBinWidth(i+1));
1118 Double_t middle = log(params[i-1] / (*fgEfficiency)(i-1) / fgUnfoldedAxis->GetBinWidth(i));
1119 Double_t left = log(params[i-2] / (*fgEfficiency)(i-2) / fgUnfoldedAxis->GetBinWidth(i-1));
1121 Double_t der1 = (right - middle) / ( log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i)) );
1122 Double_t der2 = (middle - left) /( log(fgUnfoldedAxis->GetBinCenter(i)) - log(fgUnfoldedAxis->GetBinCenter(i-1)) );
1124 double tmp = (log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i-1)))/2.;
1125 Double_t dder = (der1-der2) / tmp;
1127 chi2 += dder * dder;
1135 //____________________________________________________________________
1136 void AliUnfolding::Chi2Function(Int_t&, Double_t*, Double_t& chi2, Double_t *params, Int_t)
1139 // fit function for minuit
1140 // does: (m - Ad)W(m - Ad) where m = measured, A correlation matrix, d = guess, W = covariance matrix
1143 // TODO use static members for the variables here to speed up processing (no construction/deconstruction)
1146 TVectorD paramsVector(fgMaxParams);
1147 for (Int_t i=0; i<fgMaxParams; ++i)
1148 paramsVector[i] = params[i] * params[i];
1150 // calculate penalty factor
1151 Double_t penaltyVal = 0;
1153 switch (fgRegularizationType)
1156 case kPol0: penaltyVal = RegularizationPol0(paramsVector); break;
1157 case kPol1: penaltyVal = RegularizationPol1(paramsVector); break;
1158 case kCurvature: penaltyVal = RegularizationTotalCurvature(paramsVector); break;
1159 case kEntropy: penaltyVal = RegularizationEntropy(paramsVector); break;
1160 case kLog: penaltyVal = RegularizationLog(paramsVector); break;
1161 case kRatio: penaltyVal = RegularizationRatio(paramsVector); break;
1162 case kPowerLaw: penaltyVal = RegularizationPowerLaw(paramsVector); break;
1163 case kLogLog: penaltyVal = RegularizationLogLog(paramsVector); break;
1166 penaltyVal *= fgRegularizationWeight; //beta*PU
1169 TVectorD measGuessVector(fgMaxInput);
1170 measGuessVector = (*fgCorrelationMatrix) * paramsVector;
1173 measGuessVector -= (*fgCurrentESDVector);
1176 // new error calcuation using error on the guess
1179 TVectorD errorGuessVector(fgMaxInput);
1180 //errorGuessVector = (*fgCorrelationMatrixSquared) * paramsVector;
1181 errorGuessVector = (*fgCorrelationMatrix) * paramsVector;
1183 Double_t chi2FromFit = 0;
1184 for (Int_t i=0; i<fgMaxInput; ++i)
1187 if (errorGuessVector(i) > 0)
1188 error = errorGuessVector(i);
1189 chi2FromFit += measGuessVector(i) * measGuessVector(i) / error;
1193 // old error calcuation using the error on the measured
1194 TVectorD copy(measGuessVector);
1197 // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
1198 // normal way is like this:
1199 // measGuessVector *= (*fgCorrelationCovarianceMatrix);
1200 // optimized way like this:
1201 for (Int_t i=0; i<fgMaxInput; ++i)
1202 measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
1205 if (fgSkipBin0InChi2)
1206 measGuessVector[0] = 0;
1208 // (Ad - m) W (Ad - m)
1209 // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
1210 // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see UnfoldWithMinuit)
1211 Double_t chi2FromFit = measGuessVector * copy * 1e6;
1214 Double_t notFoundEventsConstraint = 0;
1215 Double_t currentNotFoundEvents = 0;
1216 Double_t errorNotFoundEvents = 0;
1218 if (fgNotFoundEvents > 0)
1220 for (Int_t i=0; i<fgMaxParams; ++i)
1222 Double_t eff = (1.0 / (*fgEfficiency)(i) - 1);
1224 eff = (1.0 / params[fgMaxParams] - 1);
1225 currentNotFoundEvents += eff * paramsVector(i);
1226 errorNotFoundEvents += eff * eff * paramsVector(i); // error due to guess (paramsVector)
1228 errorNotFoundEvents += (eff * 0.03) * (eff * 0.03) * paramsVector(i) * paramsVector(i); // error on eff
1229 // if ((fgCallCount % 10000) == 0)
1230 //Printf("%d %f %f %f", i, (*fgEfficiency)(i), paramsVector(i), currentNotFoundEvents);
1232 errorNotFoundEvents += fgNotFoundEvents;
1233 // TODO add error on background, background estimate
1235 notFoundEventsConstraint = (currentNotFoundEvents - fgNotFoundEvents) * (currentNotFoundEvents - fgNotFoundEvents) / errorNotFoundEvents;
1240 Float_t currentMult = 0;
1241 // try to match dn/deta
1242 for (Int_t i=0; i<fgMaxParams; ++i)
1244 avg += paramsVector(i) * currentMult;
1245 sum += paramsVector(i);
1246 currentMult += fgUnfoldedAxis->GetBinWidth(i);
1249 Float_t chi2Avg = 0; //(avg - 3.73) * (avg - 3.73) * 100;
1251 chi2 = chi2FromFit + penaltyVal + notFoundEventsConstraint + chi2Avg;
1253 if ((fgCallCount++ % 1000) == 0)
1256 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);
1258 //for (Int_t i=0; i<fgMaxInput; ++i)
1259 // Printf("%d: %f", i, measGuessVector(i) * copy(i) * 1e6);
1262 fChi2FromFit = chi2FromFit;
1263 fPenaltyVal = penaltyVal;
1266 //____________________________________________________________________
1267 void AliUnfolding::MakePads() {
1268 TPad *presult = new TPad("presult","result",0,0.4,1,1);
1269 presult->SetNumber(1);
1272 TPad *pres = new TPad("pres","residuals",0,0.2,1,0.4);
1275 TPad *ppen = new TPad("ppen","penalty",0,0.0,1,0.2);
1280 //____________________________________________________________________
1281 void AliUnfolding::DrawResults(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TCanvas *canv, Int_t reuseHists,TH1 *unfolded)
1283 // Draw histograms of
1284 // - Result folded with response
1285 // - Penalty factors
1286 // - Chisquare contributions
1287 // (Useful for debugging/sanity checks and the interactive unfolder)
1289 // If a canvas pointer is given (canv != 0), it will be used for all
1290 // plots; 3 pads are made if needed.
1293 Int_t blankCanvas = 0;
1294 TVirtualPad *presult = 0;
1295 TVirtualPad *pres = 0;
1296 TVirtualPad *ppen = 0;
1300 presult = canv->GetPad(1);
1301 pres = canv->GetPad(2);
1302 ppen = canv->GetPad(3);
1303 if (presult == 0 || pres == 0 || ppen == 0) {
1307 presult = canv->GetPad(1);
1308 pres = canv->GetPad(2);
1309 ppen = canv->GetPad(3);
1313 presult = new TCanvas;
1319 if (fgMaxInput == -1)
1321 Printf("AliUnfolding::Unfold: WARNING. Number of measured bins not set with SetNbins. Using number of bins in measured distribution");
1322 fgMaxInput = measured->GetNbinsX();
1324 if (fgMaxParams == -1)
1326 Printf("AliUnfolding::Unfold: WARNING. Number of unfolded bins not set with SetNbins. Using number of bins in measured distribution");
1327 fgMaxParams = initialConditions->GetNbinsX();
1330 if (fgOverflowBinLimit > 0)
1331 CreateOverflowBin(correlation, measured);
1333 // Uses Minuit methods
1335 SetStaticVariables(correlation, measured, efficiency);
1337 Double_t* params = new Double_t[fgMaxParams+1];
1338 for (Int_t i=0; i<fgMaxParams; ++i)
1340 params[i] = initialConditions->GetBinContent(i+1) * efficiency->GetBinContent(i+1);
1342 Bool_t fix = kFALSE;
1346 params[i] = -params[i];
1348 params[i]=TMath::Sqrt(params[i]);
1354 //fgPrintChi2Details = kTRUE;
1355 fgCallCount = 0; // To make sure that Chi2Function prints the components
1356 Chi2Function(dummy, 0, chi2, params, 0);
1359 TH1 *meas2 = (TH1*)measured->Clone("meas2");
1360 meas2->SetTitle("meas2");
1361 meas2->SetName("meas2");
1362 meas2->SetLineColor(1);
1363 meas2->SetLineWidth(2);
1364 meas2->SetMarkerColor(meas2->GetLineColor());
1365 meas2->SetMarkerStyle(20);
1366 Float_t scale = unfolded->GetXaxis()->GetBinWidth(1)/meas2->GetXaxis()->GetBinWidth(1);
1367 meas2->Scale(scale);
1371 meas2->DrawCopy("same");
1372 //meas2->GetXaxis()->SetLimits(0,fgMaxInput);
1373 meas2->SetBit(kCannotPick);
1374 DrawGuess(params, presult, pres, ppen, reuseHists,unfolded);
1377 //____________________________________________________________________
1378 void AliUnfolding::RedrawInteractive() {
1380 // Helper function for interactive unfolding
1382 DrawResults(fghCorrelation,fghEfficiency,fghMeasured,fghUnfolded,fgCanvas,1,fghUnfolded);
1384 //____________________________________________________________________
1385 void AliUnfolding::InteractiveUnfold(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions)
1388 // Function to do interactive unfolding
1389 // A canvas is drawn with the unfolding result
1390 // Change the histogram with your mouse and all histograms
1391 // will be updated automatically
1393 fgCanvas = new TCanvas("UnfoldingCanvas","Interactive unfolding",500,800);
1399 delete fghUnfolded; fghUnfolded = 0;
1402 gROOT->SetEditHistograms(kTRUE);
1404 fghUnfolded = new TH1F("AliUnfoldingfghUnfolded","Unfolded result (Interactive unfolder",efficiency->GetNbinsX(),efficiency->GetXaxis()->GetXmin(),efficiency->GetXaxis()->GetXmax());
1406 fghCorrelation = correlation;
1407 fghEfficiency = efficiency;
1408 fghMeasured = measured;
1410 if(fgMinuitMaxIterations>0)
1411 UnfoldWithMinuit(correlation,efficiency,measured,initialConditions,fghUnfolded,kFALSE);
1413 fghUnfolded = initialConditions;
1415 fghUnfolded->SetLineColor(2);
1416 fghUnfolded->SetMarkerColor(2);
1417 fghUnfolded->SetLineWidth(2);
1421 fghUnfolded->Draw();
1422 DrawResults(correlation,efficiency,measured,fghUnfolded,fgCanvas,kFALSE,fghUnfolded);
1424 TExec *execRedraw = new TExec("redraw","AliUnfolding::RedrawInteractive()");
1425 fghUnfolded->GetListOfFunctions()->Add(execRedraw);
1427 //____________________________________________________________________
1428 void AliUnfolding::DrawGuess(Double_t *params, TVirtualPad *pfolded, TVirtualPad *pres, TVirtualPad *ppen, Int_t reuseHists,TH1* unfolded)
1431 // draws residuals of solution suggested by params and effect of regularization
1435 pfolded = new TCanvas;
1442 TVectorD paramsVector(fgMaxParams);
1443 for (Int_t i=0; i<fgMaxParams; ++i)
1444 paramsVector[i] = params[i] * params[i];
1447 TVectorD measGuessVector(fgMaxInput);
1448 measGuessVector = (*fgCorrelationMatrix) * paramsVector;
1452 folded = dynamic_cast<TH1*>(gROOT->FindObject("__hfolded"));
1453 if (!reuseHists || folded == 0) {
1454 if (fgMeasuredAxis->GetXbins()->GetArray()) // variable bins
1455 folded = new TH1F("__hfolded","Folded histo from AliUnfolding",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXbins()->GetArray());
1457 folded = new TH1F("__hfolded","Folded histo from AliUnfolding",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXmin(),fgMeasuredAxis->GetXmax());
1460 folded->SetBit(kCannotPick);
1461 folded->SetLineColor(4);
1462 folded->SetLineWidth(2);
1464 for (Int_t ibin =0; ibin < fgMaxInput; ibin++)
1465 folded->SetBinContent(ibin+1, measGuessVector[ibin]);
1467 Float_t scale = unfolded->GetXaxis()->GetBinWidth(1)/folded->GetXaxis()->GetBinWidth(1);
1468 folded->Scale(scale);
1472 folded->Draw("same");
1475 measGuessVector -= (*fgCurrentESDVector);
1477 TVectorD copy(measGuessVector);
1480 // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
1481 // normal way is like this:
1482 // measGuessVector *= (*fgCorrelationCovarianceMatrix);
1483 // optimized way like this:
1484 for (Int_t i=0; i<fgMaxInput; ++i)
1485 measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
1487 // (Ad - m) W (Ad - m)
1488 // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
1489 // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see ApplyMinuitFit)
1490 //Double_t chi2FromFit = measGuessVector * copy * 1e6;
1493 // Double_t pTarray[fgMaxParams+1];
1494 // for(int i=0; i<fgMaxParams; i++) {
1495 // pTarray[i] = fgUnfoldedAxis->GetBinCenter(i)-0.5*fgUnfoldedAxis->GetBinWidth(i);
1497 // pTarray[fgMaxParams] = fgUnfoldedAxis->GetBinCenter(fgMaxParams-1)+0.5*fgUnfoldedAxis->GetBinWidth(fgMaxParams-1);
1498 // TH1* residuals = new TH1F("residuals", "residuals", fgMaxParams,pTarray);
1499 // TH1* residuals = new TH1F("residuals", "residuals", fgMaxInput, -0.5, fgMaxInput - 0.5);
1500 // for (Int_t i=0; i<fgMaxInput; ++i)
1501 // residuals->SetBinContent(i+1, measGuessVector(i) * copy(i) * 1e6);7
1502 TH1* residuals = GetResidualsPlot(params);
1503 residuals->GetXaxis()->SetTitleSize(0.06);
1504 residuals->GetXaxis()->SetTitleOffset(0.7);
1505 residuals->GetXaxis()->SetLabelSize(0.07);
1506 residuals->GetYaxis()->SetTitleSize(0.08);
1507 residuals->GetYaxis()->SetTitleOffset(0.5);
1508 residuals->GetYaxis()->SetLabelSize(0.06);
1509 pres->cd(); residuals->DrawCopy(); gPad->SetLogy();
1513 TH1* penalty = GetPenaltyPlot(params);
1514 penalty->GetXaxis()->SetTitleSize(0.06);
1515 penalty->GetXaxis()->SetTitleOffset(0.7);
1516 penalty->GetXaxis()->SetLabelSize(0.07);
1517 penalty->GetYaxis()->SetTitleSize(0.08);
1518 penalty->GetYaxis()->SetTitleOffset(0.5);
1519 penalty->GetYaxis()->SetLabelSize(0.06);
1520 ppen->cd(); penalty->DrawCopy(); gPad->SetLogy();
1523 //____________________________________________________________________
1524 TH1* AliUnfolding::GetResidualsPlot(TH1* corrected)
1527 // MvL: THIS MUST BE INCORRECT.
1528 // Need heff to calculate params from TH1 'corrected'
1532 // fill residuals histogram of solution suggested by params and effect of regularization
1535 Double_t* params = new Double_t[fgMaxParams];
1536 for (Int_t i=0; i<fgMaxParams; i++)
1539 for (Int_t i=0; i<TMath::Min(fgMaxParams, corrected->GetNbinsX()); i++)
1540 params[i] = TMath::Sqrt(TMath::Abs(corrected->GetBinContent(i+1)*(*fgEfficiency)(i)));
1542 TH1 * plot = GetResidualsPlot(params);
1547 //____________________________________________________________________
1548 TH1* AliUnfolding::GetResidualsPlot(Double_t* params)
1551 // fill residuals histogram of solution suggested by params and effect of regularization
1555 TVectorD paramsVector(fgMaxParams);
1556 for (Int_t i=0; i<fgMaxParams; ++i)
1557 paramsVector[i] = params[i] * params[i];
1560 TVectorD measGuessVector(fgMaxInput);
1561 measGuessVector = (*fgCorrelationMatrix) * paramsVector;
1564 measGuessVector -= (*fgCurrentESDVector);
1566 TVectorD copy(measGuessVector);
1569 // this step can be optimized because currently only the diagonal elements of fgCorrelationCovarianceMatrix are used
1570 // normal way is like this:
1571 // measGuessVector *= (*fgCorrelationCovarianceMatrix);
1572 // optimized way like this:
1573 for (Int_t i=0; i<fgMaxInput; ++i)
1574 measGuessVector[i] *= (*fgCorrelationCovarianceMatrix)(i, i);
1576 // (Ad - m) W (Ad - m)
1577 // the factor 1e6 prevents that a very small number (measGuessVector[i]) is multiplied with a very
1578 // big number ((*fgCorrelationCovarianceMatrix)(i, i)) (see ApplyMinuitFit)
1579 //Double_t chi2FromFit = measGuessVector * copy * 1e6;
1583 if (fgMeasuredAxis->GetXbins()->GetArray()) // variable bins
1584 residuals = new TH1F("residuals", "residuals;unfolded pos;residual",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXbins()->GetArray());
1586 residuals = new TH1F("residuals", "residuals;unfolded pos;residual",fgMeasuredAxis->GetNbins(),fgMeasuredAxis->GetXmin(), fgMeasuredAxis->GetXmax());
1587 // TH1* residuals = new TH1F("residuals", "residuals", fgMaxInput, -0.5, fgMaxInput - 0.5);
1589 Double_t sumResiduals = 0.;
1590 for (Int_t i=0; i<fgMaxInput; ++i) {
1591 residuals->SetBinContent(i+1, measGuessVector(i) * copy(i) * 1e6);
1592 sumResiduals += measGuessVector(i) * copy(i) * 1e6;
1594 fAvgResidual = sumResiduals/(double)fgMaxInput;
1599 //____________________________________________________________________
1600 TH1* AliUnfolding::GetPenaltyPlot(TH1* corrected)
1602 // draws the penalty factors as function of multiplicity of the current selected regularization
1604 Double_t* params = new Double_t[fgMaxParams];
1605 for (Int_t i=0; i<fgMaxParams; i++)
1608 for (Int_t i=0; i<TMath::Min(fgMaxParams, corrected->GetNbinsX()); i++)
1609 params[i] = (*fgEfficiency)(i)*corrected->GetBinContent(i+1);
1611 TH1* penalty = GetPenaltyPlot(params);
1618 //____________________________________________________________________
1619 TH1* AliUnfolding::GetPenaltyPlot(Double_t* params)
1621 // draws the penalty factors as function of multiplicity of the current selected regularization
1623 //TH1* penalty = new TH1F("penalty", ";unfolded multiplicity;penalty factor", fgMaxParams, -0.5, fgMaxParams - 0.5);
1624 // 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) );
1627 if (fgUnfoldedAxis->GetXbins()->GetArray())
1628 penalty = new TH1F("penalty", ";unfolded pos;penalty factor", fgUnfoldedAxis->GetNbins(),fgUnfoldedAxis->GetXbins()->GetArray());
1630 penalty = new TH1F("penalty", ";unfolded pos;penalty factor", fgUnfoldedAxis->GetNbins(),fgUnfoldedAxis->GetXmin(),fgUnfoldedAxis->GetXmax());
1632 for (Int_t i=1+fgSkipBinsBegin; i<fgMaxParams; ++i)
1635 if (fgRegularizationType == kPol0)
1637 Double_t right = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
1638 Double_t left = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
1642 Double_t diffTmp = (right - left);
1643 diff = diffTmp * diffTmp / left / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2) / 100;
1646 if (fgRegularizationType == kPol1 && i > 1+fgSkipBinsBegin)
1648 if (params[i-1] == 0)
1651 Double_t right = params[i] / fgUnfoldedAxis->GetBinWidth(i+1);
1652 Double_t middle = params[i-1] / fgUnfoldedAxis->GetBinWidth(i);
1653 Double_t left = params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1);
1655 Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
1656 Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
1658 diff = (der1 - der2) * (der1 - der2) / middle;
1661 if (fgRegularizationType == kLog && i > 1+fgSkipBinsBegin)
1663 if (params[i-1] == 0)
1666 Double_t right = log(params[i]);
1667 Double_t middle = log(params[i-1]);
1668 Double_t left = log(params[i-2]);
1670 Double_t der1 = (right - middle);
1671 Double_t der2 = (middle - left);
1673 //Double_t error = 1. / params[i] + 4. / params[i-1] + 1. / params[i-2];
1674 //Printf("%d %f %f", i, (der1 - der2) * (der1 - der2), error);
1676 diff = (der1 - der2) * (der1 - der2);// / error;
1678 if (fgRegularizationType == kCurvature && i > 1+fgSkipBinsBegin)
1680 Double_t right = params[i]; // params are sqrt
1681 Double_t middle = params[i-1];
1682 Double_t left = params[i-2];
1684 Double_t der1 = (right - middle)/0.5/(fgUnfoldedAxis->GetBinWidth(i-1) + fgUnfoldedAxis->GetBinWidth(i));
1685 Double_t der2 = (middle - left)/0.5/(fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i+1));
1687 diff = (der1 - der2)/(fgUnfoldedAxis->GetBinWidth(i-1) + fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1))*3.0;
1688 diff = 1e4*diff*diff;
1690 if (fgRegularizationType == kPowerLaw && i > 1+fgSkipBinsBegin)
1693 if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
1696 if (fgUnfoldedAxis->GetBinWidth(i+1) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8 || fgUnfoldedAxis->GetBinWidth(i) < 1e-8)
1699 double middle = TMath::Power(params[i-1] / fgUnfoldedAxis->GetBinWidth(i),fgPowern);
1702 double right = TMath::Power(params[i] / fgUnfoldedAxis->GetBinWidth(i+1),fgPowern)/middle;
1704 double left = TMath::Power(params[i-2] / fgUnfoldedAxis->GetBinWidth(i-1),fgPowern)/middle;
1708 Double_t der1 = (right - middle) / ((fgUnfoldedAxis->GetBinWidth(i+1) + fgUnfoldedAxis->GetBinWidth(i)) / 2);
1709 Double_t der2 = (middle - left) / ((fgUnfoldedAxis->GetBinWidth(i) + fgUnfoldedAxis->GetBinWidth(i-1)) / 2);
1711 diff = (der1 - der2) * (der1 - der2);// / error;
1715 if (fgRegularizationType == kLogLog && i > 1+fgSkipBinsBegin)
1718 if (params[i] < 1e-8 || params[i-1] < 1e-8 || params[i-2] < 1e-8)
1721 Double_t right = log(params[i] / (*fgEfficiency)(i) / fgUnfoldedAxis->GetBinWidth(i+1));
1722 Double_t middle = log(params[i-1] / (*fgEfficiency)(i-1) / fgUnfoldedAxis->GetBinWidth(i));
1723 Double_t left = log(params[i-2] / (*fgEfficiency)(i-2) / fgUnfoldedAxis->GetBinWidth(i-1));
1725 Double_t der1 = (right - middle) / ( log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i)) );
1726 Double_t der2 = (middle - left) /( log(fgUnfoldedAxis->GetBinCenter(i)) - log(fgUnfoldedAxis->GetBinCenter(i-1)) );
1728 double tmp = (log(fgUnfoldedAxis->GetBinCenter(i+1)) - log(fgUnfoldedAxis->GetBinCenter(i-1)))/2.;
1729 Double_t dder = (der1-der2) / tmp;
1734 penalty->SetBinContent(i, diff*fgRegularizationWeight);
1736 //Printf("%d %f %f %f %f", i-1, left, middle, right, diff);
1742 //____________________________________________________________________
1743 void AliUnfolding::TF1Function(Int_t& unused1, Double_t* unused2, Double_t& chi2, Double_t *params, Int_t unused3)
1746 // fit function for minuit
1747 // uses the TF1 stored in fgFitFunction
1750 for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
1751 fgFitFunction->SetParameter(i, params[i]);
1753 Double_t* params2 = new Double_t[fgMaxParams];
1755 for (Int_t i=0; i<fgMaxParams; ++i)
1756 params2[i] = fgFitFunction->Eval(i);
1758 Chi2Function(unused1, unused2, chi2, params2, unused3);
1766 //____________________________________________________________________
1767 Int_t AliUnfolding::UnfoldWithFunction(TH2* correlation, TH1* efficiency, TH1* measured, TH1* /* initialConditions */, TH1* aResult)
1770 // correct spectrum using minuit chi2 method applying a functional fit
1775 Printf("AliUnfolding::UnfoldWithFunction: ERROR fit function not set. Exiting.");
1779 SetChi2Regularization(kNone, 0);
1781 SetStaticVariables(correlation, measured, efficiency);
1783 // Initialize TMinuit via generic fitter interface
1784 TVirtualFitter *minuit = TVirtualFitter::Fitter(0, fgFitFunction->GetNpar());
1786 minuit->SetFCN(TF1Function);
1787 for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
1789 Double_t lower, upper;
1790 fgFitFunction->GetParLimits(i, lower, upper);
1791 minuit->SetParameter(i, Form("param%d",i), fgFitFunction->GetParameter(i), fgMinuitStepSize, lower, upper);
1794 Double_t arglist[100];
1796 minuit->ExecuteCommand("SET PRINT", arglist, 1);
1797 minuit->ExecuteCommand("SCAN", arglist, 0);
1798 minuit->ExecuteCommand("MIGRAD", arglist, 0);
1799 //minuit->ExecuteCommand("MINOS", arglist, 0);
1801 for (Int_t i=0; i<fgFitFunction->GetNpar(); i++)
1802 fgFitFunction->SetParameter(i, minuit->GetParameter(i));
1804 for (Int_t i=0; i<fgMaxParams; ++i)
1806 Double_t value = fgFitFunction->Eval(i);
1808 Printf("%d : %f", i, value);
1812 if (efficiency->GetBinContent(i+1) > 0)
1814 value /= efficiency->GetBinContent(i+1);
1819 aResult->SetBinContent(i+1, value);
1820 aResult->SetBinError(i+1, 0);
1826 //____________________________________________________________________
1827 void AliUnfolding::CreateOverflowBin(TH2* correlation, TH1* measured)
1829 // Finds the first bin where the content is below fgStatLimit and combines all values for this bin and larger bins
1830 // The same limit is applied to the correlation
1833 for (Int_t i=1; i<=measured->GetNbinsX(); ++i)
1835 if (measured->GetBinContent(i) <= fgOverflowBinLimit)
1844 Printf("AliUnfolding::CreateOverflowBin: WARNING: First bin is already below limit of %f", fgOverflowBinLimit);
1848 Printf("AliUnfolding::CreateOverflowBin: Bin limit in measured spectrum determined to be %d", lastBin);
1850 TCanvas* canvas = 0;
1854 canvas = new TCanvas("StatSolution", "StatSolution", 1000, 800);
1855 canvas->Divide(2, 2);
1858 measured->SetStats(kFALSE);
1859 measured->DrawCopy();
1863 correlation->SetStats(kFALSE);
1864 correlation->DrawCopy("COLZ");
1867 measured->SetBinContent(lastBin, measured->Integral(lastBin, measured->GetNbinsX()));
1868 for (Int_t i=lastBin+1; i<=measured->GetNbinsX(); ++i)
1870 measured->SetBinContent(i, 0);
1871 measured->SetBinError(i, 0);
1873 // the error is set to sqrt(N), better solution possible?, sum of relative errors of all contributions???
1874 measured->SetBinError(lastBin, TMath::Sqrt(measured->GetBinContent(lastBin)));
1876 Printf("AliUnfolding::CreateOverflowBin: This bin has now %f +- %f entries", measured->GetBinContent(lastBin), measured->GetBinError(lastBin));
1878 for (Int_t i=1; i<=correlation->GetNbinsX(); ++i)
1880 correlation->SetBinContent(i, lastBin, correlation->Integral(i, i, lastBin, correlation->GetNbinsY()));
1881 // the error is set to sqrt(N), better solution possible?, sum of relative errors of all contributions???
1882 correlation->SetBinError(i, lastBin, TMath::Sqrt(correlation->GetBinContent(i, lastBin)));
1884 for (Int_t j=lastBin+1; j<=correlation->GetNbinsY(); ++j)
1886 correlation->SetBinContent(i, j, 0);
1887 correlation->SetBinError(i, j, 0);
1891 Printf("AliUnfolding::CreateOverflowBin: Adjusted correlation matrix!");
1896 measured->DrawCopy();
1900 correlation->DrawCopy("COLZ");
1904 Int_t AliUnfolding::UnfoldGetBias(TH2* correlation, TH1* efficiency, TH1* measured, TH1* initialConditions, TH1* result)
1906 // unfolds and assigns bias as errors with Eq. 19 of Cowan, "a survey of unfolding methods for particle physics"
1907 // b_i = sum_j dmu_i/dn_j (nu_j - n_j) with nu_j as folded guess, n_j as data
1908 // dmu_i/dn_j is found numerically by changing the bin content and re-unfolding
1910 // return code: 0 (success) -1 (error: from Unfold(...) )
1912 if (Unfold(correlation, efficiency, measured, initialConditions, result) != 0)
1915 TMatrixD matrix(fgMaxInput, fgMaxParams);
1918 for (Int_t loop=0; loop<4; loop++)
1919 newResult[loop] = (TH1*) result->Clone(Form("newresult_%d", loop));
1921 // change bin-by-bin and built matrix of effects
1922 for (Int_t m=0; m<fgMaxInput; m++)
1924 if (measured->GetBinContent(m+1) < 1)
1927 for (Int_t loop=0; loop<4; loop++)
1929 //Printf("%d %d", i, loop);
1933 case 0: factor = 0.5; break;
1934 case 1: factor = -0.5; break;
1935 case 2: factor = 1; break;
1936 case 3: factor = -1; break;
1940 TH1* measuredClone = (TH1*) measured->Clone("measuredClone");
1942 measuredClone->SetBinContent(m+1, measured->GetBinContent(m+1) + factor * TMath::Sqrt(measured->GetBinContent(m+1)));
1943 //new TCanvas; measuredClone->Draw("COLZ");
1945 newResult[loop]->Reset();
1946 Unfold(correlation, efficiency, measuredClone, measuredClone, newResult[loop]);
1947 // WARNING if we leave here, then nothing is calculated
1950 delete measuredClone;
1953 for (Int_t t=0; t<fgMaxParams; t++)
1956 //matrix(i, j) = (result->GetBinContent(j+1) - newResult->GetBinContent(j+1)) / TMath::Sqrt(changedHist->GetBinContent(1, i+1));
1958 // four values from the derivate (procedure taken from ROOT -- suggestion by Ruben)
1959 // = 1/2D * [ 8 (f(D/2) - f(-D/2)) - (f(D)-f(-D)) ]/3
1963 measured->SetBinContent(1, m+1, 100);
1964 newResult[0]->SetBinContent(t+1, 5 * 0.5 + 10);
1965 newResult[1]->SetBinContent(t+1, 5 * -0.5 + 10);
1966 newResult[2]->SetBinContent(t+1, 5 * 1 + 10);
1967 newResult[3]->SetBinContent(t+1, 5 * -1 + 10);
1970 matrix(m, t) = 0.5 / TMath::Sqrt(measured->GetBinContent(m+1)) *
1971 ( 8. * (newResult[0]->GetBinContent(t+1) - newResult[1]->GetBinContent(t+1)) -
1972 (newResult[2]->GetBinContent(t+1) - newResult[3]->GetBinContent(t+1))
1977 for (Int_t loop=0; loop<4; loop++)
1978 delete newResult[loop];
1980 // ...calculate folded guess
1981 TH1* convoluted = (TH1*) measured->Clone("convoluted");
1982 convoluted->Reset();
1983 for (Int_t m=0; m<fgMaxInput; m++)
1986 for (Int_t t = 0; t<fgMaxParams; t++)
1988 Float_t tmp = correlation->GetBinContent(t+1, m+1) * result->GetBinContent(t+1);
1990 tmp *= efficiency->GetBinContent(t+1);
1993 convoluted->SetBinContent(m+1, value);
1997 convoluted->Scale(measured->Integral() / convoluted->Integral());
1999 //new TCanvas; convoluted->Draw(); measured->Draw("SAME"); measured->SetLineColor(2);
2003 convoluted->Add(measured, -1);
2006 for (Int_t t = 0; t<fgMaxParams; t++)
2009 for (Int_t m=0; m<fgMaxInput; m++)
2010 error += matrix(m, t) * convoluted->GetBinContent(m+1);
2011 result->SetBinError(t+1, error);
2014 //new TCanvas; result->Draw(); gPad->SetLogy();