ROOT
Version v6.32
master
v6.34
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf503_wspaceread.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
##
5
## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #503
6
##
7
## Reading and using a workspace
8
##
9
## The input file for self macro is generated by rf502_wspaceread.py
10
##
11
## \macro_image
12
## \macro_code
13
## \macro_output
14
##
15
## \date February 2018
16
## \authors Clemens Lange, Wouter Verkerke (C version)
17
18
import
ROOT
19
20
21
# Read workspace from file
22
# -----------------------------------------------
23
24
# Open input file with workspace (generated by rf503_wspacewrite)
25
f =
ROOT.TFile
(
"rf502_workspace_py.root"
)
26
27
# Retrieve workspace from file
28
w =
f.Get
(
"w"
)
29
30
# Retrieve pdf, data from workspace
31
# -----------------------------------------------------------------
32
33
# Retrieve x, and data from workspace
34
x = w[
"x"
]
35
model = w[
"model"
]
36
data = w[
"modelData"
]
37
38
# Print structure of composite p.d.f.
39
model.Print
(
"t"
)
40
41
# Fit model to data, plot model
42
# ---------------------------------------------------------
43
44
# Fit model to data
45
model.fitTo
(data, PrintLevel=-1)
46
47
# Plot data and PDF overlaid
48
xframe =
x.frame
(Title=
"Model and data read from workspace"
)
49
data.plotOn
(xframe)
50
model.plotOn
(xframe)
51
52
# Overlay the background component of model with a dashed line
53
model.plotOn
(xframe, Components=
"bkg"
, LineStyle=
"--"
)
54
55
# Overlay the background+sig2 components of model with a dotted line
56
model.plotOn
(xframe, Components=
"bkg,sig2"
, LineStyle=
":"
)
57
58
# Draw the frame on the canvas
59
c =
ROOT.TCanvas
(
"rf503_wspaceread"
,
"rf503_wspaceread"
, 600, 600)
60
ROOT.gPad.SetLeftMargin
(0.15)
61
xframe.GetYaxis
().SetTitleOffset(1.4)
62
xframe.Draw
()
63
64
c.SaveAs
(
"rf503_wspaceread.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
rf503_wspaceread.py
ROOT v6-32 - Reference Guide Generated on Fri Mar 21 2025 05:17:58 (GVA Time) using Doxygen 1.10.0