rf315_projectpdf

Multidimensional models: marginizalization of multi-dimensional pdfs through integration

Author: Clemens Lange, Wouter Verkerke (C++ version)
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Sunday, November 27, 2022 at 11:07 AM.

In [1]:
import ROOT
Welcome to JupyROOT 6.27/01

Create pdf m(x,y) = gx(x|y) * g(y)

Increase default precision of numeric integration as self exercise has high sensitivity to numeric integration precision

In [2]:
ROOT.RooAbsPdf.defaultIntegratorConfig().setEpsRel(1e-8)
ROOT.RooAbsPdf.defaultIntegratorConfig().setEpsAbs(1e-8)

Create observables

In [3]:
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -2, 2)

Create function f(y) = a0 + a1*y

In [4]:
a0 = ROOT.RooRealVar("a0", "a0", 0)
a1 = ROOT.RooRealVar("a1", "a1", -1.5, -3, 1)
fy = ROOT.RooPolyVar("fy", "fy", y, [a0, a1])

Create gaussx(x,f(y),sx)

In [5]:
sigmax = ROOT.RooRealVar("sigmax", "width of gaussian", 0.5)
gaussx = ROOT.RooGaussian("gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax)
[#0] WARNING:InputArguments -- The parameter 'sigmax' with range [-1e+30, 1e+30] of the RooGaussian 'gaussx' exceeds the safe range of (0, inf). Advise to limit its range.

Create gaussy(y,0,2)

In [6]:
gaussy = ROOT.RooGaussian("gaussy", "Gaussian in y", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))

Create gaussx(x,sx|y) * gaussy(y)

In [7]:
model = ROOT.RooProdPdf(
    "model",
    "gaussx(x|y)*gaussy(y)",
    {gaussy},
    Conditional=({gaussx}, {x}),
)

Marginalize m(x,y) to m(x)

modelx(x) = Int model(x,y) dy

In [8]:
modelx = model.createProjection({y})

Use marginalized pdf as regular 1D pdf

Sample 1000 events from modelx

In [9]:
data = modelx.generateBinned({x}, 1000)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooIntegrator1D to calculate Int(y)

Fit modelx to toy data

In [10]:
modelx.fitTo(data, Verbose=True)
Out[10]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooIntegrator1D to calculate Int(y)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --  The following expressions will be evaluated in cache-and-track mode: (gaussx)
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for a1: using 0.4
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 a1          -1.50000e+00  4.00000e-01   -3.00000e+00  1.00000e+00
 **********
 **    3 **SET ERR         0.5
 **********
 **********
 **    4 **SET PRINT           1
 **********
 **********
 **    5 **SET STR           1
 **********
 NOW USING STRATEGY  1: TRY TO BALANCE SPEED AGAINST RELIABILITY
 **********
 **    6 **MIGRAD         500           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.

prevFCN = 1900.156536   START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
a1=-1.491, 
prevFCN = 1900.000822  a1=-1.509, 
prevFCN = 1900.423606  a1=-1.499, 
prevFCN = 1900.135899  a1=-1.501, 
prevFCN = 1900.178287   FCN=1900.16 FROM MIGRAD    STATUS=INITIATE        4 CALLS           5 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a1          -1.50000e+00   4.00000e-01   2.08372e-01  -4.77802e+01
                               ERR DEF= 0.5
a1=-1.484, 
prevFCN = 1899.958651  a1=-1.483, 
prevFCN = 1899.959675  a1=-1.484, 
prevFCN = 1899.958586  a1=-1.484, 
prevFCN = 1899.958497  a1=-1.483, 
prevFCN = 1899.958935  a1=-1.485, 
prevFCN = 1899.958964   MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
a1=-1.484, 
prevFCN = 1899.958497  a1=-1.483, 
prevFCN = 1899.958935  a1=-1.485, 
prevFCN = 1899.958964  a1=-1.484, 
prevFCN = 1899.958512  a1=-1.484, 
prevFCN = 1899.958518   COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=1899.96 FROM MIGRAD    STATUS=CONVERGED      15 CALLS          16 TOTAL
                     EDM=2.41873e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a1          -1.48409e+00   2.51054e-02   3.89173e-04  -3.80132e-02
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  6.303e-04 
a1=-1.484,  **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE         500
 **********

prevFCN = 1899.958497  a1=-1.484, 
prevFCN = 1899.958512  a1=-1.484, 
prevFCN = 1899.958518  a1=-1.484, 
prevFCN = 1899.958497  a1=-1.484, 
prevFCN = 1899.958498   COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=1899.96 FROM HESSE     STATUS=OK              5 CALLS          21 TOTAL
                     EDM=2.42353e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  a1          -1.48409e+00   2.51054e-02   7.78346e-05  -2.44471e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  6.303e-04 
a1=-1.484, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot modelx over data

In [11]:
frame = x.frame(40)
data.plotOn(frame)
modelx.plotOn(frame)
Out[11]:
<cppyy.gbl.RooPlot object at 0x9fc16c0>
[#1] INFO:NumericIntegration -- RooRealIntegral::init([gaussy_NORM[y]_X_gaussx_NORM[x]]_Int[y]) using numeric integrator RooIntegrator1D to calculate Int(y)

Make 2D histogram of model(x,y)

In [12]:
hh = model.createHistogram("x,y")
hh.SetLineColor(ROOT.kBlue)

c = ROOT.TCanvas("rf315_projectpdf", "rf315_projectpdf", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.20)
hh.GetZaxis().SetTitleOffset(2.5)
hh.Draw("surf")
c.SaveAs("rf315_projectpdf.png")
Info in <TCanvas::Print>: png file rf315_projectpdf.png has been created

Draw all canvases

In [13]:
from ROOT import gROOT 
gROOT.GetListOfCanvases().Draw()