Logo ROOT  
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: creating and writing a workspace
6##
7## \macro_code
8##
9## \date February 2018
10## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14
15# Create model and dataset
16# -----------------------------------------------
17
18# Declare observable x
19x = ROOT.RooRealVar("x", "x", 0, 10)
20
21# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
22# their parameters
23mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
24sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
25sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
26
27sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
28sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
29
30# Build Chebychev polynomial p.d.f.
31a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
32a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
33bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
34
35# Sum the signal components into a composite signal p.d.f.
36sig1frac = ROOT.RooRealVar(
37 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
38sig = ROOT.RooAddPdf(
39 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
40
41# Sum the composite signal and background
42bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
43model = ROOT.RooAddPdf(
44 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
45
46# Generate a data sample of 1000 events in x from model
47data = model.generate(ROOT.RooArgSet(x), 1000)
48
49# Create workspace, import data and model
50# -----------------------------------------------------------------------------
51
52# Create a empty workspace
53w = ROOT.RooWorkspace("w", "workspace")
54
55# Import model and all its components into the workspace
56w.Import(model)
57
58# Import data into the workspace
59w.Import(data)
60
61# Print workspace contents
62w.Print()
63
64# Save workspace in file
65# -------------------------------------------
66
67# Save the workspace into a ROOT file
68w.writeToFile("rf502_workspace.root")
69