Barlow and C. Beeston, Comp. Phys. Comm. 77 (1993) 219-228,
Minuit2Minimizer: Minimize with max-calls 1345 convergence for edm < 0.01 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -10009.6795619724871
Edm = 3.20755056433236687e-07
Nfcn = 153
frac0 = 3.70672e-05 +/- 0.613661 (limited)
frac1 = 0.283223 +/- 0.0555675 (limited)
frac2 = 0.716752 +/- 0.0534245 (limited)
Status: 0
Parameter 0: true 0.010, estim. 0.000 +/- 0.614
Parameter 1: true 0.300, estim. 0.283 +/- 0.056
Parameter 2: true 0.690, estim. 0.717 +/- 0.053
#include <iostream>
void fitFraction()
{
TF1 f0(
"f0",
"[0]*(1-cos(x))/TMath::Pi()", 0.,
TMath::Pi());
f0.SetParameter(0, 1.);
f0.SetLineColor(2);
TF1 f1(
"f1",
"[0]*(1-cos(x)*cos(x))*2./TMath::Pi()", 0.,
TMath::Pi());
TF1 f2(
"f2",
"[0]*(1+cos(x))/TMath::Pi()", 0.,
TMath::Pi());
f2.SetParameter(0, 1.);
f2.SetLineColor(4);
TH1F data(
"data",
"Data angle distribution", nBins, 0,
TMath::Pi());
data.SetXTitle("x");
data.SetMarkerStyle(20);
data.SetMarkerSize(.7);
data.SetMinimum(0);
htruemc0.SetLineColor(2);
htruemc1.SetLineColor(3);
htruemc2.SetLineColor(4);
for (
Int_t i = 0; i < Ndata; i++) {
if (p < trueP0) {
} else if (p < trueP0 + trueP1) {
} else {
}
}
TH1F mc0(
"mc0",
"MC sample 0 angle distribution", nBins, 0,
TMath::Pi());
mc0.SetXTitle("x");
mc0.SetLineColor(2);
mc0.SetMarkerColor(2);
mc0.SetMarkerStyle(24);
mc0.SetMarkerSize(.7);
for (
Int_t i = 0; i < N0; i++) {
mc0.Fill(f0.GetRandom());
}
TH1F mc1(
"mc1",
"MC sample 1 angle distribution", nBins, 0,
TMath::Pi());
mc1.SetXTitle("x");
mc1.SetLineColor(3);
mc1.SetMarkerColor(3);
mc1.SetMarkerStyle(24);
mc1.SetMarkerSize(.7);
for (
Int_t i = 0; i < N1; i++) {
mc1.Fill(
f1.GetRandom());
}
TH1F mc2(
"mc2",
"MC sample 2 angle distribution", nBins, 0,
TMath::Pi());
mc2.SetXTitle("x");
mc2.SetLineColor(4);
mc2.SetMarkerColor(4);
mc2.SetMarkerStyle(24);
mc2.SetMarkerSize(.7);
for (
Int_t i = 0; i < N2; i++) {
mc2.Fill(f2.GetRandom());
}
mc.Add(&mc0);
mc.Add(&mc1);
mc.Add(&mc2);
fit.Constrain(0, 0.0, 1.0);
fit.Constrain(1, 0.0, 1.0);
fit.Constrain(2, 0.0, 1.0);
Int_t status = fit.Fit();
std::cout << "Status: " << status << std::endl;
TCanvas c(
"c",
"FractionFitter example", 700, 700);
f0.DrawClone();
f0.GetHistogram()->SetTitle("Original MC distributions");
f2.DrawClone("same");
data.SetTitle("Data distribution with true contributions");
data.DrawClone("EP");
htruemc0.Draw("same");
htruemc1.Draw("same");
htruemc2.Draw("same");
mc0.SetTitle("MC generated samples with fit predictions");
mc0.Draw("PE");
mc1.Draw("PEsame");
mc2.Draw("PEsame");
if (status == 0) {
auto mcp0 = (
TH1F *)fit.GetMCPrediction(0);
mcp0->SetLineColor(2);
mcp0->Draw("same");
auto mcp1 = (
TH1F *)fit.GetMCPrediction(1);
mcp1->SetLineColor(3);
mcp1->Draw("same");
auto mcp2 = (
TH1F *)fit.GetMCPrediction(2);
mcp2->SetLineColor(4);
mcp2->Draw("same");
}
Double_t p0, p1, p2, errP0, errP1, errP2;
if (status == 0) {
auto result = (
TH1F *)fit.GetPlot();
fit.GetResult(0, p0, errP0);
printf(" Parameter %d: true %.3f, estim. %.3f +/- %.3f\n", 0, trueP0, p0, errP0);
fit.GetResult(1, p1, errP1);
printf(" Parameter %d: true %.3f, estim. %.3f +/- %.3f\n", 1, trueP1, p1, errP1);
fit.GetResult(2, p2, errP2);
printf(" Parameter %d: true %.3f, estim. %.3f +/- %.3f\n", 2, trueP2, p2, errP2);
data.SetTitle("Data distribution with fitted contributions");
data.DrawClone("Ep");
result->Draw("same");
f0.SetParameter(0, Ndata * p0 / int0 * data.GetBinWidth(1));
f0.SetLineStyle(2);
f0.DrawClone("same");
f1.SetParameter(0, Ndata * p1 / int1 * data.GetBinWidth(1));
f2.SetParameter(0, Ndata * p2 / int2 * data.GetBinWidth(1));
f2.SetLineStyle(2);
f2.DrawClone("same");
sprintf(
text,
"%d: true %.2f, estimated %.2f +/- %.2f\n", 0, trueP0, p0, errP0);
l.DrawTextNDC(.45, .30,
text);
sprintf(
text,
"%d: true %.2f, estimated %.2f +/- %.2f\n", 1, trueP1, p1, errP1);
l.DrawTextNDC(.45, .25,
text);
sprintf(
text,
"%d: true %.2f, estimated %.2f +/- %.2f\n", 2, trueP2, p2, errP2);
l.DrawTextNDC(.45, .20,
text);
}
auto cnew =
c.DrawClone();
}
int Int_t
Signed integer 4 bytes (int).
char Char_t
Character 1 byte (char).
double Double_t
Double 8 bytes.
Fits MC fractions to data histogram.
1-D histogram with a float per channel (see TH1 documentation)