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