Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf603_multicpu.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Likelihood and minimization: setting up a multi-core parallelized unbinned maximum likelihood fit
5##
6## \macro_image
7## \macro_code
8## \macro_output
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Create 3D pdf and data
17# -------------------------------------------
18
19# Create observables
20x = ROOT.RooRealVar("x", "x", -5, 5)
21y = ROOT.RooRealVar("y", "y", -5, 5)
22z = ROOT.RooRealVar("z", "z", -5, 5)
23
24# Create signal pdf gauss(x)*gauss(y)*gauss(z)
25gx = ROOT.RooGaussian("gx", "gx", x, 0.0, 1.0)
26gy = ROOT.RooGaussian("gy", "gy", y, 0.0, 1.0)
27gz = ROOT.RooGaussian("gz", "gz", z, 0.0, 1.0)
28sig = ROOT.RooProdPdf("sig", "sig", [gx, gy, gz])
29
30# Create background pdf poly(x)*poly(y)*poly(z)
31px = ROOT.RooPolynomial("px", "px", x, [-0.1, 0.004])
32py = ROOT.RooPolynomial("py", "py", y, [0.1, -0.004])
33pz = ROOT.RooPolynomial("pz", "pz", z)
34bkg = ROOT.RooProdPdf("bkg", "bkg", [px, py, pz])
35
36# Create composite pdf sig+bkg
37fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0.0, 1.0)
38model = ROOT.RooAddPdf("model", "model", [sig, bkg], [fsig])
39
40# Generate large dataset
41data = model.generate({x, y, z}, 200000)
42
43# Parallel fitting
44# -------------------------------
45
46# In parallel mode the likelihood calculation is split in N pieces,
47# that are calculated in parallel and added a posteriori before passing
48# it back to MINUIT.
49
50# Use four processes and time results both in wall time and CPU time
51model.fitTo(data, NumCPU=4, Timer=True, PrintLevel=-1)
52
53# Parallel MC projections
54# ----------------------------------------------
55
56# Construct signal, likelihood projection on (y,z) observables and
57# likelihood ratio
58sigyz = sig.createProjection({x})
59totyz = model.createProjection({x})
60llratio_func = ROOT.RooFormulaVar("llratio", "log10(@0)-log10(@1)", [sigyz, totyz])
61
62# Calculate likelihood ratio for each event, subset of events with high
63# signal likelihood
64data.addColumn(llratio_func)
65dataSel = data.reduce(Cut="llratio>0.7")
66
67# Make plot frame and plot data
68frame = x.frame(Title="Projection on X with LLratio(y,z)>0.7", Bins=40)
69dataSel.plotOn(frame)
70
71# Perform parallel projection using MC integration of pdf using given input dataSet.
72# In self mode the data-weighted average of the pdf is calculated by splitting the
73# input dataset in N equal pieces and calculating in parallel the weighted average
74# one each subset. The N results of those calculations are then weighted into the
75# final result
76
77# Use four processes
78model.plotOn(frame, ProjWData=dataSel, NumCPU=4)
79
80c = ROOT.TCanvas("rf603_multicpu", "rf603_multicpu", 600, 600)
81ROOT.gPad.SetLeftMargin(0.15)
82frame.GetYaxis().SetTitleOffset(1.6)
83frame.Draw()
84
85c.SaveAs("rf603_multicpu.png")