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("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
30sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
31
32# Build Chebychev polynomial pdf
33a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
34a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
35bkg1 = ROOT.RooChebychev("bkg1", "Background 1", x, [a0, a1])
36
37# Build expontential pdf
38alpha = ROOT.RooRealVar("alpha", "alpha", -1)
39bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
40
41# Sum the background components into a composite background pdf
42bkg1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in background", 0.2, 0.0, 1.0)
43bkg = ROOT.RooAddPdf("bkg", "Signal", [bkg1, bkg2], [sig1frac])
44
45# Sum the composite signal and background
46bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
47model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
48
49# Make list of parameters before and after fit
50# ----------------------------------------------------------------------------------------
51
52# Make list of model parameters
53params = model.getParameters({x})
54
55# Save snapshot of prefit parameters
56initParams = params.snapshot()
57
58# Do fit to data, obtain error estimates on parameters
59data = model.generate({x}, 1000)
60model.fitTo(data)
61
62# Print LateX table of parameters of pdf
63# --------------------------------------------------------------------------
64
65# Print parameter list in LaTeX for (one column with names, column with
66# values)
67params.printLatex()
68
69# Print parameter list in LaTeX for (names values|names values)
70params.printLatex(Columns=2)
71
72# Print two parameter lists side by side (name values initvalues)
73params.printLatex(Sibling=initParams)
74
75# Print two parameter lists side by side (name values initvalues|name
76# values initvalues)
77params.printLatex(Sibling=initParams, Columns=2)
78
79# Write LaTex table to file
80params.printLatex(Sibling=initParams, OutputFile="rf407_latextables.tex")