Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf207_comptools.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #207
5## Tools and utilities for manipulation of composite objects
6##
7## \macro_code
8##
9## \date February 2018
10## \author Clemens Lange
11## \author Wouter Verkerke (C version)
12
13import ROOT
14
15# Set up composite pdf 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)
24sigma = ROOT.RooRealVar("sigma", "width of gaussians", 0.5)
25sig = ROOT.RooGaussian("sig", "Signal component 1", x, mean, sigma)
26
27# Build Chebychev polynomial p.d.f.
28a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
29a1 = ROOT.RooRealVar("a1", "a1", 0.2, 0.0, 1.0)
30bkg1 = ROOT.RooChebychev("bkg1", "Background 1", x, [a0, a1])
31
32# Build expontential pdf
33alpha = ROOT.RooRealVar("alpha", "alpha", -1)
34bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
35
36# Sum the background components into a composite background p.d.f.
37bkg1frac = ROOT.RooRealVar("bkg1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
38bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [bkg1frac])
39
40# Sum the composite signal and background
41bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
42model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
43
44# Create dummy dataset that has more observables than the above pdf
45y = ROOT.RooRealVar("y", "y", -10, 10)
46data = ROOT.RooDataSet("data", "data", {x, y})
47
48# Basic information requests
49# ---------------------------------------------
50
51
52# Get list of observables
53# ---------------------------------------------
54
55# Get list of observables of pdf in context of a dataset
56#
57# Observables are define each context as the variables
58# shared between a model and a dataset. In self case
59# that is the variable 'x'
60
61model_obs = model.getObservables(data)
62model_obs.Print("v")
63
64# Get list of parameters
65# -------------------------------------------
66
67# Get list of parameters, list of observables
68model_params = model.getParameters({x})
69model_params.Print("v")
70
71# Get list of parameters, a dataset
72# (Gives identical results to operation above)
73model_params2 = model.getParameters(data)
74model_params2.Print()
75
76# Get list of components
77# -------------------------------------------
78
79# Get list of component objects, top-level node
80model_comps = model.getComponents()
81model_comps.Print("v")
82
83# Modifications to structure of composites
84# -------------------------------------------
85
86# Create a second Gaussian
87sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
88sig2 = ROOT.RooGaussian("sig2", "Signal component 1", x, mean, sigma2)
89
90# Create a sum of the original Gaussian plus the second Gaussian
91sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
92sigsum = ROOT.RooAddPdf("sigsum", "sig+sig2", [sig, sig2], [sig1frac])
93
94# Construct a customizer utility to customize model
95cust = ROOT.RooCustomizer(model, "cust")
96
97# Instruct the customizer to replace node 'sig' with node 'sigsum'
98cust.replaceArg(sig, sigsum)
99
100# Build a clone of the input pdf according to the above customization
101# instructions. Each node that requires modified is clone so that the
102# original pdf remained untouched. The name of each cloned node is that
103# of the original node suffixed by the name of the customizer object
104#
105# The returned head node own all nodes that were cloned as part of
106# the build process so when cust_clone is deleted so will all other
107# nodes that were created in the process.
108cust_clone = cust.build(ROOT.kTRUE)
109
110# Print structure of clone of model with sig.sigsum replacement.
111cust_clone.Print("t")
112
113# The RooCustomizer has the be deleted first.
114# Otherwise, it might happen that `sig` or `sigsum` are deleted first, in which
115# case the internal TLists in the RooCustomizer will complain about deleted
116# objects.
117del cust