Logo ROOT  
Reference Guide
rf407_latextables.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4##
5## Data and categories: latex printing of lists and sets of RooArgSets
6##
7## \macro_code
8##
9## \date February 2018
10## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14
15# Setup composite pdf
16# --------------------------------------
17
18# Declare observable x
19x = ROOT.RooRealVar("x", "x", 0, 10)
20
21# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
22# their parameters
23mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
24sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
25sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
26sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
27sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
28
29# Sum the signal components into a composite signal p.d.f.
30sig1frac = ROOT.RooRealVar(
31 "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
32sig = ROOT.RooAddPdf(
33 "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
34
35# Build Chebychev polynomial p.d.f.
36a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
37a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
38bkg1 = ROOT.RooChebychev("bkg1", "Background 1",
39 x, ROOT.RooArgList(a0, a1))
40
41# Build expontential pdf
42alpha = ROOT.RooRealVar("alpha", "alpha", -1)
43bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
44
45# Sum the background components into a composite background p.d.f.
46bkg1frac = ROOT.RooRealVar(
47 "sig1frac", "fraction of component 1 in background", 0.2, 0., 1.)
48bkg = ROOT.RooAddPdf(
49 "bkg", "Signal", ROOT.RooArgList(bkg1, bkg2), ROOT.RooArgList(sig1frac))
50
51# Sum the composite signal and background
52bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
53model = ROOT.RooAddPdf(
54 "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
55
56# Make list of parameters before and after fit
57# ----------------------------------------------------------------------------------------
58
59# Make list of model parameters
60params = model.getParameters(ROOT.RooArgSet(x))
61
62# Save snapshot of prefit parameters
63initParams = params.snapshot()
64
65# Do fit to data, obtain error estimates on parameters
66data = model.generate(ROOT.RooArgSet(x), 1000)
67model.fitTo(data)
68
69# Print LateX table of parameters of pdf
70# --------------------------------------------------------------------------
71
72# Print parameter list in LaTeX for (one column with names, column with
73# values)
74params.printLatex()
75
76# Print parameter list in LaTeX for (names values|names values)
77params.printLatex(ROOT.RooFit.Columns(2))
78
79# Print two parameter lists side by side (name values initvalues)
80params.printLatex(ROOT.RooFit.Sibling(initParams))
81
82# Print two parameter lists side by side (name values initvalues|name
83# values initvalues)
84params.printLatex(ROOT.RooFit.Sibling(initParams), ROOT.RooFit.Columns(2))
85
86# Write LaTex table to file
87params.printLatex(ROOT.RooFit.Sibling(initParams),
88 ROOT.RooFit.OutputFile("rf407_latextables.tex"))