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

Detailed Description

View in nbviewer Open in SWAN
Organization and simultaneous fits: reading and writing ASCII configuration files

from __future__ import print_function
import ROOT
# Create pdf
# ------------------
# Construct gauss(x,m,s)
x = ROOT.RooRealVar("x", "x", -10, 10)
m = ROOT.RooRealVar("m", "m", 0, -10, 10)
s = ROOT.RooRealVar("s", "s", 1, -10, 10)
gauss = ROOT.RooGaussian("g", "g", x, m, s)
# Construct poly(x,p0)
p0 = ROOT.RooRealVar("p0", "p0", 0.01, 0.0, 1.0)
poly = ROOT.RooPolynomial("p", "p", x, [p0])
# model = f*gauss(x) + (1-f)*poly(x)
f = ROOT.RooRealVar("f", "f", 0.5, 0.0, 1.0)
model = ROOT.RooAddPdf("model", "model", [gauss, poly], [f])
# Fit model to toy data
# -----------------------------------------
d = model.generate({x}, 1000)
model.fitTo(d, PrintLevel=-1)
# Write parameters to ASCII file
# -----------------------------------------------------------
# Obtain set of parameters
params = model.getParameters({x})
# Write parameters to file
params.writeToFile("rf505_asciicfg_example.txt")
# Read parameters from ASCII file
# ----------------------------------------------------------------
# Read parameters from file
params.readFromFile("rf505_asciicfg_example.txt")
params.Print("v")
configFile = ROOT.gROOT.GetTutorialDir().Data() + "/roofit/rf505_asciicfg.txt"
# Read parameters from section 'Section2' of file
params.readFromFile(configFile, "", "Section2")
params.Print("v")
# Read parameters from section 'Section3' of file. Mark all
# variables that were processed with the "READ" attribute
params.readFromFile(configFile, "READ", "Section3")
# Print the list of parameters that were not read from Section3
print("The following parameters of the were _not_ read from Section3: ", params.selectByAttrib("READ", False))
# Read parameters from section 'Section4' of file, contains
# 'include file' statement of rf505_asciicfg_example.txt
# so that we effective read the same
params.readFromFile(configFile, "", "Section4")
params.Print("v")
[#0] WARNING:InputArguments -- The parameter 's' with range [-10, 10] of the RooGaussian 'g' exceeds the safe range of (0, inf). Advise to limit its range.
[#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) 0x576f0e0 RooRealVar:: f = 0.50733 +/- 0.020971 L(0 - 1) "f"
2) 0x7842470 RooRealVar:: m = 0.0064018 +/- 0.053686 L(-10 - 10) "m"
3) 0x5778900 RooRealVar:: p0 = 0.0073509 +/- 0.0078312 L(0 - 1) "p0"
4) 0x7779890 RooRealVar:: s = 0.96516 +/- 0.047052 L(-10 - 10) "s"
[#1] INFO:InputArguments -- RooArgSet::readFromStream(parameters): processing include file rf505_asciicfg_example.txt
1) 0x576f0e0 RooRealVar:: f = 0.45 +/- 0.03 L(0 - 1) "f"
2) 0x7842470 RooRealVar:: m = 0.025 +/- 0.02 L(-10 - 10) "m"
3) 0x5778900 RooRealVar:: p0 = 0.0022 +/- 0.0001 L(0 - 1) "p0"
4) 0x7779890 RooRealVar:: s = 0.98 +/- 0.03 L(-10 - 10) "s"
[#1] INFO:InputArguments -- RooArgSet::readFromStream(parameters): processing include file rf505_asciicfg_example.txt
[#1] INFO:InputArguments -- RooArgSet::readFromStream(parameters): processing include file rf505_asciicfg_example.txt
1) 0x576f0e0 RooRealVar:: f = 0.372 C L(0 - 1) "f"
2) 0x7842470 RooRealVar:: m = 0.195 C L(-10 - 10) "m"
3) 0x5778900 RooRealVar:: p0 = 0.0022 +/- 0.0001 L(0 - 1) "p0"
4) 0x7779890 RooRealVar:: s = 0.98 +/- 0.03 L(-10 - 10) "s"
The following parameters of the were _not_ read from Section3: { @0x8ca0090, @0x8ca0098 }
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf505_asciicfg.py.