rf509_wsinteractive

Organization and simultaneous fits: easy interactive access to workspace contents - CINT to CLING code migration

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 Wednesday, November 30, 2022 at 11:23 AM.

In [1]:
import ROOT


def fillWorkspace(w):
    # Create pdf and fill workspace
    # --------------------------------------------------------

    # Declare observable x
    x = ROOT.RooRealVar("x", "x", 0, 10)

    # Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and
    # their parameters
    mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
    sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
    sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)

    sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
    sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)

    # Build Chebychev polynomial pdf
    a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
    a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
    bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])

    # Sum the signal components into a composite signal pdf
    sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
    sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])

    # Sum the composite signal and background
    bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
    model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])

    w.Import(model)
Welcome to JupyROOT 6.27/01

Create and fill workspace

Create a workspace named 'w' With CINT w could exports its contents to a same-name C++ namespace in CINT 'namespace w'. but self does not work anymore in CLING. so self tutorial is an example on how to change the code

In [2]:
w = ROOT.RooWorkspace("w", True)

Fill workspace with pdf and data in a separate function

In [3]:
fillWorkspace(w)
[#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-1e+30, 1e+30] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-1e+30, 1e+30] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooChebychev::bkg
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::x
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a0
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::bkgfrac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::sig
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::mean
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sig1frac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig2
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma2

Print workspace contents

In [4]:
w.Print()
RooWorkspace(w) w contents

variables
---------
(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2,x)

p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 1
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 1/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 1/1
RooGaussian::sig1[ x=x mean=mean sigma=sigma1 ] = 1
RooGaussian::sig2[ x=x mean=mean sigma=sigma2 ] = 1

self does not work anymore with CLING use normal workspace functionality

Use workspace contents

Old syntax to use the name space prefix operator to access the workspace contents

d = w.model.generate(w.x,1000) r = w.model.fitTo(*d)

use normal workspace methods

In [5]:
model = w["model"]
x = w["x"]

d = model.generate({x}, 1000)
r = model.fitTo(d)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --  The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 a0           5.00000e-01  1.00000e-01    0.00000e+00  1.00000e+00
     2 a1           0.00000e+00  1.00000e-01    0.00000e+00  1.00000e+00
 MINUIT WARNING IN PARAM DEF
 ============== STARTING VALUE IS AT LIMIT.
 MINUIT WARNING IN PARAMETR
 ============== VARIABLE2 IS AT ITS LOWER ALLOWED LIMIT.
 MINUIT WARNING IN PARAMETR
 ============== VARIABLE2 BROUGHT BACK INSIDE LIMITS.
     3 bkgfrac      5.00000e-01  1.00000e-01    0.00000e+00  1.00000e+00
     4 mean         5.00000e+00  1.00000e+00    0.00000e+00  1.00000e+01
     5 sig1frac     8.00000e-01  1.00000e-01    0.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        2500           1
 **********
 FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
 MINUIT WARNING IN MIGrad    
 ============== VARIABLE2 IS AT ITS LOWER ALLOWED LIMIT.
 START MIGRAD MINIMIZATION.  STRATEGY  1.  CONVERGENCE WHEN EDM .LT. 1.00e-03
 FCN=1926.2 FROM MIGRAD    STATUS=INITIATE       37 CALLS          38 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a0           5.00000e-01   1.00000e-01   0.00000e+00   2.31850e+00
   2  a1           7.50936e-02   1.00000e-01   5.53213e-01  -4.13475e-03
   3  bkgfrac      5.00000e-01   1.00000e-01   0.00000e+00   2.67741e+01
   4  mean         5.00000e+00   1.00000e+00   0.00000e+00  -1.91257e+01
   5  sig1frac     8.00000e-01   1.00000e-01   0.00000e+00   3.58316e-01
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=1924.14 FROM MIGRAD    STATUS=CONVERGED     144 CALLS         145 TOTAL
                     EDM=4.87682e-06    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  a0           5.05217e-01   8.02459e-02   4.72375e-03  -6.18134e-03
   2  a1           2.66450e-01   1.34580e-01   5.73280e-03   1.14496e-03
   3  bkgfrac      4.27468e-01   3.57738e-02   1.24965e-03  -3.49141e-02
   4  mean         5.00772e+00   3.27931e-02   1.96316e-04  -2.49989e-02
   5  sig1frac     6.40178e-01   9.71885e-02   4.22823e-03   5.49851e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  5    ERR DEF=0.5
  6.495e-03  2.384e-03 -5.200e-04 -3.310e-04 -1.273e-03 
  2.384e-03  1.870e-02 -3.798e-03 -2.595e-05 -8.739e-03 
 -5.200e-04 -3.798e-03  1.282e-03 -4.006e-05  2.493e-03 
 -3.310e-04 -2.595e-05 -4.006e-05  1.075e-03 -1.593e-04 
 -1.273e-03 -8.739e-03  2.493e-03 -1.593e-04  9.578e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3      4      5
        1  0.25219   1.000  0.216 -0.180 -0.125 -0.161
        2  0.79261   0.216  1.000 -0.776 -0.006 -0.653
        3  0.82176  -0.180 -0.776  1.000 -0.034  0.711
        4  0.15132  -0.125 -0.006 -0.034  1.000 -0.050
        5  0.73023  -0.161 -0.653  0.711 -0.050  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        2500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=1924.14 FROM HESSE     STATUS=OK             31 CALLS         176 TOTAL
                     EDM=4.84164e-06    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  a0           5.05217e-01   8.02153e-02   9.44750e-04   1.04344e-02
   2  a1           2.66450e-01   1.33745e-01   2.29312e-04  -4.86007e-01
   3  bkgfrac      4.27468e-01   3.55504e-02   2.49929e-04  -1.45579e-01
   4  mean         5.00772e+00   3.27917e-02   3.92633e-05   1.54346e-03
   5  sig1frac     6.40178e-01   9.67194e-02   8.45646e-04   2.84164e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  5    ERR DEF=0.5
  6.490e-03  2.349e-03 -5.109e-04 -3.311e-04 -1.251e-03 
  2.349e-03  1.846e-02 -3.735e-03 -2.611e-05 -8.585e-03 
 -5.109e-04 -3.735e-03  1.266e-03 -3.996e-05  2.454e-03 
 -3.311e-04 -2.611e-05 -3.996e-05  1.075e-03 -1.575e-04 
 -1.251e-03 -8.585e-03  2.454e-03 -1.575e-04  9.484e-03 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3      4      5
        1  0.25076   1.000  0.215 -0.178 -0.125 -0.159
        2  0.78957   0.215  1.000 -0.773 -0.006 -0.649
        3  0.81926  -0.178 -0.773  1.000 -0.034  0.708
        4  0.15103  -0.125 -0.006 -0.034  1.000 -0.049
        5  0.72707  -0.159 -0.649  0.708 -0.049  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

old syntax to access the variable x frame = w.x.frame()

In [6]:
frame = x.frame()
d.plotOn(frame)
Out[6]:
<cppyy.gbl.RooPlot object at 0x91624d0>

OLD syntax to ommit x. NB: The 'w.' prefix can be omitted if namespace w is imported in local namespace in the usual C++ way

using namespace w model.plotOn(frame) model.plotOn(frame, Components=bkg, LineStyle="--")

correct syntax

In [7]:
bkg = w["bkg"]
model.plotOn(frame)
model.plotOn(frame, Components=bkg, LineStyle="--")
Out[7]:
<cppyy.gbl.RooPlot object at 0x91624d0>
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()

Draw the frame on the canvas

In [8]:
c = ROOT.TCanvas("rf509_wsinteractive", "rf509_wsinteractive", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()

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

Draw all canvases

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