Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf605_profilell.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Likelihood and minimization: working with the profile likelihood estimator
5///
6/// \macro_image
7/// \macro_code
8/// \macro_output
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooAddPdf.h"
17#include "RooMinimizer.h"
18#include "TCanvas.h"
19#include "TAxis.h"
20#include "RooPlot.h"
21using namespace RooFit;
22
23void rf605_profilell()
24{
25 // C r e a t e m o d e l a n d d a t a s e t
26 // -----------------------------------------------
27
28 // Observable
29 RooRealVar x("x", "x", -20, 20);
30
31 // Model (intentional strong correlations)
32 RooRealVar mean("mean", "mean of g1 and g2", 0, -10, 10);
33 RooRealVar sigma_g1("sigma_g1", "width of g1", 3);
34 RooGaussian g1("g1", "g1", x, mean, sigma_g1);
35
36 RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 6.0);
37 RooGaussian g2("g2", "g2", x, mean, sigma_g2);
38
39 RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
40 RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
41
42 // Generate 1000 events
43 std::unique_ptr<RooDataSet> data{model.generate(x, 1000)};
44
45 // C o n s t r u c t p l a i n l i k e l i h o o d
46 // ---------------------------------------------------
47
48 // Construct unbinned likelihood
49 std::unique_ptr<RooAbsReal> nll{model.createNLL(*data, NumCPU(2))};
50
51 // Minimize likelihood w.r.t all parameters before making plots
52 RooMinimizer(*nll).migrad();
53
54 // Plot likelihood scan frac
55 RooPlot *frame1 = frac.frame(Bins(10), Range(0.01, 0.95), Title("LL and profileLL in frac"));
56 nll->plotOn(frame1, ShiftToZero());
57
58 // Plot likelihood scan in sigma_g2
59 RooPlot *frame2 = sigma_g2.frame(Bins(10), Range(3.3, 5.0), Title("LL and profileLL in sigma_g2"));
60 nll->plotOn(frame2, ShiftToZero());
61
62 // 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
63 // -----------------------------------------------------------------------
64
65 // The profile likelihood estimator on nll for frac will minimize nll w.r.t
66 // all floating parameters except frac for each evaluation
67
68 std::unique_ptr<RooAbsReal> pll_frac{nll->createProfile(frac)};
69
70 // Plot the profile likelihood in frac
71 pll_frac->plotOn(frame1, LineColor(kRed));
72
73 // Adjust frame maximum for visual clarity
74 frame1->SetMinimum(0);
75 frame1->SetMaximum(3);
76
77 // 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
78 // -------------------------------------------------------------------------------
79
80 // The profile likelihood estimator on nll for sigma_g2 will minimize nll
81 // w.r.t all floating parameters except sigma_g2 for each evaluation
82 std::unique_ptr<RooAbsReal> pll_sigmag2{nll->createProfile(sigma_g2)};
83
84 // Plot the profile likelihood in sigma_g2
85 pll_sigmag2->plotOn(frame2, LineColor(kRed));
86
87 // Adjust frame maximum for visual clarity
88 frame2->SetMinimum(0);
89 frame2->SetMaximum(3);
90
91 // Make canvas and draw RooPlots
92 TCanvas *c = new TCanvas("rf605_profilell", "rf605_profilell", 800, 400);
93 c->Divide(2);
94 c->cd(1);
95 gPad->SetLeftMargin(0.15);
96 frame1->GetYaxis()->SetTitleOffset(1.4);
97 frame1->Draw();
98 c->cd(2);
99 gPad->SetLeftMargin(0.15);
100 frame2->GetYaxis()->SetTitleOffset(1.4);
101 frame2->Draw();
102}
#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:45
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:225
TAxis * GetYaxis() const
Definition RooPlot.cxx:1264
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:637
SetMaximum(ymax)
SetMinimum(ymin)
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 NumCPU(Int_t nCPU, Int_t interleave=0)
RooCmdArg ShiftToZero()
RooCmdArg LineColor(Color_t color)
Double_t x[n]
Definition legend1.C:17
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:353
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)