Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf201_composite.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ## Addition and convolution: composite p.d.f with signal and background component
5 ##
6 ## pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
7 ##
8 ## \macro_code
9 ##
10 ## \date February 2018
11 ## \author Clemens Lange, Wouter Verkerke (C++ version)
12 
13 import ROOT
14 
15 # Setup component pdfs
16 # ---------------------------------------
17 
18 # Declare observable x
19 x = 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
23 mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
24 sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
25 sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
26 
27 sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
28 sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
29 
30 # Build Chebychev polynomial p.d.f.
31 a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
32 a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
33 bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
34 
35 
36 # Method 1 - Two RooAddPdfs
37 # ------------------------------------------
38 # Add signal components
39 
40 # Sum the signal components into a composite signal p.d.f.
41 sig1frac = ROOT.RooRealVar(
42  "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
43 sig = ROOT.RooAddPdf("sig", "Signal", ROOT.RooArgList(
44  sig1, sig2), ROOT.RooArgList(sig1frac))
45 
46 # Add signal and background
47 # ------------------------------------------------
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 # Sample, fit and plot model
55 # ---------------------------------------------------
56 
57 # Generate a data sample of 1000 events in x from model
58 data = model.generate(ROOT.RooArgSet(x), 1000)
59 
60 # Fit model to data
61 model.fitTo(data)
62 
63 # Plot data and PDF overlaid
64 xframe = x.frame(ROOT.RooFit.Title(
65  "Example of composite pdf=(sig1+sig2)+bkg"))
66 data.plotOn(xframe)
67 model.plotOn(xframe)
68 
69 # Overlay the background component of model with a dashed line
70 ras_bkg = ROOT.RooArgSet(bkg)
71 model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg),
72  ROOT.RooFit.LineStyle(ROOT.kDashed))
73 
74 # Overlay the background+sig2 components of model with a dotted line
75 ras_bkg_sig2 = ROOT.RooArgSet(bkg, sig2)
76 model.plotOn(xframe, ROOT.RooFit.Components(ras_bkg_sig2),
77  ROOT.RooFit.LineStyle(ROOT.kDotted))
78 
79 # Print structure of composite p.d.f.
80 model.Print("t")
81 
82 # Method 2 - One RooAddPdf with recursive fractions
83 # ---------------------------------------------------
84 
85 # Construct sum of models on one go using recursive fraction interpretations
86 #
87 # model2 = bkg + (sig1 + sig2)
88 #
89 model2 = ROOT.RooAddPdf(
90  "model",
91  "g1+g2+a",
92  ROOT.RooArgList(
93  bkg,
94  sig1,
95  sig2),
96  ROOT.RooArgList(
97  bkgfrac,
98  sig1frac),
99  ROOT.kTRUE)
100 
101 # NB: Each coefficient is interpreted as the fraction of the
102 # left-hand component of the i-th recursive sum, i.e.
103 #
104 # sum4 = A + ( B + ( C + D) with fraction fA, and fC expands to
105 #
106 # sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
107 
108 # Plot recursive addition model
109 # ---------------------------------------------------------
110 model2.plotOn(xframe, ROOT.RooFit.LineColor(ROOT.kRed),
111  ROOT.RooFit.LineStyle(ROOT.kDashed))
112 model2.plotOn(
113  xframe,
114  ROOT.RooFit.Components(ras_bkg_sig2),
115  ROOT.RooFit.LineColor(
116  ROOT.kRed),
117  ROOT.RooFit.LineStyle(
118  ROOT.kDashed))
119 model2.Print("t")
120 
121 # Draw the frame on the canvas
122 c = ROOT.TCanvas("rf201_composite", "rf201_composite", 600, 600)
123 ROOT.gPad.SetLeftMargin(0.15)
124 xframe.GetYaxis().SetTitleOffset(1.4)
125 xframe.Draw()
126 
127 c.SaveAs("rf201_composite.png")