Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf312_multirangefit.py File Reference

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: performing fits in multiple (disjoint) ranges in one or more dimensions

import ROOT
# Create 2D pdf and data
# -------------------------------------------
# Define observables x,y
x = ROOT.RooRealVar("x", "x", -10, 10)
y = ROOT.RooRealVar("y", "y", -10, 10)
# Construct the signal pdf gauss(x)*gauss(y)
mx = ROOT.RooRealVar("mx", "mx", 1, -10, 10)
my = ROOT.RooRealVar("my", "my", 1, -10, 10)
gx = ROOT.RooGaussian("gx", "gx", x, mx, 1.0)
gy = ROOT.RooGaussian("gy", "gy", y, my, 1.0)
sig = ROOT.RooProdPdf("sig", "sig", gx, gy)
# Construct the background pdf (flat in x,y)
px = ROOT.RooPolynomial("px", "px", x)
py = ROOT.RooPolynomial("py", "py", y)
bkg = ROOT.RooProdPdf("bkg", "bkg", px, py)
# Construct the composite model sig+bkg
f = ROOT.RooRealVar("f", "f", 0.0, 1.0)
model = ROOT.RooAddPdf("model", "model", [sig, bkg], [f])
# Sample 10000 events in (x,y) from the model
modelData = model.generate({x, y}, 10000)
# Define signal and sideband regions
# -------------------------------------------------------------------
# Construct the SideBand1,SideBand2, regions
#
# |
# +-------------+-----------+
# | | |
# | Side | Sig |
# | Band1 | nal |
# | | |
# --+-------------+-----------+--
# | |
# | Side |
# | Band2 |
# | |
# +-------------+-----------+
# |
x.setRange("SB1", -10, +10)
y.setRange("SB1", -10, 0)
x.setRange("SB2", -10, 0)
y.setRange("SB2", 0, +10)
x.setRange("SIG", 0, +10)
y.setRange("SIG", 0, +10)
x.setRange("FULL", -10, +10)
y.setRange("FULL", -10, +10)
# Perform fits in individual sideband regions
# -------------------------------------------------------------------------------------
# Perform fit in SideBand1 region (ROOT.RooAddPdf coefficients will be
# interpreted in full range)
r_sb1 = model.fitTo(modelData, Range="SB1", Save=True, PrintLevel=-1)
# Perform fit in SideBand2 region (ROOT.RooAddPdf coefficients will be
# interpreted in full range)
r_sb2 = model.fitTo(modelData, Range="SB2", Save=True, PrintLevel=-1)
# Perform fits in joint sideband regions
# -----------------------------------------------------------------------------
# Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
# (ROOT.RooAddPdf coefficients will be interpreted in full range)
r_sb12 = model.fitTo(modelData, Range="SB1,SB2", Save=True, PrintLevel=-1)
# Print results for comparison
r_sb1.Print()
r_sb2.Print()
r_sb12.Print()
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'SB1' created with bounds [-10,10]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'SB1' created with bounds [-10,0]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'SB2' created with bounds [-10,0]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'SB2' created with bounds [0,10]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'SIG' created with bounds [0,10]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'SIG' created with bounds [0,10]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'FULL' created with bounds [-10,10]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'FULL' created with bounds [-10,10]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData' created with bounds [-10,10]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'fit_nll_model_modelData' created with bounds [-10,0]
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_SB1' created with bounds [-10,10]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'fit_nll_model_modelData_SB1' created with bounds [-10,0]
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'fit_nll_model_modelData_SB2' created with bounds [-10,0]
[#1] INFO:Eval -- RooRealVar::setRange(y) new range named 'fit_nll_model_modelData_SB2' created with bounds [0,10]
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_model_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
RooFitResult: minimized FCN value: 16261.4, estimated distance to minimum: 5.06415e-07
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
f 5.1135e-01 +/- 3.58e-02
mx 9.8740e-01 +/- 4.03e-02
my 9.9305e-01 +/- 9.39e-02
RooFitResult: minimized FCN value: 7578.28, estimated distance to minimum: 3.07865e-06
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
f 5.4586e-01 +/- 4.46e-02
mx 1.1276e+00 +/- 1.10e-01
my 9.6462e-01 +/- 5.60e-02
RooFitResult: minimized FCN value: 27252.6, estimated distance to minimum: 2.25202e-06
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
f 5.0082e-01 +/- 1.29e-02
mx 1.0100e+00 +/- 3.26e-02
my 9.6348e-01 +/- 3.31e-02
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf312_multirangefit.py.