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.
Before ROOT 6.24, RooFit always passed the highest function value that was encountered during the minimisation to the minimiser. If a parameter is far in a disallowed region, the minimiser has to blindly test various values of the parameters. It might find the correct values by chance, but often be unable to recover from bad starting values. Here, we use a model with such bad values.
Starting with ROOT 6.24, the minimiser receives more information. For example, when a PDF is negative, the magnitude of the "undershoot" is passed to the minimiser. The minimiser can use this to compute a gradient, which will eventually lead it out of the disallowed region. The steepness of this gradient can be chosen using RooFit::RecoverFromUndefinedRegions(double). A value of zero is equivalent to RooFit before ROOT 6.24. Positive values activate the recovery. Values between 1. and 10. were found to be a good default. If no argument is passed, RooFit uses 10.
std::unique_ptr<RooDataSet>
data(pdf.generate(
x, 10000));
a1.setVal(10.);
a2.setVal(-1.);
std::unique_ptr<RooFitResult> fitWithoutRecovery{pdf.fitTo(*
data,
RooFit::Save(),
std::cout << "\n\n\n-------------- Starting second fit ---------------\n\n" << std::endl;
a1.setVal(10.);
a2.setVal(-1.);
fitWithoutRecovery->Print();
std::cout << "Without recovery, the fitter encountered " << fitWithoutRecovery->numInvalidNLL()
<< " invalid function values. The parameters are unchanged." << std::endl;
fitWithRecovery->Print();
std::cout << "With recovery, the fitter encountered " << fitWithRecovery->numInvalidNLL()
<< " invalid function values, but the parameters are fitted." << std::endl;
legend->
AddEntry(
"noRecovery",
"Without recovery (cannot be plotted)",
"L");
legend->
AddEntry(
"recovery",
"With recovery",
"L");
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
RooArgSet is a container object that can hold multiple RooAbsArg objects.
static RooMsgService & instance()
Return reference to singleton instance.
StreamConfig & getStream(Int_t id)
A RooPlot is a plot frame and a container for graphics objects within that frame.
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
RooPolynomial implements a polynomial p.d.f of the form.
RooRealVar represents a variable that can be changed from the outside.
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
This class displays a legend box (TPaveText) containing several legend entries.
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
virtual void SetBorderSize(Int_t bordersize=4)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg RecoverFromUndefinedRegions(double strength)
When parameters are chosen such that a PDF is undefined, try to indicate to the minimiser how to leav...
RooCmdArg LineColor(Color_t color)
RooCmdArg Name(const char *name)
void removeTopic(RooFit::MsgTopic oldTopic)
[#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 = 10 +/- 0,a2 = -1 +/- 0,a3 = 0.01)
-------------- Starting second fit ---------------
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
Minuit2Minimizer: Minimize with max-calls 1000 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -863.447294814737688
Edm = 0.000768542764103495092
Nfcn = 318
a1 = -0.498558 +/- 0.0227368 (limited)
a2 = 0.198097 +/- 0.00564227 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 0, estimated distance to minimum: 0
covariance matrix quality: Not calculated at all
Status : MINIMIZE=-1 HESSE=302
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a1 1.0000e+01 +/- 0.00e+00
a2 -1.0000e+00 +/- 0.00e+00
Without recovery, the fitter encountered 23 invalid function values. The parameters are unchanged.
RooFitResult: minimized FCN value: 29650.9, estimated distance to minimum: 0.000767566
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a1 -4.9856e-01 +/- 2.27e-02
a2 1.9810e-01 +/- 5.64e-03
With recovery, the fitter encountered 65 invalid function values, but the parameters are fitted.
- Date
- 10/2020
- Author
- Stephan Hageboeck
Definition in file rf612_recoverFromInvalidParameters.C.