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
rf515_hfJSON.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Code HistFactory Models in JSON.
5
##
6
## With the HS3 standard, it is possible to code RooFit-Models of any kind as JSON files.
7
## In this tutorial, you can see how to code up a (simple) HistFactory-based model in JSON and import it into a RooWorkspace.
8
##
9
## \macro_code
10
##
11
## \date November 2021
12
## \author Carsten Burgard
13
14
import
ROOT
15
16
# start by creating an empty workspace
17
ws =
ROOT.RooWorkspace
(
"workspace"
)
18
19
# the RooJSONFactoryWSTool is responsible for importing and exporting things to and from your workspace
20
tool =
ROOT.RooJSONFactoryWSTool
(ws)
21
22
# use it to import the information from your JSON file
23
tool.importJSON
(
ROOT.gROOT.GetTutorialDir
().Data() +
"/roofit/rf515_hfJSON.json"
)
24
ws.Print
()
25
26
# now, you can easily use your workspace to run your fit (as you usually would)
27
# the model config is named after your pdf, i.e. <the pdf name>_modelConfig
28
model = ws[
"ModelConfig"
]
29
30
# for resetting the parameters after the fit
31
params =
model.GetPdf
().getParameters(ws[
"observed"
])
32
params_initial =
params.snapshot
()
33
34
# we are fitting a clone of the model now,
35
result =
model.fitTo
(ws[
"observed"
],
ROOT.RooFit.Save
(),
ROOT.RooFit.PrintLevel
(-1))
36
result.Print
()
37
# reset parameters, such that we are not double-fitting the model in the
38
# closure check.
39
params.assign
(params_initial)
40
41
# in the end, you can again write to json
42
# the result will be not completely identical to the JSON file you used as an input, but it will work just the same
43
tool.exportJSON
(
"myWorkspace.json"
)
44
45
# You can again import it if you want and check for closure
46
ws_2 =
ROOT.RooWorkspace
(
"workspace"
)
47
tool_2 =
ROOT.RooJSONFactoryWSTool
(ws_2)
48
tool_2.importJSON
(
"myWorkspace.json"
)
49
ws_2.Print
()
50
model_2 = ws_2[
"ModelConfig"
]
51
result =
model_2.fitTo
(ws_2[
"observed"
],
ROOT.RooFit.Save
(),
ROOT.RooFit.PrintLevel
(-1))
52
result.Print
()
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
rf515_hfJSON.py
ROOT v6-32 - Reference Guide Generated on Mon Apr 7 2025 05:38:39 (GVA Time) using Doxygen 1.10.0