Logo ROOT   master
Reference Guide
rf202_extendedmlfit.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ##
5 ## Addition and convolution: setting up an extended maximum likelihood fit
6 ##
7 ## \macro_code
8 ##
9 ## \date February 2018
10 ## \authors Clemens Lange, Wouter Verkerke (C++ version)
11 
12 import ROOT
13 
14 # Set up component pdfs
15 # ---------------------------------------
16 
17 # Declare observable x
18 x = 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
22 mean = ROOT.RooRealVar("mean", "mean of gaussians", 5)
23 sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
24 sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
25 
26 sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
27 sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
28 
29 # Build Chebychev polynomial p.d.f.
30 a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0., 1.)
31 a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0., 1.)
32 bkg = ROOT.RooChebychev("bkg", "Background", x, ROOT.RooArgList(a0, a1))
33 
34 # Sum the signal components into a composite signal p.d.f.
35 sig1frac = ROOT.RooRealVar(
36  "sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.)
37 sig = ROOT.RooAddPdf(
38  "sig", "Signal", ROOT.RooArgList(sig1, sig2), ROOT.RooArgList(sig1frac))
39 
40 # Method 1 - Construct extended composite model
41 # -------------------------------------------------------------------
42 
43 # Sum the composite signal and background into an extended pdf
44 # nsig*sig+nbkg*bkg
45 nsig = ROOT.RooRealVar("nsig", "number of signal events", 500, 0., 10000)
46 nbkg = ROOT.RooRealVar(
47  "nbkg", "number of background events", 500, 0, 10000)
48 model = ROOT.RooAddPdf(
49  "model",
50  "(g1+g2)+a",
51  ROOT.RooArgList(
52  bkg,
53  sig),
54  ROOT.RooArgList(
55  nbkg,
56  nsig))
57 
58 # Sample, fit and plot extended model
59 # ---------------------------------------------------------------------
60 
61 # Generate a data sample of expected number events in x from model
62 # = model.expectedEvents() = nsig+nbkg
63 data = model.generate(ROOT.RooArgSet(x))
64 
65 # Fit model to data, ML term automatically included
66 model.fitTo(data)
67 
68 # Plot data and PDF overlaid, expected number of events for p.d.f projection normalization
69 # rather than observed number of events (==data.numEntries())
70 xframe = x.frame(ROOT.RooFit.Title("extended ML fit example"))
71 data.plotOn(xframe)
72 model.plotOn(xframe, ROOT.RooFit.Normalization(
73  1.0, ROOT.RooAbsReal.RelativeExpected))
74 
75 # Overlay the background component of model with a dashed line
76 ras_bkg = ROOT.RooArgSet(bkg)
77 model.plotOn(
78  xframe, ROOT.RooFit.Components(ras_bkg), ROOT.RooFit.LineStyle(
79  ROOT.kDashed), ROOT.RooFit.Normalization(
80  1.0, ROOT.RooAbsReal.RelativeExpected))
81 
82 # Overlay the background+sig2 components of model with a dotted line
83 ras_bkg_sig2 = ROOT.RooArgSet(bkg, sig2)
84 model.plotOn(
85  xframe, ROOT.RooFit.Components(ras_bkg_sig2), ROOT.RooFit.LineStyle(
86  ROOT.kDotted), ROOT.RooFit.Normalization(
87  1.0, ROOT.RooAbsReal.RelativeExpected))
88 
89 # Print structure of composite p.d.f.
90 model.Print("t")
91 
92 
93 # Method 2 - Construct extended components first
94 # ---------------------------------------------------------------------
95 
96 # Associated nsig/nbkg as expected number of events with sig/bkg
97 esig = ROOT.RooExtendPdf("esig", "extended signal p.d.f", sig, nsig)
98 ebkg = ROOT.RooExtendPdf("ebkg", "extended background p.d.f", bkg, nbkg)
99 
100 # Sum extended components without coefs
101 # -------------------------------------------------------------------------
102 
103 # Construct sum of two extended p.d.f. (no coefficients required)
104 model2 = ROOT.RooAddPdf("model2", "(g1+g2)+a", ROOT.RooArgList(ebkg, esig))
105 
106 # Draw the frame on the canvas
107 c = ROOT.TCanvas("rf202_extendedmlfit", "rf202_extendedmlfit", 600, 600)
108 ROOT.gPad.SetLeftMargin(0.15)
109 xframe.GetYaxis().SetTitleOffset(1.4)
110 xframe.Draw()
111 
112 c.SaveAs("rf202_extendedmlfit.png")