Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf510_wsnamedsets.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Organization and simultaneous fits: working with named parameter sets and parameter snapshots in workspaces
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 04/2009 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooConstVar.h"
15 #include "RooChebychev.h"
16 #include "RooAddPdf.h"
17 #include "RooWorkspace.h"
18 #include "RooPlot.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "TFile.h"
22 #include "TH1.h"
23 
24 using namespace RooFit;
25 
26 void fillWorkspace(RooWorkspace &w);
27 
28 void rf510_wsnamedsets()
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  RooWorkspace *w = new RooWorkspace("w");
34  fillWorkspace(*w);
35 
36  // Exploit convention encoded in named set "parameters" and "observables"
37  // to use workspace contents w/o need for introspected
38  RooAbsPdf *model = w->pdf("model");
39 
40  // Generate data from p.d.f. in given observables
41  RooDataSet *data = model->generate(*w->set("observables"), 1000);
42 
43  // Fit model to data
44  model->fitTo(*data);
45 
46  // Plot fitted model and data on frame of first (only) observable
47  RooPlot *frame = ((RooRealVar *)w->set("observables")->first())->frame();
48  data->plotOn(frame);
49  model->plotOn(frame);
50 
51  // Overlay plot with model with reference parameters as stored in snapshots
52  w->loadSnapshot("reference_fit");
53  model->plotOn(frame, LineColor(kRed));
54  w->loadSnapshot("reference_fit_bkgonly");
55  model->plotOn(frame, LineColor(kRed), LineStyle(kDashed));
56 
57  // Draw the frame on the canvas
58  new TCanvas("rf510_wsnamedsets", "rf503_wsnamedsets", 600, 600);
59  gPad->SetLeftMargin(0.15);
60  frame->GetYaxis()->SetTitleOffset(1.4);
61  frame->Draw();
62 
63  // Print workspace contents
64  w->Print();
65 
66  // Workspace will remain in memory after macro finishes
67  gDirectory->Add(w);
68 }
69 
70 void fillWorkspace(RooWorkspace &w)
71 {
72  // C r e a t e m o d e l
73  // -----------------------
74 
75  // Declare observable x
76  RooRealVar x("x", "x", 0, 10);
77 
78  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
79  RooRealVar mean("mean", "mean of gaussians", 5, 0, 10);
80  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
81  RooRealVar sigma2("sigma2", "width of gaussians", 1);
82 
83  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
84  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
85 
86  // Build Chebychev polynomial p.d.f.
87  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
88  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
89  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
90 
91  // Sum the signal components into a composite signal p.d.f.
92  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
93  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
94 
95  // Sum the composite signal and background
96  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
97  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
98 
99  // Import model into p.d.f.
100  w.import(model);
101 
102  // 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
103  // ---------------------------------------------------------------------------------------
104 
105  // Define named sets "parameters" and "observables", which list which variables should be considered
106  // parameters and observables by the users convention
107  //
108  // Variables appearing in sets _must_ live in the workspace already, or the autoImport flag
109  // of defineSet must be set to import them on the fly. Named sets contain only references
110  // to the original variables, therefore the value of observables in named sets already
111  // reflect their 'current' value
112  RooArgSet *params = (RooArgSet *)model.getParameters(x);
113  w.defineSet("parameters", *params);
114  w.defineSet("observables", x);
115 
116  // 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
117  // ---------------------------------------------------------------------------------------------------
118 
119  // Define a parameter 'snapshot' in the p.d.f.
120  // Unlike a named set, a parameter snapshot stores an independent set of values for
121  // a given set of variables in the workspace. The values can be stored and reloaded
122  // into the workspace variable objects using the loadSnapshot() and saveSnapshot()
123  // methods. A snapshot saves the value of each variable, any errors that are stored
124  // with it as well as the 'Constant' flag that is used in fits to determine if a
125  // parameter is kept fixed or not.
126 
127  // Do a dummy fit to a (supposedly) reference dataset here and store the results
128  // of that fit into a snapshot
129  RooDataSet *refData = model.generate(x, 10000);
130  model.fitTo(*refData, PrintLevel(-1));
131 
132  // The kTRUE flag imports the values of the objects in (*params) into the workspace
133  // If not set, the present values of the workspace parameters objects are stored
134  w.saveSnapshot("reference_fit", *params, kTRUE);
135 
136  // Make another fit with the signal component forced to zero
137  // and save those parameters too
138 
139  bkgfrac.setVal(1);
140  bkgfrac.setConstant(kTRUE);
141  bkgfrac.removeError();
142  model.fitTo(*refData, PrintLevel(-1));
143 
144  w.saveSnapshot("reference_fit_bkgonly", *params, kTRUE);
145 }