rf703_effpdfprod

Special pdf's: using a product of an (acceptance) efficiency and a pdf as pdf

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

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

Define observables and decay pdf

Declare observables

In [2]:
t = ROOT.RooRealVar("t", "t", 0, 5)

Make pdf

In [3]:
tau = ROOT.RooRealVar("tau", "tau", -1.54, -4, -0.1)
model = ROOT.RooExponential("model", "model", t, tau)

Define efficiency function

Use error function to simulate turn-on slope

In [4]:
eff = ROOT.RooFormulaVar("eff", "0.5*(TMath::Erf((t-1)/0.5)+1)", [t])

Define decay pdf with efficiency

Multiply pdf(t) with efficiency in t

In [5]:
modelEff = ROOT.RooEffProd("modelEff", "model with efficiency", model, eff)

Plot efficiency, pdf

In [6]:
frame1 = t.frame(Title="Efficiency")
eff.plotOn(frame1, LineColor="r")

frame2 = t.frame(Title="Pdf with and without efficiency")

model.plotOn(frame2, LineStyle="--")
modelEff.plotOn(frame2)
Out[6]:
<cppyy.gbl.RooPlot object at 0x8d44190>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)

Generate toy data, fit model eff to data

Generate events. If the input pdf has an internal generator, internal generator is used and an accept/reject sampling on the efficiency is applied.

In [7]:
data = modelEff.generate({t}, 10000)

Fit pdf. The normalization integral is calculated numerically.

In [8]:
modelEff.fitTo(data)
Out[8]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#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: (eff)
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 tau         -1.54000e+00  3.90000e-01   -4.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         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=9696.74 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  tau         -1.54000e+00   3.90000e-01   2.09076e-01   1.80888e+02
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=9695.84 FROM MIGRAD    STATUS=CONVERGED      12 CALLS          13 TOTAL
                     EDM=2.41592e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  tau         -1.55887e+00   1.40931e-02   5.05936e-04   6.58163e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  1.986e-04 
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE         500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=9695.84 FROM HESSE     STATUS=OK              5 CALLS          18 TOTAL
                     EDM=2.41576e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  tau         -1.55887e+00   1.40931e-02   1.01187e-04   2.54601e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  1    ERR DEF=0.5
  1.986e-04 
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot generated data and overlay fitted pdf

In [9]:
frame3 = t.frame(Title="Fitted pdf with efficiency")
data.plotOn(frame3)
modelEff.plotOn(frame3)

c = ROOT.TCanvas("rf703_effpdfprod", "rf703_effpdfprod", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.4)
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("rf703_effpdfprod.png")
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
Info in <TCanvas::Print>: png file rf703_effpdfprod.png has been created

Draw all canvases

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