rf210_angularconv

Convolution in cyclical angular observables theta, and construction of p.d.f in terms of trasnformed angular coordinates, e.g. cos(theta), the convolution is performed in theta rather than cos(theta)

(require ROOT to be compiled with --enable-fftw3)

pdf(theta) = ROOT.T(theta) (x) gauss(theta) pdf(cosTheta) = ROOT.T(acos(cosTheta)) (x) gauss(acos(cosTheta))

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

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

Set up component pdfs

Define angle psi

In [2]:
psi = ROOT.RooRealVar("psi", "psi", 0, 3.14159268)

Define physics p.d.f T(psi)

In [3]:
Tpsi = ROOT.RooGenericPdf("Tpsi", "1+sin(2*@0)", [psi])

Define resolution R(psi)

In [4]:
gbias = ROOT.RooRealVar("gbias", "gbias", 0.2, 0.0, 1)
greso = ROOT.RooRealVar("greso", "greso", 0.3, 0.1, 1.0)
Rpsi = ROOT.RooGaussian("Rpsi", "Rpsi", psi, gbias, greso)

Define cos(psi) and function psif that calculates psi from cos(psi)

In [5]:
cpsi = ROOT.RooRealVar("cpsi", "cos(psi)", -1, 1)
psif = ROOT.RooFormulaVar("psif", "acos(cpsi)", [cpsi])

Define physics p.d.f. also as function of cos(psi): T(psif(cpsi)) = T(cpsi)

In [6]:
Tcpsi = ROOT.RooGenericPdf("T", "1+sin(2*@0)", [psif])

Construct convolution pdf in psi

Define convoluted p.d.f. as function of psi: M=T(x)R = M(psi)

In [7]:
Mpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psi, Tpsi, Rpsi)
[#1] INFO:Caching -- Changing internal binning of variable 'psi' in FFT 'Mf' from 100 to 930 to improve the precision of the numerical FFT. This can be done manually by setting an additional binning named 'cache'.

Set the buffer fraction to zero to obtain a ROOT.True cyclical convolution

In [8]:
Mpsi.setBufferFraction(0)

Sample, fit and plot convoluted pdf (psi)

Generate some events in observable psi

In [9]:
data_psi = Mpsi.generate({psi}, 10000)
[#1] INFO:Eval -- RooRealVar::setRange(psi) new range named 'refrange_fft_Mf' created with bounds [0,3.14159]
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x8418820 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x88456c0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0 from preexisting content.

Fit convoluted model as function of angle psi

In [10]:
Mpsi.fitTo(data_psi)
Out[10]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x8702e80 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0 from preexisting content.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Minimization --  The following expressions have been identified as constant and will be precalculated and cached: (Tpsi)
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 gbias        2.00000e-01  1.00000e-01    0.00000e+00  1.00000e+00
     2 greso        3.00000e-01  9.00000e-02    1.00000e-01  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        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
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
 FCN=9480.16 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  gbias        2.00000e-01   1.00000e-01   2.57889e-01   5.40962e+01
   2  greso        3.00000e-01   9.00000e-02   2.46497e-01   5.76705e+00
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=9479.64 FROM MIGRAD    STATUS=CONVERGED      32 CALLS          33 TOTAL
                     EDM=1.07647e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  gbias        1.92485e-01   7.44189e-03   1.26877e-03   1.36618e-02
   2  greso        2.98582e-01   9.19693e-03   1.65656e-03  -8.25561e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  5.539e-05  1.690e-07 
  1.690e-07  8.460e-05 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.00247   1.000  0.002
        2  0.00247   0.002  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=9479.64 FROM HESSE     STATUS=OK             10 CALLS          43 TOTAL
                     EDM=1.07707e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  gbias        1.92485e-01   7.44188e-03   5.07510e-05  -6.62425e-01
   2  greso        2.98582e-01   9.19691e-03   6.62625e-05  -5.92825e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  5.539e-05  1.248e-07 
  1.248e-07  8.460e-05 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.00182   1.000  0.002
        2  0.00182   0.002  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot cos(psi) frame with Mf(cpsi)

In [11]:
frame1 = psi.frame(Title="Cyclical convolution in angle psi")
data_psi.plotOn(frame1)
Mpsi.plotOn(frame1)
Out[11]:
<cppyy.gbl.RooPlot object at 0x8a0cac0>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x8a78730 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0

Overlay comparison to unsmeared physics p.d.f ROOT.T(psi)

In [12]:
Tpsi.plotOn(frame1, LineColor="r")
Out[12]:
<cppyy.gbl.RooPlot object at 0x8a0cac0>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)

Construct convolution pdf in cos(psi)

Define convoluted p.d.f. as function of cos(psi): M=T(x)R) = M(cpsi:

Need to give both observable psi here (for definition of convolution) and function psif here (for definition of observables, in cpsi)

In [13]:
Mcpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psif, psi, Tpsi, Rpsi)

