rf712_lagrangianmorphfit

Performing a simple fit with RooLagrangianMorphFunc

Author: Rahul Balasubramanian
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

ROOT.gStyle.SetOptStat(0)
ROOT.PyConfig.IgnoreCommandLineOptions = True
ROOT.gROOT.SetBatch(True)
Welcome to JupyROOT 6.27/01

Create functions


In [2]:
observablename = "pTV"
obsvar = ROOT.RooRealVar(observablename, "observable of pTV", 10, 600)

Setup three EFT coefficent and constant SM modifier

In [3]:
kSM = ROOT.RooRealVar("kSM", "sm modifier", 1.0)
cHq3 = ROOT.RooRealVar("cHq3", "EFT modifier", -10.0, 10.0)
cHq3.setAttribute("NewPhysics", True)
cHl3 = ROOT.RooRealVar("cHl3", "EFT modifier", -10.0, 10.0)
cHl3.setAttribute("NewPhysics", True)
cHDD = ROOT.RooRealVar("cHDD", "EFT modifier", -10.0, 10.0)
cHDD.setAttribute("NewPhysics", True)

Inputs to setup config


In [4]:
infilename = ROOT.gROOT.GetTutorialDir().Data() + "/roofit/input_histos_rf_lagrangianmorph.root"
par = "cHq3"
samplelist = [
    "SM_NPsq0",
    "cHq3_NPsq1",
    "cHq3_NPsq2",
    "cHl3_NPsq1",
    "cHl3_NPsq2",
    "cHDD_NPsq1",
    "cHDD_NPsq2",
    "cHl3_cHDD_NPsq2",
    "cHq3_cHDD_NPsq2",
    "cHl3_cHq3_NPsq2",
]

Set Config


In [5]:
config = ROOT.RooLagrangianMorphFunc.Config()
config.fileName = infilename
config.observableName = observablename
config.folderNames = samplelist
config.couplings.add(cHq3)
config.couplings.add(cHDD)
config.couplings.add(cHl3)
config.couplings.add(kSM)
Out[5]:
True

Create morphing function


