rf103_interprfuncs

Basic functionality: interpreted functions and pdfs

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:06 AM.

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

Generic interpreted pdf

Declare observable x

In [2]:
x = ROOT.RooRealVar("x", "x", -20, 20)

Construct generic pdf from interpreted expression

ROOT.To construct a proper pdf, the formula expression is explicitly normalized internally by dividing it by a numeric integral of the expresssion over x in the range [-20,20]

In [3]:
alpha = ROOT.RooRealVar("alpha", "alpha", 5, 0.1, 10)
genpdf = ROOT.RooGenericPdf("genpdf", "genpdf", "(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))", [x, alpha])

Sample, fit and plot generic pdf

Generate a toy dataset from the interpreted pdf

In [4]:
data = genpdf.generate({x}, 10000)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)

Fit the interpreted pdf to the generated data

In [5]:
genpdf.fitTo(data)
Out[5]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 alpha        5.00000e+00  9.90000e-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         500           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=35708.9 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  alpha        5.00000e+00   9.90000e-01   2.01369e-01   1.02940e+02
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=35708.6 FROM MIGRAD    STATUS=CONVERGED      12 CALLS          13 TOTAL
                     EDM=1.24247e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  alpha        4.96355e+00   4.16961e-02   1.10376e-03  -4.18389e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  1.739e-03 
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE         500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=35708.6 FROM HESSE     STATUS=OK              5 CALLS          18 TOTAL
                     EDM=1.24183e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  alpha        4.96355e+00   4.16961e-02   2.20752e-04  -1.74664e-02
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  1.739e-03 
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Make a plot of the data and the pdf overlaid

In [6]:
xframe = x.frame(Title="Interpreted expression pdf")
data.plotOn(xframe)
genpdf.plotOn(xframe)
Out[6]:
<cppyy.gbl.RooPlot object at 0xa352800>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(genpdf_Int[x]) using numeric integrator RooIntegrator1D to calculate Int(x)

Standard pdf adjust with interpreted helper function

Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian #

Construct standard pdf with formula replacing parameter

Construct parameter mean2 and sigma

In [7]:
mean2 = ROOT.RooRealVar("mean2", "mean^2", 10, 0, 200)
sigma = ROOT.RooRealVar("sigma", "sigma", 3, 0.1, 10)

Construct interpreted function mean = sqrt(mean^2)

In [8]:
mean = ROOT.RooFormulaVar("mean", "mean", "sqrt(mean2)", [mean2])

Construct a gaussian g2(x,sqrt(mean2),sigma)

In [9]:
g2 = ROOT.RooGaussian("g2", "h2", x, mean, sigma)

Generate toy data

Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian dataset with mean 10 and width 3

In [10]:
g1 = ROOT.RooGaussian("g1", "g1", x, ROOT.RooFit.RooConst(10), ROOT.RooFit.RooConst(3))
data2 = g1.generate({x}, 1000)

Fit and plot tailored standard pdf

Fit g2 to data from g1

In [11]:
r = g2.fitTo(data2, Save=True)  # ROOT.RooFitResult
r.Print()
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 mean2        1.00000e+01  5.00000e+00    0.00000e+00  2.00000e+02
     2 sigma        3.00000e+00  9.90000e-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=5148.93 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  mean2        1.00000e+01   5.00000e+00   1.18625e-01  -5.23438e+03
   2  sigma        3.00000e+00   9.90000e-01   2.22742e-01  -7.90389e+03
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=2551.39 FROM MIGRAD    STATUS=CONVERGED      59 CALLS          60 TOTAL
                     EDM=8.7852e-06    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  mean2        1.00100e+02   1.98019e+00   6.89576e-04   4.58015e-02
   2  sigma        3.11719e+00   7.12427e-02   5.29831e-04   1.79331e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  3.922e+00  2.826e-03 
  2.826e-03  5.076e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.02003   1.000  0.020
        2  0.02003   0.020  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=2551.39 FROM HESSE     STATUS=OK             10 CALLS          70 TOTAL
                     EDM=8.78617e-06    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  mean2        1.00100e+02   1.98016e+00   1.37915e-04   1.00004e-03
   2  sigma        3.11719e+00   7.12418e-02   1.05966e-04  -4.01138e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  3.922e+00  2.730e-03 
  2.730e-03  5.076e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.01935   1.000  0.019
        2  0.01935   0.019  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

  RooFitResult: minimized FCN value: 2551.39, estimated distance to minimum: 8.78617e-06
                covariance matrix quality: Full, accurate covariance matrix
                Status : MINIMIZE=0 HESSE=0 

    Floating Parameter    FinalValue +/-  Error   
  --------------------  --------------------------
                 mean2    1.0010e+02 +/-  1.98e+00
                 sigma    3.1172e+00 +/-  7.12e-02

Plot data on frame and overlay projection of g2

In [12]:
xframe2 = x.frame(Title="Tailored Gaussian pdf")
data2.plotOn(xframe2)
g2.plotOn(xframe2)
Out[12]:
<cppyy.gbl.RooPlot object at 0xa104a20>

Draw all frames on a canvas

In [13]:
c = ROOT.TCanvas("rf103_interprfuncs", "rf103_interprfuncs", 800, 400)
c.Divide(2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.4)
xframe.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
xframe2.GetYaxis().SetTitleOffset(1.4)
xframe2.Draw()

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

Draw all canvases

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