ROOT
Version v6.34
master
v6.32
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
## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #503
5
##
6
## Reading and using a workspace
7
##
8
## The input file for self macro is generated by rf502_wspaceread.py
9
##
10
## \macro_image
11
## \macro_code
12
## \macro_output
13
##
14
## \date February 2018
15
## \authors Clemens Lange, Wouter Verkerke (C version)
16
17
import
ROOT
18
19
20
# Read workspace from file
21
# -----------------------------------------------
22
23
# Open input file with workspace (generated by rf503_wspacewrite)
24
f =
ROOT.TFile
(
"rf502_workspace_py.root"
)
25
26
# Retrieve workspace from file
27
w =
f.Get
(
"w"
)
28
29
# Retrieve pdf, data from workspace
30
# -----------------------------------------------------------------
31
32
# Retrieve x, and data from workspace
33
x = w[
"x"
]
34
model = w[
"model"
]
35
data = w[
"modelData"
]
36
37
# Print structure of composite p.d.f.
38
model.Print
(
"t"
)
39
40
# Fit model to data, plot model
41
# ---------------------------------------------------------
42
43
# Fit model to data
44
model.fitTo
(data, PrintLevel=-1)
45
46
# Plot data and PDF overlaid
47
xframe =
x.frame
(Title=
"Model and data read from workspace"
)
48
data.plotOn
(xframe)
49
model.plotOn
(xframe)
50
51
# Overlay the background component of model with a dashed line
52
model.plotOn
(xframe, Components=
"bkg"
, LineStyle=
"--"
)
53
54
# Overlay the background+sig2 components of model with a dotted line
55
model.plotOn
(xframe, Components=
"bkg,sig2"
, LineStyle=
":"
)
56
57
# Draw the frame on the canvas
58
c =
ROOT.TCanvas
(
"rf503_wspaceread"
,
"rf503_wspaceread"
, 600, 600)
59
ROOT.gPad.SetLeftMargin
(0.15)
60
xframe.GetYaxis
().SetTitleOffset(1.4)
61
xframe.Draw
()
62
63
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 tags/6-34-04 - Reference Guide Generated on Wed Mar 19 2025 04:35:07 (GVA Time) using Doxygen 1.10.0