Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf504_simwstool.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4## Organization and simultaneous fits: using RooSimWSTool to construct a simultaneous pdf
5## that is built of variations of an input pdf
6##
7## \macro_code
8## \macro_output
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Create master pdf
17# ---------------------------------
18
19# Construct gauss(x,m,s)
20x = ROOT.RooRealVar("x", "x", -10, 10)
21m = ROOT.RooRealVar("m", "m", 0, -10, 10)
22s = ROOT.RooRealVar("s", "s", 1, -10, 10)
23gauss = ROOT.RooGaussian("g", "g", x, m, s)
24
25# Construct poly(x,p0)
26p0 = ROOT.RooRealVar("p0", "p0", 0.01, 0.0, 1.0)
27poly = ROOT.RooPolynomial("p", "p", x, [p0])
28
29# model = f*gauss(x) + (1-f)*poly(x)
30f = ROOT.RooRealVar("f", "f", 0.5, 0.0, 1.0)
31model = ROOT.RooAddPdf("model", "model", [gauss, poly], [f])
32
33# Create category observables for splitting
34# ----------------------------------------------------------------------------------
35
36# Define two categories that can be used for splitting
37c = ROOT.RooCategory("c", "c")
38c.defineType("run1")
39c.defineType("run2")
40
41d = ROOT.RooCategory("d", "d")
42d.defineType("foo")
43d.defineType("bar")
44
45# Set up SimWSTool
46# -----------------------------
47
48# Import ingredients in a workspace
49w = ROOT.RooWorkspace("w", "w")
50w.Import({model, c, d})
51
52# Make Sim builder tool
53sct = ROOT.RooSimWSTool(w)
54
55# Build a simultaneous model with one split
56# ---------------------------------------------------------------------------------
57
58# Construct a simultaneous pdf with the following form
59#
60# model_run1(x) = f*gauss_run1(x,m_run1,s) + (1-f)*poly
61# model_run2(x) = f*gauss_run2(x,m_run2,s) + (1-f)*poly
62# simpdf(x,c) = model_run1(x) if c=="run1"
63# = model_run2(x) if c=="run2"
64#
65# Returned pdf is owned by the workspace
66model_sim = sct.build("model_sim", "model", SplitParam=("m", "c"))
67
68# Print tree structure of model
69model_sim.Print("t")
70
71# Adjust model_sim parameters in workspace
72w.var("m_run1").setVal(-3)
73w.var("m_run2").setVal(+3)
74
75# Print contents of workspace
76w.Print("v")
77
78# Build a simultaneous model with product split
79# -----------------------------------------------------------------------------------------
80
81# Build another simultaneous pdf using a composite split in states c X d
82model_sim2 = sct.build("model_sim2", "model", SplitParam=("p0", "c,d"))
83
84# Print tree structure of self model
85model_sim2.Print("t")