rf105_funcbinding

'BASIC FUNCTIONALITY' RooFit tutorial macro #105 Demonstration of binding ROOT Math functions as RooFit 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

Bind ROOT TMath::Erf C function

Bind one-dimensional ROOT.TMath.Erf function as ROOT.RooAbsReal function

In [2]:
x = ROOT.RooRealVar("x", "x", -3, 3)
erf = ROOT.RooFit.bindFunction("erf", ROOT.TMath.Erf, x)

Print erf definition

In [3]:
erf.Print()
RooCFunction1Binding<double,double>::erf[ function=TMath::Erf x=x ] = 0

Plot erf on frame

In [4]:
frame1 = x.frame(Title="TMath.Erf bound as ROOT.RooFit function")
erf.plotOn(frame1)
Out[4]:
<cppyy.gbl.RooPlot object at 0x8eaf770>

Bind ROOT::Math::beta_pdf C function

Bind pdf ROOT.Math.Beta with three variables as ROOT.RooAbsPdf function

In [5]:
x2 = ROOT.RooRealVar("x2", "x2", 0, 0.999)
a = ROOT.RooRealVar("a", "a", 5, 0, 10)
b = ROOT.RooRealVar("b", "b", 2, 0, 10)
beta = ROOT.RooFit.bindPdf("beta", ROOT.Math.beta_pdf, x2, a, b)

Perf beta definition

In [6]:
beta.Print()
RooCFunction3PdfBinding<double,double,double,double>::beta[ function=ROOT::Math::beta_pdf x=x2 a=a b=b ] = 0.934689

Generate some events and fit

In [7]:
data = beta.generate({x2}, 10000)
beta.fitTo(data)
Out[7]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(beta_Int[x2]) using numeric integrator RooIntegrator1D to calculate Int(x2)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(beta_Int[x2]) using numeric integrator RooIntegrator1D to calculate Int(x2)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(beta_Int[x2]) using numeric integrator RooIntegrator1D to calculate Int(x2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 a            5.00000e+00  1.00000e+00    0.00000e+00  1.00000e+01
     2 b            2.00000e+00  1.00000e+00    0.00000e+00  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=-4851.7 FROM MIGRAD    STATUS=INITIATE       10 CALLS          11 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a            5.00000e+00   1.00000e+00   2.01358e-01  -5.86406e+01
   2  b            2.00000e+00   1.00000e+00   2.57889e-01   3.98974e+02
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-4851.9 FROM MIGRAD    STATUS=CONVERGED      39 CALLS          40 TOTAL
                     EDM=4.45661e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a            4.99036e+00   7.12139e-02   3.69873e-04   8.65173e-02
   2  b            1.98812e+00   2.63908e-02   1.71677e-04  -1.47761e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  5.072e-03  1.582e-03 
  1.582e-03  6.965e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.84198   1.000  0.842
        2  0.84198   0.842  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=-4851.9 FROM HESSE     STATUS=OK             10 CALLS          50 TOTAL
                     EDM=4.46009e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  a            4.99036e+00   7.12317e-02   7.39745e-05  -1.92762e-03
   2  b            1.98812e+00   2.63974e-02   3.43355e-05  -6.46475e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  5.074e-03  1.583e-03 
  1.583e-03  6.968e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.84207   1.000  0.842
        2  0.84207   0.842  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot data and pdf on frame

In [8]:
frame2 = x2.frame(Title="ROOT.Math.Beta bound as ROOT.RooFit pdf")
data.plotOn(frame2)
beta.plotOn(frame2)
Out[8]:
<cppyy.gbl.RooPlot object at 0x96c6d10>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(beta_Int[x2]) using numeric integrator RooIntegrator1D to calculate Int(x2)

Bind ROOT TF1 as RooFit function

Create a ROOT TF1 function

In [9]:
fa1 = ROOT.TF1("fa1", "sin(x)/x", 0, 10)

Create an observable

In [10]:
x3 = ROOT.RooRealVar("x3", "x3", 0.01, 20)

Create binding of TF1 object to above observable

In [11]:
rfa1 = ROOT.RooFit.bindFunction(fa1, x3)

Print rfa1 definition

In [12]:
rfa1.Print()
RooTFnBinding::fa1[ TFn={fa1=sin(x)/x} obs=(x3) ] = -0.0547936

Make plot frame in observable, TF1 binding function

In [13]:
frame3 = x3.frame(Title="TF1 bound as ROOT.RooFit function")
rfa1.plotOn(frame3)

c = ROOT.TCanvas("rf105_funcbinding", "rf105_funcbinding", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.6)
frame1.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.6)
frame2.Draw()
c.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
frame3.GetYaxis().SetTitleOffset(1.6)
frame3.Draw()

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

Draw all canvases

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