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