ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf605_profilell.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #605
4 ///
5 /// Working with the profile likelihood estimator
6 ///
7 ///
8 ///
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooConstVar.h"
20 #include "RooAddPdf.h"
21 #include "RooMinuit.h"
22 #include "TCanvas.h"
23 #include "TAxis.h"
24 #include "RooPlot.h"
25 using namespace RooFit ;
26 
27 
28 void rf605_profilell()
29 {
30  // C r e a t e m o d e l a n d d a t a s e t
31  // -----------------------------------------------
32 
33  // Observable
34  RooRealVar x("x","x",-20,20) ;
35 
36  // Model (intentional strong correlations)
37  RooRealVar mean("mean","mean of g1 and g2",0,-10,10) ;
38  RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
39  RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
40 
41  RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
42  RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
43 
44  RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
45  RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
46 
47  // Generate 1000 events
48  RooDataSet* data = model.generate(x,1000) ;
49 
50 
51 
52  // C o n s t r u c t p l a i n l i k e l i h o o d
53  // ---------------------------------------------------
54 
55  // Construct unbinned likelihood
56  RooAbsReal* nll = model.createNLL(*data,NumCPU(2)) ;
57 
58  // Minimize likelihood w.r.t all parameters before making plots
59  RooMinuit(*nll).migrad() ;
60 
61  // Plot likelihood scan frac
62  RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
63  nll->plotOn(frame1,ShiftToZero()) ;
64 
65  // Plot likelihood scan in sigma_g2
66  RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
67  nll->plotOn(frame2,ShiftToZero()) ;
68 
69 
70 
71  // 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
72  // -----------------------------------------------------------------------
73 
74  // The profile likelihood estimator on nll for frac will minimize nll w.r.t
75  // all floating parameters except frac for each evaluation
76 
77  RooAbsReal* pll_frac = nll->createProfile(frac) ;
78 
79  // Plot the profile likelihood in frac
80  pll_frac->plotOn(frame1,LineColor(kRed)) ;
81 
82  // Adjust frame maximum for visual clarity
83  frame1->SetMinimum(0) ;
84  frame1->SetMaximum(3) ;
85 
86 
87 
88  // 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
89  // -------------------------------------------------------------------------------
90 
91  // The profile likelihood estimator on nll for sigma_g2 will minimize nll
92  // w.r.t all floating parameters except sigma_g2 for each evaluation
93  RooAbsReal* pll_sigmag2 = nll->createProfile(sigma_g2) ;
94 
95  // Plot the profile likelihood in sigma_g2
96  pll_sigmag2->plotOn(frame2,LineColor(kRed)) ;
97 
98  // Adjust frame maximum for visual clarity
99  frame2->SetMinimum(0) ;
100  frame2->SetMaximum(3) ;
101 
102 
103 
104  // Make canvas and draw RooPlots
105  TCanvas *c = new TCanvas("rf605_profilell","rf605_profilell",800, 400);
106  c->Divide(2);
107  c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
108  c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
109 
110  delete pll_frac ;
111  delete pll_sigmag2 ;
112  delete nll ;
113 }
114 
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:245
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.
friend class RooMinuit
Definition: RooAbsArg.h:302
RooCmdArg LineColor(Color_t color)
tuple g2
Definition: multifit.py:39
Definition: Rtypes.h:61
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
tuple g1
Definition: multifit.py:38
RooCmdArg Title(const char *name)
RooCmdArg Range(const char *rangeName, Bool_t adjustNorm=kTRUE)
virtual void SetMinimum(Double_t minimum=-1111)
Set minimum value of Y axis.
Definition: RooPlot.cxx:959
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooCmdArg NumCPU(Int_t nCPU, Int_t interleave=0)
RooCmdArg ShiftToZero()
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual void SetMaximum(Double_t maximum=-1111)
Set maximum value of Y axis.
Definition: RooPlot.cxx:949
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
virtual RooAbsReal * createProfile(const RooArgSet &paramsOfInterest)
Create a RooProfileLL object that eliminates all nuisance parameters in the present function...
Definition: RooAbsReal.cxx:463
RooCmdArg Bins(Int_t nbin)
#define gPad
Definition: TVirtualPad.h:288
RooAbsMoment * mean(RooRealVar &obs)
Definition: RooAbsReal.h:297
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559