Logo ROOT   6.16/01
Reference Guide
rf502_wspacewrite.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4##
5## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #502
6##
7## Creating and writing a workspace
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15import ROOT
16
17
18# Create model and dataset
19# -----------------------------------------------
20
21# Declare observable x
22x = ROOT.RooRealVar("x", "x", 0, 10)
23
24# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
25# their parameters
26mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
27sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
28sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
29
30sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
31sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
32
33# Build Chebychev polynomial p.d.f.
34a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
35a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
36bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
37
38# Sum the signal components into a composite signal p.d.f.
39sig1frac = ROOT.RooRealVar(
40 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
41sig = ROOT.RooAddPdf(
42 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
43
44# Sum the composite signal and background
45bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
46model = ROOT.RooAddPdf(
47 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
48
49# Generate a data sample of 1000 events in x from model
50data = model.generate(ROOT.RooArgSet(x), 1000)
51
52# Create workspace, import data and model
53# -----------------------------------------------------------------------------
54
55# Create a empty workspace
56w = ROOT.RooWorkspace("w", "workspace")
57
58# Import model and all its components into the workspace
59getattr(w, 'import')(model)
60
61# Import data into the workspace
62getattr(w, 'import')(data)
63
64# Print workspace contents
65w.Print()
66
67# Save workspace in file
68# -------------------------------------------
69
70# Save the workspace into a ROOT file
71w.writeToFile("rf502_workspace.root")
72
73# Workspace will remain in memory after macro finishes
74ROOT.gDirectory.Add(w)