Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf204a_extendedLikelihood.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Extended maximum likelihood fit in multiple ranges.

When an extended pdf and multiple ranges are used, the RooExtendPdf cannot correctly interpret the coefficients used for extension. This can be solved by using a RooAddPdf for extending the model.

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooChebychev.h"
#include "RooAddPdf.h"
#include "RooExtendPdf.h"
#include "RooFitResult.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;
{
// S e t u p c o m p o n e n t p d f s
// ---------------------------------------
// Declare observable x
RooRealVar x("x","x",0,11) ;
// Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
RooRealVar mean("mean","mean of gaussians",5) ;
RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
RooRealVar sigma2("sigma2","width of gaussians",1) ;
RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
// Build Chebychev polynomial pdf
RooRealVar a0("a0","a0",0.5,0.,1.) ;
RooRealVar a1("a1","a1",0.2,0.,1.) ;
RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
// Sum the signal components into a composite signal pdf
RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
// E x t e n d t h e p d f s
// -----------------------------
// Define signal range in which events counts are to be defined
x.setRange("signalRange",4,6) ;
// Associated nsig/nbkg as expected number of events with sig/bkg _in_the_range_ "signalRange"
RooRealVar nsig("nsig","number of signal events in signalRange",500,0.,10000) ;
RooRealVar nbkg("nbkg","number of background events in signalRange",500,0,10000) ;
// Use AddPdf to extend the model. Giving as many coefficients as pdfs switches
// on extension.
RooAddPdf model("model","(g1+g2)+a", RooArgList(bkg,sig), RooArgList(nbkg,nsig)) ;
// S a m p l e d a t a , f i t m o d e l
// -------------------------------------------
// Generate 1000 events from model so that nsig,nbkg come out to numbers <<500 in fit
std::unique_ptr<RooDataSet> data{model.generate(x,1000)};
auto canv = new TCanvas("Canvas", "Canvas", 1500, 600);
canv->Divide(3,1);
// Fit full range
// -------------------------------------------
canv->cd(1);
// Perform unbinned ML fit to data, full range
// IMPORTANT:
// The model needs to be copied when fitting with different ranges because
// the interpretation of the coefficients is tied to the fit range
// that's used in the first fit
RooAddPdf model1(model);
std::unique_ptr<RooFitResult> r{model1.fitTo(*data,Save(), PrintLevel(-1))};
r->Print() ;
RooPlot * frame = x.frame(Title("Full range fitted"));
data->plotOn(frame);
model1.plotOn(frame, VisualizeError(*r));
model1.plotOn(frame);
model1.paramOn(frame);
frame->Draw();
// Fit in two regions
// -------------------------------------------
canv->cd(2);
x.setRange("left", 0., 4.);
x.setRange("right", 6., 10.);
RooAddPdf model2(model);
std::unique_ptr<RooFitResult> r2{model2.fitTo(*data, Range("left,right"), Save(), PrintLevel(-1))};
r2->Print();
RooPlot * frame2 = x.frame(Title("Fit in left/right sideband"));
data->plotOn(frame2);
model2.plotOn(frame2, VisualizeError(*r2));
model2.plotOn(frame2);
model2.paramOn(frame2);
frame2->Draw();
// Fit in one region
// -------------------------------------------
// Note how restricting the region to only the left tail increases
// the fit uncertainty
canv->cd(3);
x.setRange("leftToMiddle", 0., 5.);
RooAddPdf model3(model);
std::unique_ptr<RooFitResult> r3{model3.fitTo(*data, Range("leftToMiddle"), Save(), PrintLevel(-1))};
r3->Print();
RooPlot * frame3 = x.frame(Title("Fit from left to middle"));
data->plotOn(frame3);
model3.plotOn(frame3, VisualizeError(*r3));
model3.plotOn(frame3);
model3.paramOn(frame3);
frame3->Draw();
canv->Draw();
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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 r
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Chebychev polynomial p.d.f.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:237
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg VisualizeError(const RooDataSet &paramData, double Z=1)
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68
Ta Range(0, 0, 1, 1)
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-inf, inf] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-inf, inf] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'signalRange' created with bounds [4,6]
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: -3872.49, estimated distance to minimum: 4.30406e-05
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a0 4.2647e-01 +/- 7.59e-02
a1 1.7594e-01 +/- 1.10e-01
nbkg 5.1101e+02 +/- 3.60e+01
nsig 4.8899e+02 +/- 3.57e+01
sig1frac 8.6392e-01 +/- 1.08e-01
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'left' created with bounds [0,4]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'right' created with bounds [6,10]
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_left' created with bounds [0,4]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_right' created with bounds [6,10]
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: -1134.15, estimated distance to minimum: 3.61209e-05
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a0 3.2415e-01 +/- 1.09e-01
a1 3.0373e-02 +/- 2.12e-01
nbkg 5.0182e+02 +/- 3.94e+01
nsig 4.1091e+02 +/- 2.63e+02
sig1frac 8.5838e-01 +/- 2.74e-01
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData_left,fit_nll_model_modelData_right'
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'leftToMiddle' created with bounds [0,5]
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData' created with bounds [0,5]
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: -1620.17, estimated distance to minimum: 0.000658484
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
a0 7.0685e-01 +/- 7.00e-01
a1 2.2397e-01 +/- 6.03e-01
nbkg 7.1175e+02 +/- 1.22e+03
nsig 4.4359e+02 +/- 1.28e+02
sig1frac 9.7243e-01 +/- 8.92e-01
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f was fitted in a subrange and no explicit Range() and NormRange() was specified. Plotting / normalising in fit range. To override, do one of the following
- Clear the automatic fit range attribute: <pdf>.removeStringAttribute("fitrange");
- Explicitly specify the plotting range: Range("<rangeName>").
- Explicitly specify where to compute the normalisation: NormRange("<rangeName>").
The default (full) range can be denoted with Range("") / NormRange("").
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) only plotting range 'fit_nll_model_modelData'
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) p.d.f. curve is normalized using explicit choice of ranges 'fit_nll_model_modelData'
Date
December 2018
Author
Stephan Hageboeck

Definition in file rf204a_extendedLikelihood.C.