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

Detailed Description

View in nbviewer Open in SWAN
'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #503

Reading and using a workspace

The input file for self macro is generated by rf502_wspaceread.py

import ROOT
# Read workspace from file
# -----------------------------------------------
# Open input file with workspace (generated by rf503_wspacewrite)
f = ROOT.TFile("rf502_workspace_py.root")
# Retrieve workspace from file
w = f.Get("w")
# Retrieve pdf, data from workspace
# -----------------------------------------------------------------
# Retrieve x, and data from workspace
x = w["x"]
model = w["model"]
data = w["modelData"]
# Print structure of composite p.d.f.
model.Print("t")
# Fit model to data, plot model
# ---------------------------------------------------------
# Fit model to data
model.fitTo(data, PrintLevel=-1)
# Plot data and PDF overlaid
xframe = x.frame(Title="Model and data read from workspace")
data.plotOn(xframe)
model.plotOn(xframe)
# Overlay the background component of model with a dashed line
model.plotOn(xframe, Components="bkg", LineStyle="--")
# Overlay the background+sig2 components of model with a dotted line
model.plotOn(xframe, Components="bkg,sig2", LineStyle=":")
# Draw the frame on the canvas
c = ROOT.TCanvas("rf503_wspaceread", "rf503_wspaceread", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.SaveAs("rf503_wspaceread.png")
0x8fdb600 RooAddPdf::model = 1/1 [Auto,Clean]
0x91994b0/V- RooChebychev::bkg = 1 [Auto,Dirty]
0x880af10/V- RooRealVar::x = 5
0x27be9c0/V- RooRealVar::a0 = 0.5
0x2765460/V- RooRealVar::a1 = 0
0x9279200/V- RooRealVar::bkgfrac = 0.5
0x926e430/V- RooAddPdf::sig = 1/1 [Auto,Clean]
0x91f8100/V- RooGaussian::sig1 = 1 [Auto,Dirty]
0x880af10/V- RooRealVar::x = 5
0x87e28d0/V- RooRealVar::mean = 5
0x4735f60/V- RooRealVar::sigma1 = 0.5
0x4735b70/V- RooRealVar::sig1frac = 0.8
0x9205290/V- RooGaussian::sig2 = 1 [Auto,Dirty]
0x880af10/V- RooRealVar::x = 5
0x87e28d0/V- RooRealVar::mean = 5
0x873e5a0/V- RooRealVar::sigma2 = 1
[#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:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg,sig2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: (sig)
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C version)

Definition in file rf503_wspaceread.py.