Multidimensional models: multi-dimensional pdfs through composition, e.g.
substituting a pdf parameter with a function that depends on other observables
pdf = gauss(x,f(y),s)
with f(y) = a0 + a1*y
import ROOT
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -5, 5)
a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
fy = ROOT.RooPolyVar("fy", "fy", y, [a0, a1])
sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
model = ROOT.RooGaussian("model", "Gaussian with shifting mean", x, fy, sigma)
data = model.generate({x, y}, 10000)
xframe = x.frame()
data.plotOn(xframe)
model.plotOn(xframe)
yframe = y.frame()
data.plotOn(yframe)
model.plotOn(yframe)
hh_model = model.createHistogram("hh_model", x, Binning=50, YVar=dict(var=y, Binning=50))
hh_model.SetLineColor(ROOT.kBlue)
c = ROOT.TCanvas("rf301_composition", "rf301_composition", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
yframe.GetYaxis().SetTitleOffset(1.4)
yframe.Draw()
c.cd(3)
ROOT.gPad.SetLeftMargin(0.20)
hh_model.GetZaxis().SetTitleOffset(2.5)
hh_model.Draw("surf")
c.SaveAs("rf301_composition.png")
[#0] WARNING:InputArguments -- The parameter 'sigma' with range [-inf, inf] of the RooGaussian 'model' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on x integrates over variables (y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[y]_Norm[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on y integrates over variables (x)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_Int[x,y]) using numeric integrator RooRombergIntegrator to calculate Int(y)
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf301_composition.py.