Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf511_wsfactory_basic.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// Organization and simultaneous fits: basic use of the 'object factory' associated with a workspace to rapidly build
5 /// p.d.f.s functions and their parameter components
6 ///
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 using namespace RooFit;
22 
23 void rf511_wsfactory_basic(Bool_t compact = kFALSE)
24 {
25  RooWorkspace *w = new RooWorkspace("w");
26 
27  // C r e a t i n g a n d a d d i n g b a s i c p . d . f . s
28  // ----------------------------------------------------------------
29 
30  // Remake example p.d.f. of tutorial rs502_wspacewrite.C:
31  //
32  // Basic p.d.f. construction: ClassName::ObjectName(constructor arguments)
33  // Variable construction : VarName[x,xlo,xhi], VarName[xlo,xhi], VarName[x]
34  // P.d.f. addition : SUM::ObjectName(coef1*pdf1,...coefM*pdfM,pdfN)
35  //
36 
37  if (!compact) {
38 
39  // Use object factory to build p.d.f. of tutorial rs502_wspacewrite
40  w->factory("Gaussian::sig1(x[-10,10],mean[5,0,10],0.5)");
41  w->factory("Gaussian::sig2(x,mean,1)");
42  w->factory("Chebychev::bkg(x,{a0[0.5,0.,1],a1[0.2,0.,1.]})");
43  w->factory("SUM::sig(sig1frac[0.8,0.,1.]*sig1,sig2)");
44  w->factory("SUM::model(bkgfrac[0.5,0.,1.]*bkg,sig)");
45 
46  } else {
47 
48  // Use object factory to build p.d.f. of tutorial rs502_wspacewrite but
49  // - Contracted to a single line recursive expression,
50  // - Omitting explicit names for components that are not referred to explicitly later
51 
52  w->factory("SUM::model(bkgfrac[0.5,0.,1.]*Chebychev::bkg(x[-10,10],{a0[0.5,0.,1],a1[0.2,0.,1.]}),"
53  "SUM(sig1frac[0.8,0.,1.]*Gaussian(x,mean[5,0,10],0.5), Gaussian(x,mean,1)))");
54  }
55 
56  // A d v a n c e d p . d . f . c o n s t r u c t o r a r g u m e n t s
57  // ----------------------------------------------------------------
58  //
59  // P.d.f. constructor arguments may by any type of RooAbsArg, but also
60  //
61  // Double_t --> converted to RooConst(...)
62  // {a,b,c} --> converted to RooArgSet() or RooArgList() depending on required ctor arg
63  // dataset name --> converted to RooAbsData reference for any dataset residing in the workspace
64  // enum --> Any enum label that belongs to an enum defined in the (base) class
65 
66  // Make a dummy dataset p.d.f. 'model' and import it in the workspace
67  RooDataSet *data = w->pdf("model")->generate(*w->var("x"), 1000);
68  w->import(*data, Rename("data"));
69 
70  // Construct a KEYS p.d.f. passing a dataset name and an enum type defining the
71  // mirroring strategy
72  w->factory("KeysPdf::k(x,data,NoMirror,0.2)");
73 
74  // Print workspace contents
75  w->Print();
76 }