Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf612_recoverFromInvalidParameters.py File Reference

Namespaces

namespace  rf612_recoverFromInvalidParameters
 

Detailed Description

View in nbviewer Open in SWAN
Likelihood and minimization: Recover from regions where the function is not defined.

We demonstrate improved recovery from disallowed parameters. For this, we use a polynomial PDF of the form

\[ \mathrm{Pol2} = \mathcal{N} \left( c + a_1 \cdot x + a_2 \cdot x^2 + 0.01 \cdot x^3 \right), \]

where \( \mathcal{N} \) is a normalisation factor. Unless the parameters are chosen carefully, this function can be negative, and hence, it cannot be used as a PDF. In this case, RooFit passes an error to the minimiser, which might try to recover.

import ROOT
# Create a fit model:
# The polynomial is notoriously unstable, because it can quickly go negative.
# Since PDFs need to be positive, one often ends up with an unstable fit model.
x = ROOT.RooRealVar("x", "x", -15, 15)
a1 = ROOT.RooRealVar("a1", "a1", -0.5, -10.0, 20.0)
a2 = ROOT.RooRealVar("a2", "a2", 0.2, -10.0, 20.0)
a3 = ROOT.RooRealVar("a3", "a3", 0.01)
pdf = ROOT.RooPolynomial("pol3", "c + a1 * x + a2 * x*x + 0.01 * x*x*x", x, [a1, a2, a3])
# Create toy data with all-positive coefficients:
data = pdf.generate(x, 10000)
# For plotting.
# We create pointers to the plotted objects. We want these objects to leak out of the function,
# so we can still see them after it returns.
c = ROOT.TCanvas()
frame = x.frame()
data.plotOn(frame, Name="data")
# Plotting a PDF with disallowed parameters doesn't work. We would get a lot of error messages.
# Therefore, we disable plotting messages in RooFit's message streams:
ROOT.RooMsgService.instance().getStream(0).removeTopic(ROOT.RooFit.Plotting)
ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.Plotting)
# RooFit before ROOT 6.24
# --------------------------------
# Before 6.24, RooFit wasn't able to recover from invalid parameters. The minimiser just errs around
# the starting values of the parameters without finding any improvement.
# Set up the parameters such that the PDF would come out negative. The PDF is now undefined.
a1.setVal(10.0)
a2.setVal(-1.0)
# Perform a fit:
fitWithoutRecovery = pdf.fitTo(
data,
Save=True,
RecoverFromUndefinedRegions=0.0, # This is how RooFit behaved prior to ROOT 6.24
PrintEvalErrors=-1, # We are expecting a lot of evaluation errors. -1 switches off printing.
PrintLevel=-1,
)
pdf.plotOn(frame, LineColor="r", Name="noRecovery")
# RooFit since ROOT 6.24
# --------------------------------
# The minimiser gets information about the "badness" of the violation of the function definition. It uses this
# to find its way out of the disallowed parameter regions.
print("\n\n\n-------------- Starting second fit ---------------\n\n")
# Reset the parameters such that the PDF is again undefined.
a1.setVal(10.0)
a2.setVal(-1.0)
# Fit again, but pass recovery information to the minimiser:
fitWithRecovery = pdf.fitTo(
data,
Save=True,
RecoverFromUndefinedRegions=1.0, # The magnitude of the recovery information can be chosen here.
# Higher values mean more aggressive recovery.
PrintEvalErrors=-1, # We are still expecting a few evaluation errors.
PrintLevel=0,
)
pdf.plotOn(frame, LineColor="b", Name="recovery")
# Collect results and plot.
# --------------------------------
# We print the two fit results, and plot the fitted curves.
# The curve of the fit without recovery cannot be plotted, because the PDF is undefined if a2 < 0.
fitWithoutRecovery.Print()
print(
"Without recovery, the fitter encountered {}".format(fitWithoutRecovery.numInvalidNLL())
+ " invalid function values. The parameters are unchanged.\n"
)
fitWithRecovery.Print()
print(
"With recovery, the fitter encountered {}".format(fitWithoutRecovery.numInvalidNLL())
+ " invalid function values, but the parameters are fitted.\n"
)
legend = ROOT.TLegend(0.5, 0.7, 0.9, 0.9)
legend.SetBorderSize(0)
legend.SetFillStyle(0)
legend.AddEntry("data", "Data", "P")
legend.AddEntry("noRecovery", "Without recovery (cannot be plotted)", "L")
legend.AddEntry("recovery", "With recovery", "L")
frame.Draw()
legend.Draw()
c.Draw()
c.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
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#0] ERROR:Eval -- RooAbsReal::logEvalError(pol3) evaluation error,
origin : RooPolynomial::pol3[ x=x coefList=(a1,a2,a3) ]
message : p.d.f normalization integral is zero or negative: -2220.000000
server values: x=x=0, coefList=(a1 = 2.60781 +/- 11.9002,a2 = -1 +/- 11.5683,a3 = 0.01)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
**********
** 15 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 a1 1.00000e+01 1.19002e+01 -1.00000e+01 2.00000e+01
2 a2 -1.00000e+00 1.15683e+01 -1.00000e+01 2.00000e+01
**********
** 16 **SET ERR 0.5
**********
**********
** 17 **SET PRINT 0
**********
**********
** 18 **SET STR 1
**********
**********
** 19 **MIGRAD 1000 1
**********
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
FCN=-858.564 FROM MIGRAD STATUS=CONVERGED 243 CALLS 244 TOTAL
EDM=7.33131e-05 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a1 -4.98209e-01 2.27025e-02 2.77075e-05 -2.18163e+00
2 a2 1.98271e-01 5.64128e-03 1.48249e-05 -2.54212e+01
ERR DEF= 0.5
**********
** 20 **SET ERR 0.5
**********
**********
** 21 **SET PRINT 0
**********
**********
** 22 **HESSE 1000
**********
FCN=-858.564 FROM HESSE STATUS=OK 10 CALLS 254 TOTAL
EDM=7.3377e-05 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 a1 -4.98209e-01 2.27254e-02 8.14968e-06 -3.41822e+01
2 a2 1.98271e-01 5.64697e-03 7.41245e-06 3.10901e+01
ERR DEF= 0.5
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 0, estimated distance to minimum: 0
covariance matrix quality: Approximation only, not accurate
Status : MINIMIZE=200 HESSE=200
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a1 2.6078e+00 +/- 1.19e+01
a2 -1.0000e+00 +/- 1.16e+01
RooFitResult: minimized FCN value: 29650.9, estimated distance to minimum: 7.3377e-05
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a1 -4.9821e-01 +/- 2.27e-02
a2 1.9827e-01 +/- 5.65e-03
-------------- Starting second fit ---------------
Without recovery, the fitter encountered 66 invalid function values. The parameters are unchanged.
With recovery, the fitter encountered 66 invalid function values, but the parameters are fitted.
Date
June 2021
Author
Harshal Shende, Stephan Hageboeck (C++ version)

Definition in file rf612_recoverFromInvalidParameters.py.