In [6]:
morphfunc = ROOT.RooLagrangianMorphFunc("morphfunc", "morphed dist. of pTV", config)
[#0] PROGRESS:InputArguments -- initializing physics inputs from file /home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/src/master.build/tutorials/roofit/input_histos_rf_lagrangianmorph.root with object name(s) 'pTV'

Create pseudo data histogram to fit at cHq3 = 0.01, cHl3 = 1.0, cHDD = 0.2


In [7]:
morphfunc.setParameter("cHq3", 0.01)
morphfunc.setParameter("cHl3", 1.0)
morphfunc.setParameter("cHDD", 0.2)

pseudo_hist = morphfunc.createTH1("pseudo_hist")
pseudo_dh = ROOT.RooDataHist("pseudo_dh", "pseudo_dh", [obsvar], pseudo_hist)
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x95ff690
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
[#1] INFO:DataHandling -- RooDataHist::adjustBinning(pseudo_dh): fit range of variable pTV expanded to nearest bin boundaries: [10,600] --> [0,600]

reset parameters to zeros before fit

In [8]:
morphfunc.setParameter("cHq3", 0.0)
morphfunc.setParameter("cHl3", 0.0)
morphfunc.setParameter("cHDD", 0.0)

set error to set initial step size in fit

In [9]:
cHq3.setError(0.1)
cHl3.setError(0.1)
cHDD.setError(0.1)

Wrap pdf on morphfunc and fit to data histogram


wrapper pdf to normalise morphing function to a morphing pdf

In [10]:
model = ROOT.RooWrapperPdf("wrap_pdf", "wrap_pdf", morphfunc)
fitres = model.fitTo(pseudo_dh, SumW2Error=True, Optimize=False, Save=True)
[#0] PROGRESS:Caching -- creating cache from getCache function for 0xa182030
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
 **********
 **    1 **SET PRINT           1
 **********
 **********
 **    2 **SET NOGRAD
 **********
 PARAMETER DEFINITIONS:
    NO.   NAME         VALUE      STEP SIZE      LIMITS
     1 cHDD         0.00000e+00  1.00000e-01   -1.00000e+01  1.00000e+01
     2 cHl3         0.00000e+00  1.00000e-01   -1.00000e+01  1.00000e+01
     3 cHq3         0.00000e+00  1.00000e-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        1500           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=39609.1 FROM MIGRAD    STATUS=INITIATE       14 CALLS          15 TOTAL
                     EDM= unknown      STRATEGY= 1      NO ERROR MATRIX       
  EXT PARAMETER               CURRENT GUESS       STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  cHDD         0.00000e+00   1.00000e-01   1.00002e-02  -2.72888e-01
   2  cHl3         0.00000e+00   1.00000e-01   1.00002e-02  -1.83697e+01
   3  cHq3         0.00000e+00   1.00000e-01   1.00002e-02  -3.57539e+03
                               ERR DEF= 0.5
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39615.1) to force MIGRAD to back out of this region. Error log follows.
Parameter values: 	cHDD=0.223088	cHl3=2.02101	cHq3=0.00780325
RooWrapperPdf::wrap_pdf[ inputFunction=morphfunc ]
     p.d.f value is less than zero (-0.041358), trying to recover @ inputFunction=morphfunc=-0.0413585
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0413585
     p.d.f value is less than zero (-0.043016), trying to recover @ inputFunction=morphfunc=-0.043016
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.043016
RooNLLVar::nll_wrap_pdf_pseudo_dh[ parameters=(binWidth_pTV,cHDD,cHl3,cHq3,kSM,nNP0,nNP1,nNP2,nNP3,nNP4) ]
     function value is NAN @ parameters=(binWidth_pTV = 0.05,cHDD = 0.223088 +/- 0.1,cHl3 = 2.02101 +/- 0.1,cHq3 = 0.00780325 +/- 0.1,kSM = 1,nNP0 = 1,nNP1 = 1,nNP2 = 1,nNP3 = 1,nNP4 = 1)

RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39615.1) to force MIGRAD to back out of this region. Error log follows.
Parameter values: 	cHDD=0.56869	cHl3=4.96014	cHq3=0.0199009
RooWrapperPdf::wrap_pdf[ inputFunction=morphfunc ]
     p.d.f value is less than zero (-0.074541), trying to recover @ inputFunction=morphfunc=-0.0745407
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0745407
     p.d.f value is less than zero (-1.173323), trying to recover @ inputFunction=morphfunc=-1.17332
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-1.17332
     p.d.f value is less than zero (-1.063175), trying to recover @ inputFunction=morphfunc=-1.06317
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-1.06317
     p.d.f value is less than zero (-1.234351), trying to recover @ inputFunction=morphfunc=-1.23435
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-1.23435
     p.d.f value is less than zero (-1.674381), trying to recover @ inputFunction=morphfunc=-1.67438
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-1.67438
     p.d.f value is less than zero (-0.075079), trying to recover @ inputFunction=morphfunc=-0.0750793
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0750793
    ... (remaining 10 messages suppressed)
RooNLLVar::nll_wrap_pdf_pseudo_dh[ parameters=(binWidth_pTV,cHDD,cHl3,cHq3,kSM,nNP0,nNP1,nNP2,nNP3,nNP4) ]
     function value is NAN @ parameters=(binWidth_pTV = 0.05,cHDD = 0.56869 +/- 0.1,cHl3 = 4.96014 +/- 0.1,cHq3 = 0.0199009 +/- 0.1,kSM = 1,nNP0 = 1,nNP1 = 1,nNP2 = 1,nNP3 = 1,nNP4 = 1)

RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39615.1) to force MIGRAD to back out of this region. Error log follows.
Parameter values: 	cHDD=0.318765	cHl3=2.86725	cHq3=0.0111509
RooWrapperPdf::wrap_pdf[ inputFunction=morphfunc ]
     p.d.f value is less than zero (-0.070498), trying to recover @ inputFunction=morphfunc=-0.0704979
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0704979
     p.d.f value is less than zero (-0.063999), trying to recover @ inputFunction=morphfunc=-0.0639991
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0639991
     p.d.f value is less than zero (-0.140607), trying to recover @ inputFunction=morphfunc=-0.140607
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.140607
RooNLLVar::nll_wrap_pdf_pseudo_dh[ parameters=(binWidth_pTV,cHDD,cHl3,cHq3,kSM,nNP0,nNP1,nNP2,nNP3,nNP4) ]
     function value is NAN @ parameters=(binWidth_pTV = 0.05,cHDD = 0.318765 +/- 0.1,cHl3 = 2.86725 +/- 0.1,cHq3 = 0.0111509 +/- 0.1,kSM = 1,nNP0 = 1,nNP1 = 1,nNP2 = 1,nNP3 = 1,nNP4 = 1)

RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39615.1) to force MIGRAD to back out of this region. Error log follows.
Parameter values: 	cHDD=0.193718	cHl3=1.7579	cHq3=0.00677582
RooWrapperPdf::wrap_pdf[ inputFunction=morphfunc ]
     p.d.f value is less than zero (-0.029110), trying to recover @ inputFunction=morphfunc=-0.0291104
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0291104
     p.d.f value is less than zero (-0.011532), trying to recover @ inputFunction=morphfunc=-0.0115319
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0115319
RooNLLVar::nll_wrap_pdf_pseudo_dh[ parameters=(binWidth_pTV,cHDD,cHl3,cHq3,kSM,nNP0,nNP1,nNP2,nNP3,nNP4) ]
     function value is NAN @ parameters=(binWidth_pTV = 0.05,cHDD = 0.193718 +/- 0.1,cHl3 = 1.7579 +/- 0.1,cHq3 = 0.00677582 +/- 0.1,kSM = 1,nNP0 = 1,nNP1 = 1,nNP2 = 1,nNP3 = 1,nNP4 = 1)

RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39615.1) to force MIGRAD to back out of this region. Error log follows.
Parameter values: 	cHDD=0.820784	cHl3=1.67112	cHq3=0.0140305
RooWrapperPdf::wrap_pdf[ inputFunction=morphfunc ]
     p.d.f value is less than zero (-0.010720), trying to recover @ inputFunction=morphfunc=-0.0107196
     getLogVal() top-level p.d.f evaluates to NaN @ inputFunction=morphfunc=-0.0107196
RooNLLVar::nll_wrap_pdf_pseudo_dh[ parameters=(binWidth_pTV,cHDD,cHl3,cHq3,kSM,nNP0,nNP1,nNP2,nNP3,nNP4) ]
     function value is NAN @ parameters=(binWidth_pTV = 0.05,cHDD = 0.820784 +/- 0.1,cHl3 = 1.67112 +/- 0.1,cHq3 = 0.0140305 +/- 0.1,kSM = 1,nNP0 = 1,nNP1 = 1,nNP2 = 1,nNP3 = 1,nNP4 = 1)

 MIGRAD MINIMIZATION HAS CONVERGED.
 MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=39604.9 FROM MIGRAD    STATUS=CONVERGED      97 CALLS          98 TOTAL
                     EDM=4.90989e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                   STEP         FIRST   
  NO.   NAME      VALUE            ERROR          SIZE      DERIVATIVE 
   1  cHDD         1.98078e-01   1.34076e+00   1.76827e-02  -1.24714e-02
   2  cHl3         9.98258e-01   2.82061e-01   3.84123e-03  -2.35191e-01
   3  cHq3         1.00089e-02   4.07058e-03   5.39370e-05   6.27101e+00
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  3    ERR DEF=0.5
  1.809e+00 -5.353e-02 -1.409e-03 
 -5.353e-02  7.958e-02  1.262e-04 
 -1.409e-03  1.262e-04  1.657e-05 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3
        1  0.28129   1.000 -0.141 -0.257
        2  0.16034  -0.141  1.000  0.110
        3  0.26790  -0.257  0.110  1.000
 **********
 **    7 **SET ERR         0.5
 **********
 **********
 **    8 **SET PRINT           1
 **********
 **********
 **    9 **HESSE        1500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=39604.9 FROM HESSE     STATUS=OK             16 CALLS         114 TOTAL
                     EDM=4.93538e-05    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  cHDD         1.98078e-01   1.32824e+00   7.07309e-04   1.98091e-02
   2  cHl3         9.98258e-01   2.81948e-01   1.53649e-04   9.99923e-02
   3  cHq3         1.00089e-02   4.04044e-03   2.15748e-06   1.00089e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  3    ERR DEF=0.5
  1.775e+00 -4.527e-02 -1.221e-03 
 -4.527e-02  7.952e-02  1.208e-04 
 -1.221e-03  1.208e-04  1.633e-05 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3
        1  0.24668   1.000 -0.121 -0.227
        2  0.14509  -0.121  1.000  0.106
        3  0.24026  -0.227  0.106  1.000

