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_output
8/// \macro_code
9///
10/// \date July 2008
11/// \author Wouter Verkerke
12
13#include "RooRealVar.h"
14#include "RooDataSet.h"
15#include "RooGaussian.h"
16#include "RooConstVar.h"
17#include "RooAddPdf.h"
18#include "RooMinimizer.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
24void rf605_profilell()
25{
26 // C r e a t e m o d e l a n d d a t a s e t
27 // -----------------------------------------------
28
29 // Observable
30 RooRealVar x("x", "x", -20, 20);
31
32 // Model (intentional strong correlations)
33 RooRealVar mean("mean", "mean of g1 and g2", 0, -10, 10);
34 RooRealVar sigma_g1("sigma_g1", "width of g1", 3);
35 RooGaussian g1("g1", "g1", x, mean, sigma_g1);
36
37 RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 6.0);
38 RooGaussian g2("g2", "g2", x, mean, sigma_g2);
39
40 RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
41 RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
42
43 // Generate 1000 events
44 RooDataSet *data = model.generate(x, 1000);
45
46 // C o n s t r u c t p l a i n l i k e l i h o o d
47 // ---------------------------------------------------
48
49 // Construct unbinned likelihood
50 RooAbsReal *nll = model.createNLL(*data, NumCPU(2));
51
52 // Minimize likelihood w.r.t all parameters before making plots
53 RooMinimizer(*nll).migrad();
54
55 // Plot likelihood scan frac
56 RooPlot *frame1 = frac.frame(Bins(10), Range(0.01, 0.95), Title("LL and profileLL in frac"));
57 nll->plotOn(frame1, ShiftToZero());
58
59 // Plot likelihood scan in sigma_g2
60 RooPlot *frame2 = sigma_g2.frame(Bins(10), Range(3.3, 5.0), Title("LL and profileLL in sigma_g2"));
61 nll->plotOn(frame2, ShiftToZero());
62
63 // 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
64 // -----------------------------------------------------------------------
65
66 // The profile likelihood estimator on nll for frac will minimize nll w.r.t
67 // all floating parameters except frac for each evaluation
68
69 RooAbsReal *pll_frac = nll->createProfile(frac);
70
71 // Plot the profile likelihood in frac
72 pll_frac->plotOn(frame1, LineColor(kRed));
73
74 // Adjust frame maximum for visual clarity
75 frame1->SetMinimum(0);
76 frame1->SetMaximum(3);
77
78 // 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
79 // -------------------------------------------------------------------------------
80
81 // The profile likelihood estimator on nll for sigma_g2 will minimize nll
82 // w.r.t all floating parameters except sigma_g2 for each evaluation
83 RooAbsReal *pll_sigmag2 = nll->createProfile(sigma_g2);
84
85 // Plot the profile likelihood in sigma_g2
86 pll_sigmag2->plotOn(frame2, LineColor(kRed));
87
88 // Adjust frame maximum for visual clarity
89 frame2->SetMinimum(0);
90 frame2->SetMaximum(3);
91
92 // Make canvas and draw RooPlots
93 TCanvas *c = new TCanvas("rf605_profilell", "rf605_profilell", 800, 400);
94 c->Divide(2);
95 c->cd(1);
96 gPad->SetLeftMargin(0.15);
97 frame1->GetYaxis()->SetTitleOffset(1.4);
98 frame1->Draw();
99 c->cd(2);
100 gPad->SetLeftMargin(0.15);
101 frame2->GetYaxis()->SetTitleOffset(1.4);
102 frame2->Draw();
103
104 delete pll_frac;
105 delete pll_sigmag2;
106 delete nll;
107}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
#define gPad
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg(), const RooCmdArg &arg10=RooCmdArg()) const
Plot (project) PDF on specified frame.
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function.
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Int_t migrad()
Execute MIGRAD.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
virtual void SetMinimum(Double_t minimum=-1111)
Set minimum value of Y axis.
Definition RooPlot.cxx:1091
TAxis * GetYaxis() const
Definition RooPlot.cxx:1263
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.
Definition RooPlot.cxx:249
virtual void SetMaximum(Double_t maximum=-1111)
Set maximum value of Y axis.
Definition RooPlot.cxx:1081
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:293
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
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition TXMLSetup.cxx:68
Ta Range(0, 0, 1, 1)