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