[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf) Calculating sum-of-weights-squared correction matrix for covariance matrix
 **********
 **   10 **SET ERR         0.5
 **********
 **********
 **   11 **SET PRINT           1
 **********
 **********
 **   12 **HESSE        1500
 **********
 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
 FCN=40071 FROM HESSE     STATUS=OK             22 CALLS         136 TOTAL
                     EDM=4.93731    STRATEGY= 1      ERROR MATRIX ACCURATE 
  EXT PARAMETER                                INTERNAL      INTERNAL  
  NO.   NAME      VALUE            ERROR       STEP SIZE       VALUE   
   1  cHDD         1.98078e-01   1.36717e+00   1.44291e-02   1.98091e-02
   2  cHl3         9.98258e-01   1.81304e-01   3.13445e-03   9.99923e-02
   3  cHq3         1.00089e-02   4.47315e-03   4.40126e-05   1.00089e-03
                               ERR DEF= 0.5
 EXTERNAL ERROR MATRIX.    NDIM=  25    NPAR=  3    ERR DEF=0.5
  1.881e+00  6.681e-02 -2.510e-03 
  6.681e-02  3.287e-02 -3.029e-06 
 -2.510e-03 -3.029e-06  2.001e-05 
 PARAMETER  CORRELATION COEFFICIENTS  
       NO.  GLOBAL      1      2      3
        1  0.48867   1.000  0.269 -0.409
        2  0.29280   0.269  1.000 -0.004
        3  0.42378  -0.409 -0.004  1.000

run the fit Get the correlation matrix

In [11]:
hcorr = fitres.correlationHist()

Extract postfit distribution and plot with initial histogram


In [12]:
postfit_hist = morphfunc.createTH1("morphing_postfit_hist")
postfit_dh = ROOT.RooDataHist("morphing_postfit_dh", "morphing_postfit_dh", [obsvar], postfit_hist)

frame0 = obsvar.frame(Title="Input templates for p_{T}^{V}")
postfit_dh.plotOn(
    frame0,
    Name="postfit_dist",
    DrawOption="C",
    LineColor="b",
    DataError=None,
    XErrorSize=0,
)
pseudo_dh.plotOn(frame0, Name="input")
Out[12]:
<cppyy.gbl.RooPlot object at 0xa617930>
[#1] INFO:DataHandling -- RooDataHist::adjustBinning(morphing_postfit_dh): fit range of variable pTV expanded to nearest bin boundaries: [0,600] --> [0,600]
[#1] INFO:InputArguments -- RooAbsData::plotOn(pseudo_dh) INFO: dataset has non-integer weights, auto-selecting SumW2 errors instead of Poisson errors
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 7389.24 will supercede previous event count of 7396.07 for normalization of PDF projections

Draw plots on canvas


In [13]:
c1 = ROOT.TCanvas("fig3", "fig3", 800, 400)
c1.Divide(2, 1)

c1.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
ROOT.gPad.SetRightMargin(0.05)

model.paramOn(frame0, ROOT.RooFit.Layout(0.50, 0.75, 0.9))
frame0.GetXaxis().SetTitle("p_{T}^{V}")
frame0.Draw()

c1.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
ROOT.gPad.SetRightMargin(0.15)
ROOT.gStyle.SetPaintTextFormat("4.1f")
ROOT.gStyle.SetOptStat(0)
hcorr.SetMarkerSize(3.0)
hcorr.SetTitle("correlation matrix")
hcorr.GetYaxis().SetTitleOffset(1.4)
hcorr.GetYaxis().SetLabelSize(0.1)
hcorr.GetXaxis().SetLabelSize(0.1)
hcorr.GetYaxis().SetBinLabel(1, "c_{HDD}")
hcorr.GetYaxis().SetBinLabel(2, "c_{Hl^{(3)}}")
hcorr.GetYaxis().SetBinLabel(3, "c_{Hq^{(3)}}")
hcorr.GetXaxis().SetBinLabel(3, "c_{HDD}")
hcorr.GetXaxis().SetBinLabel(2, "c_{Hl^{(3)}}")
hcorr.GetXaxis().SetBinLabel(1, "c_{Hq^{(3)}}")
hcorr.GetYaxis().SetTitleOffset(1.4)
hcorr.Draw("colz text")

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

Draw all canvases

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