Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf510_wsnamedsets.py File Reference

Detailed Description

View in nbviewer Open in SWAN

'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #510

Working with named parameter sets and parameter snapshots in workspaces

import ROOT
def fillWorkspace(w):
# Create model
# -----------------------
# Declare observable x
x = ROOT.RooRealVar("x", "x", 0, 10)
# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
# their parameters
mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
# Build Chebychev polynomial p.d.f.
a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
# Sum the signal components into a composite signal p.d.f.
sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
# Sum the composite signal and background
bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
# Import model into p.d.f.
w.Import(model)
# Encode definition of parameters in workspace
# ---------------------------------------------------------------------------------------
# Define named sets "parameters" and "observables", list which variables should be considered
# parameters and observables by the users convention
#
# Variables appearing in sets _must_ live in the workspace already, the autoImport flag
# of defineSet must be set to import them on the fly. Named sets contain only references
# to the original variables, the value of observables in named sets already
# reflect their 'current' value
params = model.getParameters({x})
w.defineSet("parameters", params)
w.defineSet("observables", {x})
# Encode reference value for parameters in workspace
# ---------------------------------------------------------------------------------------------------
# Define a parameter 'snapshot' in the p.d.f.
# Unlike a named set, parameter snapshot stores an independent set of values for
# a given set of variables in the workspace. The values can be stored and reloaded
# into the workspace variable objects using the loadSnapshot() and saveSnapshot()
# methods. A snapshot saves the value of each variable, errors that are stored
# with it as well as the 'Constant' flag that is used in fits to determine if a
# parameter is kept fixed or not.
# Do a dummy fit to a (supposedly) reference dataset here and store the results
# of that fit into a snapshot
refData = model.generate({x}, 10000)
model.fitTo(refData, PrintLevel=-1)
# The kTRUE flag imports the values of the objects in (*params) into the workspace
# If not set, present values of the workspace parameters objects are stored
w.saveSnapshot("reference_fit", params, True)
# Make another fit with the signal componentforced to zero
# and save those parameters too
bkgfrac.setVal(1)
bkgfrac.setConstant(True)
bkgfrac.removeError()
model.fitTo(refData, PrintLevel=-1)
w.saveSnapshot("reference_fit_bkgonly", params, True)
# Create model and dataset
# -----------------------------------------------
w = ROOT.RooWorkspace("w")
fillWorkspace(w)
# Exploit convention encoded in named set "parameters" and "observables"
# to use workspace contents w/o need for introspected
model = w["model"]
# Generate data from p.d.f. in given observables
data = model.generate(w.set("observables"), 1000)
# Fit model to data
model.fitTo(data, PrintLevel=-1)
# Plot fitted model and data on frame of first (only) observable
frame = (w.set("observables").first()).frame()
data.plotOn(frame)
model.plotOn(frame)
# Overlay plot with model with reference parameters as stored in snapshots
w.loadSnapshot("reference_fit")
model.plotOn(frame, LineColor="r")
w.loadSnapshot("reference_fit_bkgonly")
model.plotOn(frame, LineColor="r", LineStyle="--")
# Draw the frame on the canvas
c = ROOT.TCanvas("rf510_wsnamedsets", "rf503_wsnamedsets", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()
c.SaveAs("rf510_wsnamedsets.png")
# Print workspace contents
w.Print()
# Workspace will remain in memory after macro finishes
ROOT.gDirectory.Add(w)
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-inf, inf] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-inf, inf] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooChebychev::bkg
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::x
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a0
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::bkgfrac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::sig
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::mean
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sig1frac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig2
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma2
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooWorkspace(w) w contents
variables
---------
(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2,x)
p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 1
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 1/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 0.999807/1
RooGaussian::sig1[ x=x mean=mean sigma=sigma1 ] = 0.999765
RooGaussian::sig2[ x=x mean=mean sigma=sigma2 ] = 0.999941
parameter snapshots
-------------------
reference_fit = (a0=0.488363 +/- 0.0241765,a1=0.0249649 +/- 0.0384439,bkgfrac=0.484424 +/- 0.0115194,mean=5.01085 +/- 0.0103629,sigma1=0.5[C],sig1frac=0.757184 +/- 0.0357496,sigma2=1[C])
reference_fit_bkgonly = (a0=0.469143 +/- 0,a1=7.71653e-08 +/- 0,bkgfrac=1[C],mean=5.01085 +/- 0,sigma1=0.5[C],sig1frac=0.757184 +/- 0,sigma2=1[C])
named sets
----------
observables:(x)
parameters:(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2)
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C version)

Definition in file rf510_wsnamedsets.py.