Logo ROOT   6.16/01
Reference Guide
rf605_profilell.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605
5///
6/// Working with the profile likelihood estimator
7///
8/// \macro_image
9/// \macro_output
10/// \macro_code
11/// \author 07/2008 - Wouter Verkerke
12
13
14#include "RooRealVar.h"
15#include "RooDataSet.h"
16#include "RooGaussian.h"
17#include "RooConstVar.h"
18#include "RooAddPdf.h"
19#include "RooMinimizer.h"
20#include "TCanvas.h"
21#include "TAxis.h"
22#include "RooPlot.h"
23using namespace RooFit ;
24
25
26void rf605_profilell()
27{
28 // C r e a t e m o d e l a n d d a t a s e t
29 // -----------------------------------------------
30
31 // Observable
32 RooRealVar x("x","x",-20,20) ;
33
34 // Model (intentional strong correlations)
35 RooRealVar mean("mean","mean of g1 and g2",0,-10,10) ;
36 RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
37 RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
38
39 RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
40 RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
41
42 RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
43 RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
44
45 // Generate 1000 events
46 RooDataSet* data = model.generate(x,1000) ;
47
48
49
50 // C o n s t r u c t p l a i n l i k e l i h o o d
51 // ---------------------------------------------------
52
53 // Construct unbinned likelihood
54 RooAbsReal* nll = model.createNLL(*data,NumCPU(2)) ;
55
56 // Minimize likelihood w.r.t all parameters before making plots
57 RooMinimizer(*nll).migrad() ;
58
59 // Plot likelihood scan frac
60 RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
61 nll->plotOn(frame1,ShiftToZero()) ;
62
63 // Plot likelihood scan in sigma_g2
64 RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
65 nll->plotOn(frame2,ShiftToZero()) ;
66
67
68
69 // 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
70 // -----------------------------------------------------------------------
71
72 // The profile likelihood estimator on nll for frac will minimize nll w.r.t
73 // all floating parameters except frac for each evaluation
74
75 RooAbsReal* pll_frac = nll->createProfile(frac) ;
76
77 // Plot the profile likelihood in frac
78 pll_frac->plotOn(frame1,LineColor(kRed)) ;
79
80 // Adjust frame maximum for visual clarity
81 frame1->SetMinimum(0) ;
82 frame1->SetMaximum(3) ;
83
84
85
86 // 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
87 // -------------------------------------------------------------------------------
88
89 // The profile likelihood estimator on nll for sigma_g2 will minimize nll
90 // w.r.t all floating parameters except sigma_g2 for each evaluation
91 RooAbsReal* pll_sigmag2 = nll->createProfile(sigma_g2) ;
92
93 // Plot the profile likelihood in sigma_g2
94 pll_sigmag2->plotOn(frame2,LineColor(kRed)) ;
95
96 // Adjust frame maximum for visual clarity
97 frame2->SetMinimum(0) ;
98 frame2->SetMaximum(3) ;
99
100
101
102 // Make canvas and draw RooPlots
103 TCanvas *c = new TCanvas("rf605_profilell","rf605_profilell",800, 400);
104 c->Divide(2);
105 c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
106 c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
107
108 delete pll_frac ;
109 delete pll_sigmag2 ;
110 delete nll ;
111}
112
#define c(i)
Definition: RSha256.hxx:101
@ kRed
Definition: Rtypes.h:63
#define gPad
Definition: TVirtualPad.h:286
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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.
Definition: RooAbsReal.cxx:462
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooMinimizer is a wrapper class around ROOT::Fit:Fitter that provides a seamless interface between th...
Definition: RooMinimizer.h:38
Int_t migrad()
Execute MIGRAD.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
virtual void SetMinimum(Double_t minimum=-1111)
Set minimum value of Y axis.
Definition: RooPlot.cxx:958
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void SetMaximum(Double_t maximum=-1111)
Set maximum value of Y axis.
Definition: RooPlot.cxx:948
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
Double_t x[n]
Definition: legend1.C:17
RooCmdArg ShiftToZero()
RooCmdArg NumCPU(Int_t nCPU, Int_t interleave=0)
RooCmdArg LineColor(Color_t color)
RooCmdArg Bins(Int_t nbin)
const char * Title
Definition: TXMLSetup.cxx:67
Ta Range(0, 0, 1, 1)