Logo ROOT  
Reference Guide
 
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
14import ROOT
15
16# start by creating an empty workspace
17ws = ROOT.RooWorkspace("workspace")
18
19# the RooJSONFactoryWSTool is responsible for importing and exporting things to and from your workspace
20tool = ROOT.RooJSONFactoryWSTool(ws)
21
22# use it to import the information from your JSON file
23tool.importJSON(ROOT.gROOT.GetTutorialDir().Data() + "/roofit/rf515_hfJSON.json")
24ws.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
28model = ws["ModelConfig"]
29
30# for resetting the parameters after the fit
31params = model.GetPdf().getParameters(ws["observed"])
32ROOT.SetOwnership(params, True)
33params_initial = params.snapshot()
34ROOT.SetOwnership(params_initial, True)
35
36# we are fitting a clone of the model now,
37result = model.fitTo(ws["observed"], ROOT.RooFit.Save(), ROOT.RooFit.PrintLevel(-1))
38ROOT.SetOwnership(result, True)
39result.Print()
40# reset parameters, such that we are not double-fitting the model in the
41# closure check.
42params.assign(params_initial)
43
44# in the end, you can again write to json
45# the result will be not completely identical to the JSON file you used as an input, but it will work just the same
46tool.exportJSON("myWorkspace.json")
47
48# You can again import it if you want and check for closure
49ws_2 = ROOT.RooWorkspace("workspace")
50tool_2 = ROOT.RooJSONFactoryWSTool(ws_2)
51tool_2.importJSON("myWorkspace.json")
52ws_2.Print()
53model_2 = ws_2["ModelConfig"]
54result = model_2.fitTo(ws_2["observed"], ROOT.RooFit.Save(), ROOT.RooFit.PrintLevel(-1))
55ROOT.SetOwnership(result, True)
56result.Print()