Implementing the Barlow-Beeston method for taking into account the statistical uncertainty of a Monte-Carlo fit template.
import ROOT
x = ROOT.RooRealVar("x", "x", -20, 20)
x.setBins(25)
meanG = ROOT.RooRealVar("meanG", "meanG", 1, -10, 10)
sigG = ROOT.RooRealVar("sigG", "sigG", 1.5, -10, 10)
g = ROOT.RooGaussian("g", "Gauss", x, meanG, sigG)
u = ROOT.RooUniform("u", "Uniform", x)
sigData = g.generate(x, 50)
bkgData = u.generate(x, 1000)
sumData = ROOT.RooDataSet("sumData", "Gauss + Uniform", x)
sumData.append(sigData)
sumData.append(bkgData)
dh_sig = g.generateBinned(x, 50)
dh_bkg = u.generateBinned(x, 10000)
p_h_sig = ROOT.RooHistFunc("p_h_sig", "p_h_sig", x, dh_sig)
p_h_bkg = ROOT.RooHistFunc("p_h_bkg", "p_h_bkg", x, dh_bkg)
Asig0 = ROOT.RooRealVar("Asig", "Asig", 1, 0.01, 5000)
Abkg0 = ROOT.RooRealVar("Abkg", "Abkg", 1, 0.01, 5000)
model0 = ROOT.RooRealSumPdf("model0", "model0", [p_h_sig, p_h_bkg], [Asig0, Abkg0], True)
p_ph_sig1 = ROOT.RooParamHistFunc("p_ph_sig", "p_ph_sig", dh_sig)
p_ph_bkg1 = ROOT.RooParamHistFunc("p_ph_bkg", "p_ph_bkg", dh_bkg)
Asig1 = ROOT.RooRealVar("Asig", "Asig", 1, 0.01, 5000)
Abkg1 = ROOT.RooRealVar("Abkg", "Abkg", 1, 0.01, 5000)
model_tmp = ROOT.RooRealSumPdf("sp_ph", "sp_ph", [p_ph_sig1, p_ph_bkg1], [Asig1, Abkg1], True)
hc_sig = ROOT.RooHistConstraint("hc_sig", "hc_sig", p_ph_sig1)
hc_bkg = ROOT.RooHistConstraint("hc_bkg", "hc_bkg", p_ph_bkg1)
model1 = ROOT.RooProdPdf("model1", "model1", {hc_sig, hc_bkg}, Conditional=(model_tmp, x))
p_ph_sig2 = ROOT.RooParamHistFunc("p_ph_sig2", "p_ph_sig2", dh_sig)
p_ph_bkg2 = ROOT.RooParamHistFunc("p_ph_bkg2", "p_ph_bkg2", dh_bkg, p_ph_sig2, True)
Asig2 = ROOT.RooRealVar("Asig", "Asig", 1, 0.01, 5000)
Abkg2 = ROOT.RooRealVar("Abkg", "Abkg", 1, 0.01, 5000)
model2_tmp = ROOT.RooRealSumPdf("sp_ph", "sp_ph", [p_ph_sig2, p_ph_bkg2], [Asig2, Abkg2], True)
hc_sigbkg = ROOT.RooHistConstraint("hc_sigbkg", "hc_sigbkg", {p_ph_sig2, p_ph_bkg2})
model2 = ROOT.RooProdPdf("model2", "model2", hc_sigbkg, Conditional=(model2_tmp, x))
result0 = model0.fitTo(sumData, PrintLevel=0, Save=True)
result1 = model1.fitTo(sumData, PrintLevel=0, Save=True)
result2 = model2.fitTo(sumData, PrintLevel=0, Save=True)
can = ROOT.TCanvas("can", "", 1500, 600)
can.Divide(3, 1)
pt = ROOT.TPaveText(-19.5, 1, -2, 25)
pt.SetFillStyle(0)
pt.SetBorderSize(0)
can.cd(1)
frame = x.frame(Title="No template uncertainties")
sumData.plotOn(frame)
model0.plotOn(frame, LineColor="b", VisualizeError=result0)
sumData.plotOn(frame)
model0.plotOn(frame, LineColor="b")
p_ph_sig_set = {p_h_sig}
p_ph_bkg_set = {p_h_bkg}
model0.plotOn(frame, Components=p_ph_sig_set, LineColor="kAzure")
model0.plotOn(frame, Components=p_ph_bkg_set, LineColor="r")
model0.paramOn(frame)
sigData.plotOn(frame, MarkerColor="b")
frame.Draw()
pt_text1 = [
"No template uncertainties",
"are taken into account.",
"This leads to low errors",
"for the parameters A, since",
"the only source of errors",
"are the data statistics.",
]
for text in pt_text1:
pt.AddText(text)
pt.DrawClone()
can.cd(2)
frame = x.frame(Title="Barlow Beeston for Sig & Bkg separately")
sumData.plotOn(frame)
model1.plotOn(frame, LineColor="b", VisualizeError=result1)
sumData.plotOn(frame)
model1.plotOn(frame, LineColor="b")
p_ph_sig1_set = {p_ph_sig1}
p_ph_bkg1_set = {p_ph_bkg1}
model1.plotOn(frame, Components=p_ph_sig1_set, LineColor="kAzure")
model1.plotOn(frame, Components=p_ph_bkg1_set, LineColor="r")
model1.paramOn(frame, Parameters={Asig1, Abkg1})
sigData.plotOn(frame, MarkerColor="b")
frame.Draw()
pt.Clear()
pt_text2 = [
"With gamma parameters, the",
"signal & background templates",
"can adapt to the data.",
"Note how the blue signal",
"template changes its shape.",
"This leads to higher errors",
"of the scale parameters A.",
]
for text in pt_text2:
pt.AddText(text)
pt.DrawClone()
can.cd(3)
frame = x.frame(Title="Barlow Beeston light for (Sig+Bkg)")
sumData.plotOn(frame)
model2.plotOn(frame, LineColor="b", VisualizeError=result2)
sumData.plotOn(frame)
model2.plotOn(frame, LineColor="b")
p_ph_sig2_set = {p_ph_sig2}
p_ph_bkg2_set = {p_ph_bkg2}
model2.plotOn(frame, Components=p_ph_sig2_set, LineColor="kAzure")
model2.plotOn(frame, Components=p_ph_bkg2_set, LineColor="r")
model2.paramOn(frame, Parameters={Asig2, Abkg2})
sigData.plotOn(frame, MarkerColor="b")
frame.Draw()
pt.Clear()
pt_text3 = [
"When signal and background",
"template share one gamma para-",
"meter per bin, they adapt less.",
"The errors of the A parameters",
"also shrink slightly.",
]
for text in pt_text3:
pt.AddText(text)
pt.DrawClone()
print(
"Asig [normal ] = {} +/- {}".
format(Asig0.getVal(), Asig0.getError()))
print(
"Asig [BB ] = {} +/- {}".
format(Asig1.getVal(), Asig1.getError()))
print(
"Asig [BBlight] = {} +/- {}".
format(Asig2.getVal(), Asig2.getError()))
can.SaveAs("rf709_BarlowBeeston.png")
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
[#0] WARNING:InputArguments -- The parameter 'sigG' with range [-10, 10] of the RooGaussian 'g' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (p_h_sig,p_h_bkg)
Minuit2Minimizer: Minimize with max-calls 1000 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -2388.31039421315518
Edm = 3.25299757922590944e-06
Nfcn = 60
Abkg = 0.0614547 +/- 0.00211669 (limited)
Asig = 0.833778 +/- 0.1898 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- Including the following constraint terms in minimization: (hc_bkg,hc_sig)
[#1] INFO:Minimization -- The global observables are not defined , normalize constraints with respect to the parameters (Abkg,Asig,p_ph_bkg_gamma_bin_0,p_ph_bkg_gamma_bin_1,p_ph_bkg_gamma_bin_10,p_ph_bkg_gamma_bin_11,p_ph_bkg_gamma_bin_12,p_ph_bkg_gamma_bin_13,p_ph_bkg_gamma_bin_14,p_ph_bkg_gamma_bin_15,p_ph_bkg_gamma_bin_16,p_ph_bkg_gamma_bin_17,p_ph_bkg_gamma_bin_18,p_ph_bkg_gamma_bin_19,p_ph_bkg_gamma_bin_2,p_ph_bkg_gamma_bin_20,p_ph_bkg_gamma_bin_21,p_ph_bkg_gamma_bin_22,p_ph_bkg_gamma_bin_23,p_ph_bkg_gamma_bin_24,p_ph_bkg_gamma_bin_3,p_ph_bkg_gamma_bin_4,p_ph_bkg_gamma_bin_5,p_ph_bkg_gamma_bin_6,p_ph_bkg_gamma_bin_7,p_ph_bkg_gamma_bin_8,p_ph_bkg_gamma_bin_9,p_ph_sig_gamma_bin_0,p_ph_sig_gamma_bin_1,p_ph_sig_gamma_bin_10,p_ph_sig_gamma_bin_11,p_ph_sig_gamma_bin_12,p_ph_sig_gamma_bin_13,p_ph_sig_gamma_bin_14,p_ph_sig_gamma_bin_15,p_ph_sig_gamma_bin_16,p_ph_sig_gamma_bin_17,p_ph_sig_gamma_bin_18,p_ph_sig_gamma_bin_19,p_ph_sig_gamma_bin_2,p_ph_sig_gamma_bin_20,p_ph_sig_gamma_bin_21,p_ph_sig_gamma_bin_22,p_ph_sig_gamma_bin_23,p_ph_sig_gamma_bin_24,p_ph_sig_gamma_bin_3,p_ph_sig_gamma_bin_4,p_ph_sig_gamma_bin_5,p_ph_sig_gamma_bin_6,p_ph_sig_gamma_bin_7,p_ph_sig_gamma_bin_8,p_ph_sig_gamma_bin_9)
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model1_sumData_with_constr) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (p_ph_sig,p_ph_bkg)
Minuit2Minimizer: Minimize with max-calls 16000 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -2282.22801982066039
Edm = 9.75046667922085591e-08
Nfcn = 1365
Abkg = 0.061372 +/- 0.00222727 (limited)
Asig = 0.850643 +/- 0.234541 (limited)
p_ph_bkg_gamma_bin_0 = 0.998055 +/- 0.0474272 (limited)
p_ph_bkg_gamma_bin_1 = 1.0003 +/- 0.0474806 (limited)
p_ph_bkg_gamma_bin_10 = 0.980305 +/- 0.0456254 (limited)
p_ph_bkg_gamma_bin_11 = 1.00631 +/- 0.0508601 (limited)
p_ph_bkg_gamma_bin_12 = 1.00477 +/- 0.0510228 (limited)
p_ph_bkg_gamma_bin_13 = 0.991229 +/- 0.0489466 (limited)
p_ph_bkg_gamma_bin_14 = 0.99605 +/- 0.0497199 (limited)
p_ph_bkg_gamma_bin_15 = 1.00899 +/- 0.0476474 (limited)
p_ph_bkg_gamma_bin_16 = 1.02997 +/- 0.0492181 (limited)
p_ph_bkg_gamma_bin_17 = 1.00769 +/- 0.0467483 (limited)
p_ph_bkg_gamma_bin_18 = 1.01965 +/- 0.0509664 (limited)
p_ph_bkg_gamma_bin_19 = 0.975013 +/- 0.045874 (limited)
p_ph_bkg_gamma_bin_2 = 0.996158 +/- 0.0468678 (limited)
p_ph_bkg_gamma_bin_20 = 0.999858 +/- 0.0473542 (limited)
p_ph_bkg_gamma_bin_21 = 0.998271 +/- 0.0474906 (limited)
p_ph_bkg_gamma_bin_22 = 0.998324 +/- 0.0487672 (limited)
p_ph_bkg_gamma_bin_23 = 0.98395 +/- 0.0461929 (limited)
p_ph_bkg_gamma_bin_24 = 1.02664 +/- 0.047986 (limited)
p_ph_bkg_gamma_bin_3 = 1.00554 +/- 0.0495906 (limited)
p_ph_bkg_gamma_bin_4 = 0.990585 +/- 0.0483895 (limited)
p_ph_bkg_gamma_bin_5 = 0.992304 +/- 0.0482458 (limited)
p_ph_bkg_gamma_bin_6 = 1.01041 +/- 0.0497109 (limited)
p_ph_bkg_gamma_bin_7 = 0.997207 +/- 0.0460182 (limited)
p_ph_bkg_gamma_bin_8 = 0.980295 +/- 0.0463258 (limited)
p_ph_bkg_gamma_bin_9 = 1.00976 +/- 0.0478226 (limited)
p_ph_sig_gamma_bin_11 = 1.09508 +/- 0.337672 (limited)
p_ph_sig_gamma_bin_12 = 1.07036 +/- 0.220991 (limited)
p_ph_sig_gamma_bin_13 = 0.890729 +/- 0.20285 (limited)
p_ph_sig_gamma_bin_14 = 0.947873 +/- 0.35799 (limited)
p_ph_sig_gamma_bin_15 = 1.14083 +/- 1.11115 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- Including the following constraint terms in minimization: (hc_sigbkg)
[#1] INFO:Minimization -- The global observables are not defined , normalize constraints with respect to the parameters (Abkg,Asig,p_ph_sig2_gamma_bin_0,p_ph_sig2_gamma_bin_1,p_ph_sig2_gamma_bin_10,p_ph_sig2_gamma_bin_11,p_ph_sig2_gamma_bin_12,p_ph_sig2_gamma_bin_13,p_ph_sig2_gamma_bin_14,p_ph_sig2_gamma_bin_15,p_ph_sig2_gamma_bin_16,p_ph_sig2_gamma_bin_17,p_ph_sig2_gamma_bin_18,p_ph_sig2_gamma_bin_19,p_ph_sig2_gamma_bin_2,p_ph_sig2_gamma_bin_20,p_ph_sig2_gamma_bin_21,p_ph_sig2_gamma_bin_22,p_ph_sig2_gamma_bin_23,p_ph_sig2_gamma_bin_24,p_ph_sig2_gamma_bin_3,p_ph_sig2_gamma_bin_4,p_ph_sig2_gamma_bin_5,p_ph_sig2_gamma_bin_6,p_ph_sig2_gamma_bin_7,p_ph_sig2_gamma_bin_8,p_ph_sig2_gamma_bin_9)
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model2_sumData_with_constr) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (p_ph_sig2,p_ph_bkg2)
Minuit2Minimizer: Minimize with max-calls 13500 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -2291.45127145603146
Edm = 2.76473089713644995e-05
Nfcn = 1119
Abkg = 0.0614393 +/- 0.00221881 (limited)
Asig = 0.834401 +/- 0.203008 (limited)
p_ph_sig2_gamma_bin_0 = 0.997977 +/- 0.0474184 (limited)
p_ph_sig2_gamma_bin_1 = 1.00022 +/- 0.0474718 (limited)
p_ph_sig2_gamma_bin_10 = 0.980256 +/- 0.045617 (limited)
p_ph_sig2_gamma_bin_11 = 1.01067 +/- 0.0488754 (limited)
p_ph_sig2_gamma_bin_12 = 1.01175 +/- 0.0483914 (limited)
p_ph_sig2_gamma_bin_13 = 0.983606 +/- 0.0466986 (limited)
p_ph_sig2_gamma_bin_14 = 0.9948 +/- 0.0483051 (limited)
p_ph_sig2_gamma_bin_15 = 1.00965 +/- 0.0473307 (limited)
p_ph_sig2_gamma_bin_16 = 1.02983 +/- 0.0492085 (limited)
p_ph_sig2_gamma_bin_17 = 1.00759 +/- 0.0467394 (limited)
p_ph_sig2_gamma_bin_18 = 1.01954 +/- 0.0509568 (limited)
p_ph_sig2_gamma_bin_19 = 0.974972 +/- 0.0458656 (limited)
p_ph_sig2_gamma_bin_2 = 0.996083 +/- 0.0468591 (limited)
p_ph_sig2_gamma_bin_20 = 0.999777 +/- 0.0473454 (limited)
p_ph_sig2_gamma_bin_21 = 0.998193 +/- 0.0474818 (limited)
p_ph_sig2_gamma_bin_22 = 0.998245 +/- 0.0487583 (limited)
p_ph_sig2_gamma_bin_23 = 0.983895 +/- 0.0461844 (limited)
p_ph_sig2_gamma_bin_24 = 1.02651 +/- 0.0479766 (limited)
p_ph_sig2_gamma_bin_3 = 1.00545 +/- 0.0495815 (limited)
p_ph_sig2_gamma_bin_4 = 0.990519 +/- 0.0483808 (limited)
p_ph_sig2_gamma_bin_5 = 0.992236 +/- 0.048237 (limited)
p_ph_sig2_gamma_bin_6 = 1.01031 +/- 0.0497017 (limited)
p_ph_sig2_gamma_bin_7 = 0.99713 +/- 0.0460096 (limited)
p_ph_sig2_gamma_bin_8 = 0.980245 +/- 0.0463173 (limited)
p_ph_sig2_gamma_bin_9 = 1.00966 +/- 0.0478136 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model0) directly selected PDF components: (p_h_sig)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model0) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model0) directly selected PDF components: (p_h_bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model0) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 50 will supersede previous event count of 1050 for normalization of PDF projections
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model1) directly selected PDF components: (p_ph_sig)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model1) indirectly selected PDF components: (sp_ph)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model1) directly selected PDF components: (p_ph_bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model1) indirectly selected PDF components: (sp_ph)
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 50 will supersede previous event count of 1050 for normalization of PDF projections
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model2) directly selected PDF components: (p_ph_sig2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model2) indirectly selected PDF components: (sp_ph)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model2) directly selected PDF components: (p_ph_bkg2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model2) indirectly selected PDF components: (sp_ph)
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 50 will supersede previous event count of 1050 for normalization of PDF projections
Asig [normal ] = 0.8337778709310433 +/- 0.1898141852388937
Asig [BB ] = 0.8506430697754669 +/- 0.23584537196554434
Asig [BBlight] = 0.834401491499811 +/- 0.20297209693561152
Based on a demo by Wouter Verkerke
- Date
- June 2021
- Author
- Harshal Shende, Stephan Hageboeck (C++ version)
Definition in file rf709_BarlowBeeston.py.