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