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

Detailed Description

View in nbviewer Open in SWAN
Likelihood and minimization: working with the profile likelihood estimator

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooAddPdf.h"
#include "RooMinimizer.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;
{
// C r e a t e m o d e l a n d d a t a s e t
// -----------------------------------------------
// Observable
RooRealVar x("x", "x", -20, 20);
// Model (intentional strong correlations)
RooRealVar mean("mean", "mean of g1 and g2", 0, -10, 10);
RooRealVar sigma_g1("sigma_g1", "width of g1", 3);
RooGaussian g1("g1", "g1", x, mean, sigma_g1);
RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 6.0);
RooGaussian g2("g2", "g2", x, mean, sigma_g2);
RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
// Generate 1000 events
std::unique_ptr<RooDataSet> data{model.generate(x, 1000)};
// C o n s t r u c t p l a i n l i k e l i h o o d
// ---------------------------------------------------
// Construct unbinned likelihood
std::unique_ptr<RooAbsReal> nll{model.createNLL(*data, NumCPU(2))};
// Minimize likelihood w.r.t all parameters before making plots
// Plot likelihood scan frac
RooPlot *frame1 = frac.frame(Bins(10), Range(0.01, 0.95), Title("LL and profileLL in frac"));
nll->plotOn(frame1, ShiftToZero());
// Plot likelihood scan in sigma_g2
RooPlot *frame2 = sigma_g2.frame(Bins(10), Range(3.3, 5.0), Title("LL and profileLL in sigma_g2"));
nll->plotOn(frame2, ShiftToZero());
// C o n s t r u c t p r o f i l e l i k e l i h o o d i n f r a c
// -----------------------------------------------------------------------
// The profile likelihood estimator on nll for frac will minimize nll w.r.t
// all floating parameters except frac for each evaluation
std::unique_ptr<RooAbsReal> pll_frac{nll->createProfile(frac)};
// Plot the profile likelihood in frac
pll_frac->plotOn(frame1, LineColor(kRed));
// Adjust frame maximum for visual clarity
frame1->SetMinimum(0);
frame1->SetMaximum(3);
// C o n s t r u c t p r o f i l e l i k e l i h o o d i n s i g m a _ g 2
// -------------------------------------------------------------------------------
// The profile likelihood estimator on nll for sigma_g2 will minimize nll
// w.r.t all floating parameters except sigma_g2 for each evaluation
std::unique_ptr<RooAbsReal> pll_sigmag2{nll->createProfile(sigma_g2)};
// Plot the profile likelihood in sigma_g2
pll_sigmag2->plotOn(frame2, LineColor(kRed));
// Adjust frame maximum for visual clarity
frame2->SetMinimum(0);
frame2->SetMaximum(3);
// Make canvas and draw RooPlots
TCanvas *c = new TCanvas("rf605_profilell", "rf605_profilell", 800, 400);
c->Divide(2);
c->cd(1);
gPad->SetLeftMargin(0.15);
frame1->GetYaxis()->SetTitleOffset(1.4);
frame1->Draw();
c->cd(2);
gPad->SetLeftMargin(0.15);
frame2->GetYaxis()->SetTitleOffset(1.4);
frame2->Draw();
}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
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
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Wrapper class around ROOT::Fit:Fitter that provides a seamless interface between the minimizer functi...
int migrad()
Execute MIGRAD.
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
virtual void SetMinimum(double minimum=-1111)
Set minimum value of Y axis.
Definition RooPlot.cxx:1059
virtual void SetMaximum(double maximum=-1111)
Set maximum value of Y axis.
Definition RooPlot.cxx:1049
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
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
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg ShiftToZero()
RooCmdArg LineColor(Color_t color)
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 'sigma_g1' with range [-inf, inf] of the RooGaussian 'g1' exceeds the safe range of (0, inf). Advise to limit its range.
[#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
Minuit2Minimizer: Minimize with max-calls 1500 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = 2659.73712858695399
Edm = 0.000190395763129910388
Nfcn = 60
frac = 0.62118 +/- 0.165788 (limited)
mean = 0.00442366 +/- 0.109372 (limited)
sigma_g2 = 4.10789 +/- 0.405468 (limited)
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[frac]) Creating instance of MINUIT
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[frac]) determining minimum likelihood for current configurations w.r.t all observable
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[frac]) minimum found at (frac=0.62104)
..................................................................................
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[sigma_g2]) Creating instance of MINUIT
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[sigma_g2]) determining minimum likelihood for current configurations w.r.t all observable
[#1] INFO:Minimization -- RooProfileLL::evaluate(RooEvaluatorWrapper_Profile[sigma_g2]) minimum found at (sigma_g2=4.11258)
....................................................................................
Date
July 2008
Author
Wouter Verkerke

Definition in file rf605_profilell.C.