Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf712_lagrangianmorphfit.py File Reference

Detailed Description

View in nbviewer Open in SWAN
Performing a simple fit with RooLagrangianMorphFunc

import ROOT
ROOT.gStyle.SetOptStat(0)
ROOT.PyConfig.IgnoreCommandLineOptions = True
ROOT.gROOT.SetBatch(True)
# Create functions
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
observablename = "pTV"
obsvar = ROOT.RooRealVar(observablename, "observable of pTV", 10, 600)
# Setup three EFT coefficient and constant SM modifier
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
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
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
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
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)
# Create morphing function
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
morphfunc = ROOT.RooLagrangianMorphFunc("morphfunc", "morphed dist. of pTV", config)
# Create pseudo data histogram to fit at cHq3 = 0.01, cHl3 = 1.0, cHDD = 0.2
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
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)
# reset parameters to zeros before fit
morphfunc.setParameter("cHq3", 0.0)
morphfunc.setParameter("cHl3", 0.0)
morphfunc.setParameter("cHDD", 0.0)
# set error to set initial step size in fit
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
model = ROOT.RooWrapperPdf("wrap_pdf", "wrap_pdf", morphfunc)
fitres = model.fitTo(pseudo_dh, SumW2Error=True, Optimize=False, Save=True, PrintLevel=-1)
# run the fit
# Get the correlation matrix
hcorr = fitres.correlationHist()
# Extract postfit distribution and plot with initial histogram
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
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")
# Draw plots on canvas
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
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")
[#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'
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x860fda0
[#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]
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x9371280
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf_over_wrap_pdf_Int[pTV]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_wrap_pdf_over_wrap_pdf_Int[pTV]_pseudo_dh) Summation contains a RooNLLVar, using its error level
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0 cHl3=0 cHq3=-0.0202918
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.213672 cHl3=1.97898 cHq3=0.00773174
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.0861399 cHl3=-9.50561 cHq3=0.0801661
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.186765 cHl3=8.8591 cHq3=-0.971282
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-6.32705, denominator=wrap_pdf_Int[pTV]=46316
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.218731 cHl3=0.37397 cHq3=-2.08166
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=3.22694 cHl3=-7.04051 cHq3=0.54016
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-2.21831, denominator=wrap_pdf_Int[pTV]=89722.5
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=3.44258 cHl3=4.96668 cHq3=0.0273884
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
... (remaining 14 messages suppressed)
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=-3.35967 cHl3=-9.58704 cHq3=-6.27461
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.154263 cHl3=2.95902 cHq3=-2.78828
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=59.1285, denominator=wrap_pdf_Int[pTV]=200921
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.338546 cHl3=0.879879 cHq3=-1.35856
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.372361 cHl3=0.491134 cHq3=-0.886807
RooAbsPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=5.8312, denominator=wrap_pdf_Int[pTV]=12183.6
[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf) Calculating sum-of-weights-squared correction matrix for covariance matrix
[#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 supersede previous event count of 9313.81 for normalization of PDF projections
Date
January 2022
Author
Rahul Balasubramanian

Definition in file rf712_lagrangianmorphfit.py.