Logo ROOT   6.16/01
Reference Guide
rf510_wsnamedsets.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #510
5///
6/// Working with named parameter sets and parameter snapshots in
7/// workspaces
8///
9/// \macro_image
10/// \macro_output
11/// \macro_code
12/// \author 04/2009 - Wouter Verkerke
13
14
15#include "RooRealVar.h"
16#include "RooDataSet.h"
17#include "RooGaussian.h"
18#include "RooConstVar.h"
19#include "RooChebychev.h"
20#include "RooAddPdf.h"
21#include "RooWorkspace.h"
22#include "RooPlot.h"
23#include "TCanvas.h"
24#include "TAxis.h"
25#include "TFile.h"
26#include "TH1.h"
27
28using namespace RooFit;
29
30
31void fillWorkspace(RooWorkspace& w) ;
32
33void rf510_wsnamedsets()
34{
35 // C r e a t e m o d e l a n d d a t a s e t
36 // -----------------------------------------------
37
38 RooWorkspace* w = new RooWorkspace("w") ;
39 fillWorkspace(*w) ;
40
41 // Exploit convention encoded in named set "parameters" and "observables"
42 // to use workspace contents w/o need for introspected
43 RooAbsPdf* model = w->pdf("model") ;
44
45 // Generate data from p.d.f. in given observables
46 RooDataSet* data = model->generate(*w->set("observables"),1000) ;
47
48 // Fit model to data
49 model->fitTo(*data) ;
50
51 // Plot fitted model and data on frame of first (only) observable
52 RooPlot* frame = ((RooRealVar*)w->set("observables")->first())->frame() ;
53 data->plotOn(frame) ;
54 model->plotOn(frame) ;
55
56 // Overlay plot with model with reference parameters as stored in snapshots
57 w->loadSnapshot("reference_fit") ;
58 model->plotOn(frame,LineColor(kRed)) ;
59 w->loadSnapshot("reference_fit_bkgonly") ;
60 model->plotOn(frame,LineColor(kRed),LineStyle(kDashed)) ;
61
62
63 // Draw the frame on the canvas
64 new TCanvas("rf510_wsnamedsets","rf503_wsnamedsets",600,600) ;
65 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
66
67
68 // Print workspace contents
69 w->Print() ;
70
71
72 // Workspace will remain in memory after macro finishes
73 gDirectory->Add(w) ;
74
75}
76
77
78
79void fillWorkspace(RooWorkspace& w)
80{
81 // C r e a t e m o d e l
82 // -----------------------
83
84 // Declare observable x
85 RooRealVar x("x","x",0,10) ;
86
87 // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
88 RooRealVar mean("mean","mean of gaussians",5,0,10) ;
89 RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
90 RooRealVar sigma2("sigma2","width of gaussians",1) ;
91
92 RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
93 RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
94
95 // Build Chebychev polynomial p.d.f.
96 RooRealVar a0("a0","a0",0.5,0.,1.) ;
97 RooRealVar a1("a1","a1",0.2,0.,1.) ;
98 RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
99
100 // Sum the signal components into a composite signal p.d.f.
101 RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
102 RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
103
104 // Sum the composite signal and background
105 RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
106 RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
107
108 // Import model into p.d.f.
109 w.import(model) ;
110
111
112 // E n c o d e d e f i n i t i o n o f p a r a m e t e r s i n w o r k s p a c e
113 // ---------------------------------------------------------------------------------------
114
115
116 // Define named sets "parameters" and "observables", which list which variables should be considered
117 // parameters and observables by the users convention
118 //
119 // Variables appearing in sets _must_ live in the workspace already, or the autoImport flag
120 // of defineSet must be set to import them on the fly. Named sets contain only references
121 // to the original variables, therefore the value of observables in named sets already
122 // reflect their 'current' value
123 RooArgSet* params = (RooArgSet*) model.getParameters(x) ;
124 w.defineSet("parameters",*params) ;
125 w.defineSet("observables",x) ;
126
127
128 // E n c o d e r e f e r e n c e v a l u e f o r p a r a m e t e r s i n w o r k s p a c e
129 // ---------------------------------------------------------------------------------------------------
130
131
132 // Define a parameter 'snapshot' in the p.d.f.
133 // Unlike a named set, a parameter snapshot stores an independent set of values for
134 // a given set of variables in the workspace. The values can be stored and reloaded
135 // into the workspace variable objects using the loadSnapshot() and saveSnapshot()
136 // methods. A snapshot saves the value of each variable, any errors that are stored
137 // with it as well as the 'Constant' flag that is used in fits to determine if a
138 // parameter is kept fixed or not.
139
140 // Do a dummy fit to a (supposedly) reference dataset here and store the results
141 // of that fit into a snapshot
142 RooDataSet* refData = model.generate(x,10000) ;
143 model.fitTo(*refData,PrintLevel(-1)) ;
144
145 // The kTRUE flag imports the values of the objects in (*params) into the workspace
146 // If not set, the present values of the workspace parameters objects are stored
147 w.saveSnapshot("reference_fit",*params,kTRUE) ;
148
149 // Make another fit with the signal component forced to zero
150 // and save those parameters too
151
152 bkgfrac.setVal(1) ;
153 bkgfrac.setConstant(kTRUE) ;
154 bkgfrac.removeError() ;
155 model.fitTo(*refData,PrintLevel(-1)) ;
156
157 w.saveSnapshot("reference_fit_bkgonly",*params,kTRUE) ;
158
159
160}
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:63
@ kDashed
Definition: TAttLine.h:48
#define gDirectory
Definition: TDirectory.h:213
#define gPad
Definition: TVirtualPad.h:286
RooAbsArg * first() const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
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
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
void Print(Option_t *opts=0) const
Print contents of the workspace.
Bool_t defineSet(const char *name, const RooArgSet &aset, Bool_t importMissing=kFALSE)
Define a named RooArgSet with given constituents.
Bool_t saveSnapshot(const char *name, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of parameters 'params' If importValues ...
Bool_t loadSnapshot(const char *name)
Load the values and attributes of the parameters in the snapshot saved with the given name.
Bool_t import(const RooAbsArg &arg, 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())
Import a RooAbsArg object, e.g.
const RooArgSet * set(const char *name)
Return pointer to previously defined named set with given nmame If no such set is found a null pointe...
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
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 PrintLevel(Int_t code)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)