Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf801_mcstudy.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Validation and MC studies: toy Monte Carlo study that perform cycles of event generation and fitting
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13
14# Create model
15# -----------------------
16
17# Declare observable x
18x = ROOT.RooRealVar("x", "x", 0, 10)
19x.setBins(40)
20
21# Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
22# their parameters
23mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
24sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
25sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
26
27sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
28sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
29
30# Build Chebychev polynomial pdf
31a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
32a1 = ROOT.RooRealVar("a1", "a1", -0.2, -1, 1.0)
33bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
34
35# Sum the signal components into a composite signal pdf
36sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
37sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
38
39# Sum the composite signal and background
40nbkg = ROOT.RooRealVar("nbkg", "number of background events, ", 150, 0, 1000)
41nsig = ROOT.RooRealVar("nsig", "number of signal events", 150, 0, 1000)
42model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [nbkg, nsig])
43
44# Create manager
45# ---------------------------
46
47# Instantiate ROOT.RooMCStudy manager on model with x as observable and given choice of fit options
48#
49# The Silence() option kills all messages below the PROGRESS level, only a single message
50# per sample executed, any error message that occur during fitting
51#
52# The Extended() option has two effects:
53# 1) The extended ML term is included in the likelihood and
54# 2) A poisson fluctuation is introduced on the number of generated events
55#
56# The FitOptions() given here are passed to the fitting stage of each toy experiment.
57# If Save() is specified, fit result of each experiment is saved by the manager
58#
59# A Binned() option is added in self example to bin the data between generation and fitting
60# to speed up the study at the expemse of some precision
61
62mcstudy = ROOT.RooMCStudy(
63 model,
64 {x},
65 Binned=True,
66 Silence=True,
67 Extended=True,
68 FitOptions=dict(Save=True, PrintEvalErrors=0),
69)
70
71# Generate and fit events
72# ---------------------------------------------
73
74# Generate and fit 1000 samples of Poisson(nExpected) events
75mcstudy.generateAndFit(1000)
76
77# Explore results of study
78# ------------------------------------------------
79
80# Make plots of the distributions of mean, error on mean and the pull of
81# mean
82frame1 = mcstudy.plotParam(mean, Bins=40)
83frame2 = mcstudy.plotError(mean, Bins=40)
84frame3 = mcstudy.plotPull(mean, Bins=40, FitGauss=True)
85
86# Plot distribution of minimized likelihood
87frame4 = mcstudy.plotNLL(Bins=40)
88
89# Make some histograms from the parameter dataset
90hh_cor_a0_s1f = ROOT.RooAbsData.createHistogram(mcstudy.fitParDataSet(), "hh", a1, YVar=sig1frac)
91hh_cor_a0_a1 = ROOT.RooAbsData.createHistogram(mcstudy.fitParDataSet(), "hh", a0, YVar=a1)
92
93# Access some of the saved fit results from individual toys
94corrHist000 = mcstudy.fitResult(0).correlationHist("c000")
95corrHist127 = mcstudy.fitResult(127).correlationHist("c127")
96corrHist953 = mcstudy.fitResult(953).correlationHist("c953")
97
98# Draw all plots on a canvas
99ROOT.gStyle.SetPalette(1)
100ROOT.gStyle.SetOptStat(0)
101c = ROOT.TCanvas("rf801_mcstudy", "rf801_mcstudy", 900, 900)
102c.Divide(3, 3)
103c.cd(1)
104ROOT.gPad.SetLeftMargin(0.15)
105frame1.GetYaxis().SetTitleOffset(1.4)
106frame1.Draw()
107c.cd(2)
108ROOT.gPad.SetLeftMargin(0.15)
109frame2.GetYaxis().SetTitleOffset(1.4)
110frame2.Draw()
111c.cd(3)
112ROOT.gPad.SetLeftMargin(0.15)
113frame3.GetYaxis().SetTitleOffset(1.4)
114frame3.Draw()
115c.cd(4)
116ROOT.gPad.SetLeftMargin(0.15)
117frame4.GetYaxis().SetTitleOffset(1.4)
118frame4.Draw()
119c.cd(5)
120ROOT.gPad.SetLeftMargin(0.15)
121hh_cor_a0_s1f.GetYaxis().SetTitleOffset(1.4)
122hh_cor_a0_s1f.Draw("box")
123c.cd(6)
124ROOT.gPad.SetLeftMargin(0.15)
125hh_cor_a0_a1.GetYaxis().SetTitleOffset(1.4)
126hh_cor_a0_a1.Draw("box")
127c.cd(7)
128ROOT.gPad.SetLeftMargin(0.15)
129corrHist000.GetYaxis().SetTitleOffset(1.4)
130corrHist000.Draw("colz")
131c.cd(8)
132ROOT.gPad.SetLeftMargin(0.15)
133corrHist127.GetYaxis().SetTitleOffset(1.4)
134corrHist127.Draw("colz")
135c.cd(9)
136ROOT.gPad.SetLeftMargin(0.15)
137corrHist953.GetYaxis().SetTitleOffset(1.4)
138corrHist953.Draw("colz")
139
140c.SaveAs("rf801_mcstudy.png")
141
142# Make ROOT.RooMCStudy object available on command line after
143# macro finishes
144ROOT.gDirectory.Add(mcstudy)