Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf612_recoverFromInvalidParameters.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Likelihood and minimization: Recover from regions where the function is not defined.
5##
6## We demonstrate improved recovery from disallowed parameters. For this, we use a polynomial PDF of the form
7## \f[
8## \mathrm{Pol2} = \mathcal{N} \left( c + a_1 \cdot x + a_2 \cdot x^2 + 0.01 \cdot x^3 \right),
9## \f]
10## where \f$ \mathcal{N} \f$ is a normalisation factor. Unless the parameters are chosen carefully,
11## this function can be negative, and hence, it cannot be used as a PDF. In this case, RooFit passes
12## an error to the minimiser, which might try to recover.
13##
14## \macro_code
15## \macro_output
16##
17## \date June 2021
18## \author Harshal Shende, Stephan Hageboeck (C++ version)
19
20import ROOT
21
22
23# Create a fit model:
24# The polynomial is notoriously unstable, because it can quickly go negative.
25# Since PDFs need to be positive, one often ends up with an unstable fit model.
26x = ROOT.RooRealVar("x", "x", -15, 15)
27a1 = ROOT.RooRealVar("a1", "a1", -0.5, -10.0, 20.0)
28a2 = ROOT.RooRealVar("a2", "a2", 0.2, -10.0, 20.0)
29a3 = ROOT.RooRealVar("a3", "a3", 0.01)
30pdf = ROOT.RooPolynomial("pol3", "c + a1 * x + a2 * x*x + 0.01 * x*x*x", x, [a1, a2, a3])
31
32# Create toy data with all-positive coefficients:
33data = pdf.generate(x, 10000)
34
35# For plotting.
36# We create pointers to the plotted objects. We want these objects to leak out of the function,
37# so we can still see them after it returns.
38c = ROOT.TCanvas()
39frame = x.frame()
40data.plotOn(frame, Name="data")
41
42# Plotting a PDF with disallowed parameters doesn't work. We would get a lot of error messages.
43# Therefore, we disable plotting messages in RooFit's message streams:
44ROOT.RooMsgService.instance().getStream(0).removeTopic(ROOT.RooFit.Plotting)
45ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.Plotting)
46
47
48# RooFit before ROOT 6.24
49# --------------------------------
50# Before 6.24, RooFit wasn't able to recover from invalid parameters. The minimiser just errs around
51# the starting values of the parameters without finding any improvement.
52
53# Set up the parameters such that the PDF would come out negative. The PDF is now undefined.
54a1.setVal(10.0)
55a2.setVal(-1.0)
56
57# Perform a fit:
58fitWithoutRecovery = pdf.fitTo(
59 data,
60 Save=True,
61 RecoverFromUndefinedRegions=0.0, # This is how RooFit behaved prior to ROOT 6.24
62 PrintEvalErrors=-1, # We are expecting a lot of evaluation errors. -1 switches off printing.
63 PrintLevel=-1,
64)
65
66pdf.plotOn(frame, LineColor="r", Name="noRecovery")
67
68
69# RooFit since ROOT 6.24
70# --------------------------------
71# The minimiser gets information about the "badness" of the violation of the function definition. It uses this
72# to find its way out of the disallowed parameter regions.
73print("\n\n\n-------------- Starting second fit ---------------\n\n")
74
75# Reset the parameters such that the PDF is again undefined.
76a1.setVal(10.0)
77a2.setVal(-1.0)
78
79# Fit again, but pass recovery information to the minimiser:
80fitWithRecovery = pdf.fitTo(
81 data,
82 Save=True,
83 RecoverFromUndefinedRegions=1.0, # The magnitude of the recovery information can be chosen here.
84 # Higher values mean more aggressive recovery.
85 PrintEvalErrors=-1, # We are still expecting a few evaluation errors.
86 PrintLevel=0,
87)
88
89pdf.plotOn(frame, LineColor="b", Name="recovery")
90
91
92# Collect results and plot.
93# --------------------------------
94# We print the two fit results, and plot the fitted curves.
95# The curve of the fit without recovery cannot be plotted, because the PDF is undefined if a2 < 0.
96fitWithoutRecovery.Print()
97print(
98 "Without recovery, the fitter encountered {}".format(fitWithoutRecovery.numInvalidNLL())
99 + " invalid function values. The parameters are unchanged.\n"
100)
101
102fitWithRecovery.Print()
103print(
104 "With recovery, the fitter encountered {}".format(fitWithoutRecovery.numInvalidNLL())
105 + " invalid function values, but the parameters are fitted.\n"
106)
107
108legend = ROOT.TLegend(0.5, 0.7, 0.9, 0.9)
109legend.SetBorderSize(0)
110legend.SetFillStyle(0)
111legend.AddEntry("data", "Data", "P")
112legend.AddEntry("noRecovery", "Without recovery (cannot be plotted)", "L")
113legend.AddEntry("recovery", "With recovery", "L")
114frame.Draw()
115legend.Draw()
116c.Draw()
117
118c.SaveAs("rf612_recoverFromInvalidParameters.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