Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf206_treevistools.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook -nodraw
4 ## Addition and convolution: tools for visualization of ROOT.RooAbsArg expression trees
5 ##
6 ## \macro_code
7 ##
8 ## \date February 2018
9 ## \author Clemens Lange, Wouter Verkerke (C++ version)
10 
11 import ROOT
12 
13 # Set up composite pdf
14 # --------------------------------------
15 
16 # Declare observable x
17 x = ROOT.RooRealVar("x", "x", 0, 10)
18 
19 # Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
20 # their parameters
21 mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
22 sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
23 sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
24 sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
25 sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
26 
27 # Sum the signal components into a composite signal p.d.f.
28 sig1frac = ROOT.RooRealVar(
29  "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
30 sig = ROOT.RooAddPdf(
31  "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
32 
33 # Build Chebychev polynomial p.d.f.
34 a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
35 a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
36 bkg1 = ROOT.RooChebychev("bkg1", "Background 1",
37  x, ROOT.RooArgList(a0, a1))
38 
39 # Build expontential pdf
40 alpha = ROOT.RooRealVar("alpha", "alpha", -1)
41 bkg2 = ROOT.RooExponential("bkg2", "Background 2", x, alpha)
42 
43 # Sum the background components into a composite background p.d.f.
44 bkg1frac = ROOT.RooRealVar(
45  "bkg1frac", "fraction of component 1 in background", 0.2, 0., 1.)
46 bkg = ROOT.RooAddPdf(
47  "bkg", "Signal", ROOT.RooArgList(bkg1, bkg2), ROOT.RooArgList(bkg1frac))
48 
49 # Sum the composite signal and background
50 bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0., 1.)
51 model = ROOT.RooAddPdf(
52  "model", "g1+g2+a", ROOT.RooArgList(bkg, sig), ROOT.RooArgList(bkgfrac))
53 
54 # Print composite tree in ASCII
55 # -----------------------------------------------------------
56 
57 # Print tree to stdout
58 model.Print("t")
59 
60 # Print tree to file
61 model.printCompactTree("", "rf206_asciitree.txt")
62 
63 # Draw composite tree graphically
64 # -------------------------------------------------------------
65 
66 # Print GraphViz DOT file with representation of tree
67 model.graphVizTree("rf206_model.dot")
68 
69 # Make graphic output file with one of the GraphViz tools
70 # (freely available from www.graphviz.org)
71 #
72 # 'Top-to-bottom graph'
73 # unix> dot -Tgif -o rf207_model_dot.gif rf207_model.dot
74 #
75 # 'Spring-model graph'
76 # unix> fdp -Tgif -o rf207_model_fdp.gif rf207_model.dot