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.
␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- RooMinimizer::optimizeConst: 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)
-------------- Starting second fit ---------------
[#1] INFO:Minization -- RooMinimizer::optimizeConst: 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:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 0, estimated distance to minimum: 0
covariance matrix quality: Approximation only, not accurate
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a1 2.6078e+00 +/- 1.19e+01
a2 -1.0000e+00 +/- 1.16e+01
Without recovery, the fitter encountered 66 invalid function values. The parameters are unchanged.
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
With recovery, the fitter encountered 72 invalid function values, but the parameters are fitted.
void rf612_recoverFromInvalidParameters() {
std::unique_ptr<RooDataSet> data(pdf.generate(
x, 10000));
a1.setVal(10.);
a2.setVal(-1.);
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(frame->
findObject(
"noRecovery"),
"Without recovery (cannot be plotted)",
"L");
}
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Int_t numInvalidNLL() const
Return number of NLL evaluations with problems.
RooPlot * plotOn(RooPlot *frame, const RooAbsArg &par1, const RooAbsArg &par2, const char *options="ME") const
Add objects to a 2D plot.
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_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
TObject * findObject(const char *name, const TClass *clas=0) const
Find the named object in our list of items and return a pointer to it.
virtual void Draw(Option_t *options=0)
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.
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
virtual void SetBorderSize(Int_t bordersize=4)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg Save(Bool_t flag=kTRUE)
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)