ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf601_intminuit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #601
4 ///
5 /// Interactive minimization with MINUIT
6 ///
7 ///
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 
12 #ifndef __CINT__
13 #include "RooGlobalFunc.h"
14 #endif
15 #include "RooRealVar.h"
16 #include "RooDataSet.h"
17 #include "RooGaussian.h"
18 #include "RooConstVar.h"
19 #include "RooProdPdf.h"
20 #include "RooAddPdf.h"
21 #include "RooMinuit.h"
22 #include "RooFitResult.h"
23 #include "RooPlot.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "TH1.h"
27 using namespace RooFit ;
28 
29 
30 void rf601_intminuit()
31 {
32  // S e t u p p d f a n d l i k e l i h o o d
33  // -----------------------------------------------
34 
35  // Observable
36  RooRealVar x("x","x",-20,20) ;
37 
38  // Model (intentional strong correlations)
39  RooRealVar mean("mean","mean of g1 and g2",0) ;
40  RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
41  RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
42 
43  RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
44  RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
45 
46  RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
47  RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
48 
49  // Generate 1000 events
50  RooDataSet* data = model.generate(x,1000) ;
51 
52  // Construct unbinned likelihood of model w.r.t. data
53  RooAbsReal* nll = model.createNLL(*data) ;
54 
55  // I n t e r a c t i v e m i n i m i z a t i o n , e r r o r a n a l y s i s
56  // -------------------------------------------------------------------------------
57 
58  // Create MINUIT interface object
59  RooMinuit m(*nll) ;
60 
61  // Activate verbose logging of MINUIT parameter space stepping
62  m.setVerbose(kTRUE) ;
63 
64  // Call MIGRAD to minimize the likelihood
65  m.migrad() ;
66 
67  // Print values of all parameters, that reflect values (and error estimates)
68  // that are back propagated from MINUIT
69  model.getParameters(x)->Print("s") ;
70 
71  // Disable verbose logging
72  m.setVerbose(kFALSE) ;
73 
74  // Run HESSE to calculate errors from d2L/dp2
75  m.hesse() ;
76 
77  // Print value (and error) of sigma_g2 parameter, that reflects
78  // value and error back propagated from MINUIT
79  sigma_g2.Print() ;
80 
81  // Run MINOS on sigma_g2 parameter only
82  m.minos(sigma_g2) ;
83 
84  // Print value (and error) of sigma_g2 parameter, that reflects
85  // value and error back propagated from MINUIT
86  sigma_g2.Print() ;
87 
88 
89 
90  // S a v i n g r e s u l t s , c o n t o u r p l o t s
91  // ---------------------------------------------------------
92 
93  // Save a snapshot of the fit result. This object contains the initial
94  // fit parameters, the final fit parameters, the complete correlation
95  // matrix, the EDM, the minimized FCN , the last MINUIT status code and
96  // the number of times the RooFit function object has indicated evaluation
97  // problems (e.g. zero probabilities during likelihood evaluation)
98  RooFitResult* r = m.save() ;
99 
100  // Make contour plot of mx vs sx at 1,2,3 sigma
101  RooPlot* frame = m.contour(frac,sigma_g2,1,2,3) ;
102  frame->SetTitle("RooMinuit contour plot") ;
103 
104  // Print the fit result snapshot
105  r->Print("v") ;
106 
107 
108 
109  // C h a n g e p a r a m e t e r v a l u e s , f l o a t i n g
110  // -----------------------------------------------------------------
111 
112  // At any moment you can manually change the value of a (constant)
113  // parameter
114  mean = 0.3 ;
115 
116  // Rerun MIGRAD,HESSE
117  m.migrad() ;
118  m.hesse() ;
119  frac.Print() ;
120 
121  // Now fix sigma_g2
122  sigma_g2.setConstant(kTRUE) ;
123 
124  // Rerun MIGRAD,HESSE
125  m.migrad() ;
126  m.hesse() ;
127  frac.Print() ;
128 
129 
130 
131  new TCanvas("rf601_intminuit","rf601_intminuit",600,600) ;
132  gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
133 
134 }
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:65
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
tuple g2
Definition: multifit.py:39
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
const Bool_t kFALSE
Definition: Rtypes.h:92
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1099
tuple g1
Definition: multifit.py:38
Double_t x[n]
Definition: legend1.C:17
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
TMarker * m
Definition: textangle.C:8
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
virtual void Print(Option_t *option="") const
Dump this marker with its attributes.
Definition: TMarker.cxx:277
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
#define gPad
Definition: TVirtualPad.h:288
const Bool_t kTRUE
Definition: Rtypes.h:91
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
RooMinuit is a wrapper class around TFitter/TMinuit that provides a seamless interface between the MI...
Definition: RooMinuit.h:39