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