Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf509_wsinteractive.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Organization and simultaneous fits: easy interactive access to workspace contents - CINT to CLING code migration
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 04/2009 - Wouter Verkerke
10 
11 using namespace RooFit;
12 
13 void fillWorkspace(RooWorkspace &w);
14 
15 void rf509_wsinteractive()
16 {
17  // C r e a t e a n d f i l l w o r k s p a c e
18  // ------------------------------------------------
19 
20  // Create a workspace named 'w'
21  // With CINT w could exports its contents to
22  // a same-name C++ namespace in CINT 'namespace w'.
23  // but this does not work anymore in CLING.
24  // so this tutorial is an example on how to
25  // change the code
26  RooWorkspace *w1 = new RooWorkspace("w", kTRUE);
27 
28  // Fill workspace with p.d.f. and data in a separate function
29  fillWorkspace(*w1);
30 
31  // Print workspace contents
32  w1->Print();
33 
34  // this does not work anymore with CLING
35  // use normal workspace functionality
36 
37  // U s e w o r k s p a c e c o n t e n t s
38  // ----------------------------------------------
39 
40  // Old syntax to use the name space prefix operator to access the workspace contents
41  //
42  // RooDataSet* d = w::model.generate(w::x,1000) ;
43  // RooFitResult* r = w::model.fitTo(*d) ;
44 
45  // use normal workspace methods
46  RooAbsPdf *model = w1->pdf("model");
47  RooRealVar *x = w1->var("x");
48 
49  RooDataSet *d = model->generate(*x, 1000);
50  RooFitResult *r = model->fitTo(*d);
51 
52  // old syntax to access the variable x
53  // RooPlot* frame = w::x.frame() ;
54 
55  RooPlot *frame = x->frame();
56  d->plotOn(frame);
57 
58  // OLD syntax to omit x::
59  // NB: The 'w::' prefix can be omitted if namespace w is imported in local namespace
60  // in the usual C++ way
61  //
62  // using namespace w;
63  // model.plotOn(frame) ;
64  // model.plotOn(frame,Components(bkg),LineStyle(kDashed)) ;
65 
66  // new correct syntax
67  RooAbsPdf *bkg = w1->pdf("bkg");
68  model->plotOn(frame);
69  model->plotOn(frame, Components(*bkg), LineStyle(kDashed));
70 
71  // Draw the frame on the canvas
72  new TCanvas("rf509_wsinteractive", "rf509_wsinteractive", 600, 600);
73  gPad->SetLeftMargin(0.15);
74  frame->GetYaxis()->SetTitleOffset(1.4);
75  frame->Draw();
76 }
77 
78 void fillWorkspace(RooWorkspace &w)
79 {
80  // C r e a t e p d f a n d f i l l w o r k s p a c e
81  // --------------------------------------------------------
82 
83  // Declare observable x
84  RooRealVar x("x", "x", 0, 10);
85 
86  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
87  RooRealVar mean("mean", "mean of gaussians", 5, 0, 10);
88  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
89  RooRealVar sigma2("sigma2", "width of gaussians", 1);
90 
91  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
92  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
93 
94  // Build Chebychev polynomial p.d.f.
95  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
96  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
97  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
98 
99  // Sum the signal components into a composite signal p.d.f.
100  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
101  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
102 
103  // Sum the composite signal and background
104  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
105  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
106 
107  w.import(model);
108 }