Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf601_intminuit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Likelihood and minimization: interactive minimization with MINUIT
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 "RooProdPdf.h"
18#include "RooAddPdf.h"
19#include "RooMinimizer.h"
20#include "RooFitResult.h"
21#include "RooPlot.h"
22#include "TCanvas.h"
23#include "TAxis.h"
24#include "TH1.h"
25using namespace RooFit;
26
27void rf601_intminuit()
28{
29 // S e t u p p d f a n d l i k e l i h o o d
30 // -----------------------------------------------
31
32 // Observable
33 RooRealVar x("x", "x", -20, 20);
34
35 // Model (intentional strong correlations)
36 RooRealVar mean("mean", "mean of g1 and g2", 0);
37 RooRealVar sigma_g1("sigma_g1", "width of g1", 3);
38 RooGaussian g1("g1", "g1", x, mean, sigma_g1);
39
40 RooRealVar sigma_g2("sigma_g2", "width of g2", 4, 3.0, 6.0);
41 RooGaussian g2("g2", "g2", x, mean, sigma_g2);
42
43 RooRealVar frac("frac", "frac", 0.5, 0.0, 1.0);
44 RooAddPdf model("model", "model", RooArgList(g1, g2), frac);
45
46 // Generate 1000 events
47 RooDataSet *data = model.generate(x, 1000);
48
49 // Construct unbinned likelihood of model w.r.t. data
50 RooAbsReal *nll = model.createNLL(*data);
51
52 // 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
53 // -------------------------------------------------------------------------------
54
55 // Create MINUIT interface object
56 RooMinimizer m(*nll);
57
58 // Activate verbose logging of MINUIT parameter space stepping
59 m.setVerbose(kTRUE);
60
61 // Call MIGRAD to minimize the likelihood
62 m.migrad();
63
64 // Print values of all parameters, that reflect values (and error estimates)
65 // that are back propagated from MINUIT
66 model.getParameters(x)->Print("s");
67
68 // Disable verbose logging
69 m.setVerbose(kFALSE);
70
71 // Run HESSE to calculate errors from d2L/dp2
72 m.hesse();
73
74 // Print value (and error) of sigma_g2 parameter, that reflects
75 // value and error back propagated from MINUIT
76 sigma_g2.Print();
77
78 // Run MINOS on sigma_g2 parameter only
79 m.minos(sigma_g2);
80
81 // Print value (and error) of sigma_g2 parameter, that reflects
82 // value and error back propagated from MINUIT
83 sigma_g2.Print();
84
85 // 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
86 // ---------------------------------------------------------
87
88 // Save a snapshot of the fit result. This object contains the initial
89 // fit parameters, the final fit parameters, the complete correlation
90 // matrix, the EDM, the minimized FCN , the last MINUIT status code and
91 // the number of times the RooFit function object has indicated evaluation
92 // problems (e.g. zero probabilities during likelihood evaluation)
93 RooFitResult *r = m.save();
94
95 // Make contour plot of mx vs sx at 1,2,3 sigma
96 RooPlot *frame = m.contour(frac, sigma_g2, 1, 2, 3);
97 frame->SetTitle("Minuit contour plot");
98
99 // Print the fit result snapshot
100 r->Print("v");
101
102 // 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
103 // -----------------------------------------------------------------
104
105 // At any moment you can manually change the value of a (constant)
106 // parameter
107 mean = 0.3;
108
109 // Rerun MIGRAD,HESSE
110 m.migrad();
111 m.hesse();
112 frac.Print();
113
114 // Now fix sigma_g2
115 sigma_g2.setConstant(kTRUE);
116
117 // Rerun MIGRAD,HESSE
118 m.migrad();
119 m.hesse();
120 frac.Print();
121
122 new TCanvas("rf601_intminuit", "rf601_intminuit", 600, 600);
123 gPad->SetLeftMargin(0.15);
124 frame->GetYaxis()->SetTitleOffset(1.4);
125 frame->Draw();
126}
ROOT::R::TRInterface & r
Definition Object.C:4
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#define gPad
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
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
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
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...
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition RooPlot.cxx:1242
TAxis * GetYaxis() const
Definition RooPlot.cxx:1263
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
virtual void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition TObject.cxx:552
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
auto * m
Definition textangle.C:8