Rf 3 0 7_Fullpereventerrors

Multidimensional models: usage of full pdf with per-event errors

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 Friday, May 13, 2022 at 09:26 AM.

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

B-physics pdf with per-event Gaussian resolution

Observables

In [2]:
dt = ROOT.RooRealVar("dt", "dt", -10, 10)
dterr = ROOT.RooRealVar("dterr", "per-event error on dt", 0.01, 10)

Build a gaussian resolution model scaled by the per-error = gauss(dt,bias,sigma*dterr)

In [3]:
bias = ROOT.RooRealVar("bias", "bias", 0, -10, 10)
sigma = ROOT.RooRealVar("sigma", "per-event error scale factor", 1, 0.1, 10)
gm = ROOT.RooGaussModel("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr)

Construct decay(dt) (x) gauss1(dt|dterr)

In [4]:
tau = ROOT.RooRealVar("tau", "tau", 1.548)
decay_gm = ROOT.RooDecay("decay_gm", "decay", dt, tau, gm, type="DoubleSided")

Construct empirical pdf for per-event error

Use landau pdf to get empirical distribution with long tail

In [5]:
pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, ROOT.RooFit.RooConst(1), ROOT.RooFit.RooConst(0.25))
expDataDterr = pdfDtErr.generate({dterr}, 10000)

Construct a histogram pdf to describe the shape of the dtErr distribution

In [6]:
expHistDterr = expDataDterr.binnedClone()
pdfErr = ROOT.RooHistPdf("pdfErr", "pdfErr", {dterr}, expHistDterr)

Construct conditional product decay_dm(dt|dterr)*pdf(dterr)

Construct production of conditional decay_dm(dt|dterr) with empirical pdfErr(dterr)

In [7]:
model = ROOT.RooProdPdf("model", "model", {pdfErr}, Conditional=({decay_gm}, {dt}))

(Alternatively you could also use the landau shape pdfDtErr) ROOT.RooProdPdf model("model", "model",pdfDtErr, ROOT.RooFit.Conditional(decay_gm,dt))

Sample, fit and plot product model

Specify external dataset with dterr values to use model_dm as conditional pdf

In [8]:
data = model.generate({dt, dterr}, 10000)

Fit conditional decay_dm(dt|dterr)

Specify dterr as conditional observable

In [9]:
model.fitTo(data)
Out[9]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --  The following expressions have been identified as constant and will be precalculated and cached: (pdfErr)
[#1] INFO:Minimization --  The following expressions will be evaluated in cache-and-track mode: (decay_gm)
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 bias         0.00000e+00  2.00000e+00   -1.00000e+01  1.00000e+01
     2 sigma        1.00000e+00  4.50000e-01    1.00000e-01  1.00000e+01
 **********
 **    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        1000           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
 START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
 FCN=36045.2 FROM MIGRAD    STATUS=INITIATE        8 CALLS           9 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  bias         0.00000e+00   2.00000e+00   2.01358e-01   3.32978e+02
   2  sigma        1.00000e+00   4.50000e-01   1.63378e-01   6.12093e+01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=36045 FROM MIGRAD    STATUS=CONVERGED      23 CALLS          24 TOTAL
                     EDM=2.04993e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  bias        -1.00111e-02   1.72603e-02   2.27324e-04  -2.30921e+00
   2  sigma        9.91174e-01   2.01704e-02   9.34197e-04  -3.06107e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  2.979e-04 -2.715e-06 
 -2.715e-06  4.069e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.00780   1.000 -0.008
        2  0.00780  -0.008  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=36045 FROM HESSE     STATUS=OK             10 CALLS          34 TOTAL
                     EDM=2.04714e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  bias        -1.00111e-02   1.72602e-02   4.54647e-05  -1.00111e-03
   2  sigma        9.91174e-01   2.01706e-02   1.86839e-04  -9.61350e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  2.979e-04 -3.117e-06 
 -3.117e-06  4.069e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.00895   1.000 -0.009
        2  0.00895  -0.009  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot conditional decay_dm(dt|dterr)

Make two-dimensional plot of conditional pdf in (dt,dterr)

In [10]:
hh_model = model.createHistogram("hh_model", dt, Binning=50, YVar=dict(var=dterr, Binning=50))
hh_model.SetLineColor(ROOT.kBlue)

Make projection of data an dt

In [11]:
frame = dt.frame(Title="Projection of model(dt|dterr) on dt")
data.plotOn(frame)
model.plotOn(frame)
Out[11]:
<cppyy.gbl.RooPlot object at 0x88663a0>
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on dt integrates over variables (dterr)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([pdfErr_NORM[dterr]_X_decay_gm_NORM[dt]]_Int[dterr]) using numeric integrator RooIntegrator1D to calculate Int(dterr)

Draw all frames on canvas

In [12]:
c = ROOT.TCanvas("rf307_fullpereventerrors", "rf307_fullpereventerrors", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.20)
hh_model.GetZaxis().SetTitleOffset(2.5)
hh_model.Draw("surf")
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.6)
frame.Draw()

c.SaveAs("rf307_fullpereventerrors.png")
Info in <TCanvas::Print>: png file rf307_fullpereventerrors.png has been created

Draw all canvases

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