Set the buffer fraction to zero to obtain a ROOT.True cyclical convolution

In [14]:
Mcpsi.setBufferFraction(0)

Sample, fit and plot convoluted pdf (cospsi)

Generate some events

In [15]:
data_cpsi = Mcpsi.generate({cpsi}, 10000)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x8a78e30 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x8a415e0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)

set psi constant to exclude to be a parameter of the fit

In [16]:
psi.setConstant(True)

Fit convoluted model as function of cos(psi)

In [17]:
Mcpsi.fitTo(data_cpsi)
Out[17]:
<cppyy.gbl.RooFitResult object at 0x(nil)>
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x86b9f70 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#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: (psif)
 **********
 **   10 **SET PRINT           1
 **********
 **********
 **   11 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 gbias        1.92485e-01  7.44188e-03    0.00000e+00  1.00000e+00
     2 greso        2.98582e-01  9.19691e-03    1.00000e-01  1.00000e+00
 **********
 **   12 **SET ERR         0.5
 **********
 **********
 **   13 **SET PRINT           1
 **********
 **********
 **   14 **SET STR           1
 **********
 NOW USING STRATEGY  1: TRY TO BALANCE SPEED AGAINST RELIABILITY
 **********
 **   15 **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
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
 FCN=5178.91 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  gbias        1.92485e-01   7.44188e-03   1.88791e-02   2.79436e+01
   2  greso        2.98582e-01   9.19691e-03   2.46483e-02  -8.41273e+00
                               ERR DEF= 0.5
 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=5178.63 FROM MIGRAD    STATUS=CONVERGED      28 CALLS          29 TOTAL
                     EDM=3.37796e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  gbias        1.85814e-01   9.36461e-03   1.16493e-03   2.44691e-02
   2  greso        3.02469e-01   1.02714e-02   1.32671e-03   1.52048e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  8.771e-05 -2.174e-05 
 -2.174e-05  1.055e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.22594   1.000 -0.226
        2  0.22594  -0.226  1.000
 **********
 **   16 **SET ERR         0.5
 **********
 **********
 **   17 **SET PRINT           1
 **********
 **********
 **   18 **HESSE        1000
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=5178.63 FROM HESSE     STATUS=OK             10 CALLS          39 TOTAL
                     EDM=3.36905e-07    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  gbias        1.85814e-01   9.36561e-03   2.32986e-04  -6.79459e-01
   2  greso        3.02469e-01   1.02725e-02   5.30684e-05  -5.82446e-01
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  2    ERR DEF=0.5
  8.773e-05 -2.179e-05 
 -2.179e-05  1.056e-04 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2
        1  0.22639   1.000 -0.226
        2  0.22639  -0.226  1.000
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization

Plot cos(psi) frame with Mf(cpsi)

In [18]:
frame2 = cpsi.frame(Title="Same convolution in psi, in cos(psi)")
data_cpsi.plotOn(frame2)
Mcpsi.plotOn(frame2)
Out[18]:
<cppyy.gbl.RooPlot object at 0x8b63b20>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x86b9f70 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)

Overlay comparison to unsmeared physics p.d.f ROOT.Tf(cpsi)

In [19]:
Tcpsi.plotOn(frame2, LineColor="r")
Out[19]:
<cppyy.gbl.RooPlot object at 0x8b63b20>
[#1] INFO:NumericIntegration -- RooRealIntegral::init(T_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)

Draw frame on canvas

In [20]:
c = ROOT.TCanvas("rf210_angularconv", "rf210_angularconv", 800, 400)
c.Divide(2)
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.4)
frame2.Draw()

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

Draw all canvases

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