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
# 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/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.fileName = infilename
config.observableName = observablename
config.folderNames = samplelist
# 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
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
pseudo_hist = morphfunc.createTH1("pseudo_hist")
pseudo_dh = ROOT.RooDataHist("pseudo_dh", "pseudo_dh", [obsvar], pseudo_hist)
# reset parameters to zeros before fit
# set error to set initial step size in fit
# 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
# 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}")
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)
model.paramOn(frame0, ROOT.RooFit.Layout(0.50, 0.75, 0.9))
frame0.GetXaxis().SetTitle("p_{T}^{V}")
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")
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
[#0] PROGRESS:InputArguments -- initializing physics inputs from file /github/home/ROOT-CI/build/tutorials/roofit/roofit/input_histos_rf_lagrangianmorph.root with object name(s) 'pTV'
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x5560325d7580
[#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 0x55603339af20
[#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 generic CPU library compiled with no vectorizations
[#1] INFO:Fitting -- Creation of NLL object took 8.80133 ms
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_wrap_pdf_over_wrap_pdf_Int[pTV]_pseudo_dh) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- [fitFCN] No discrete parameters, performing continuous minimization only
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.22238 cHl3=2.00088 cHq3=0.00780325
RooFit::Detail::RooNormalizedPdf::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.0889525, denominator=wrap_pdf_Int[pTV]=4479.59
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0889525, denominator=wrap_pdf_Int[pTV]=4479.59
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0889525, denominator=wrap_pdf_Int[pTV]=4479.59
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0889525, denominator=wrap_pdf_Int[pTV]=4479.59
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=1.16153 cHl3=8.6965 cHq3=0.0408466
RooFit::Detail::RooNormalizedPdf::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.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.894368, denominator=wrap_pdf_Int[pTV]=3980.01
... (remaining 4 messages suppressed)
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.607791 cHl3=5.23432 cHq3=0.0213386
RooFit::Detail::RooNormalizedPdf::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.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.315769, denominator=wrap_pdf_Int[pTV]=237.773
... (remaining 16 messages suppressed)
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.330108 cHl3=2.94627 cHq3=0.0115846
RooFit::Detail::RooNormalizedPdf::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.131032, denominator=wrap_pdf_Int[pTV]=2465.37
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.131032, denominator=wrap_pdf_Int[pTV]=2465.37
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.131032, denominator=wrap_pdf_Int[pTV]=2465.37
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.131032, denominator=wrap_pdf_Int[pTV]=2465.37
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.131032, denominator=wrap_pdf_Int[pTV]=2465.37
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.131032, denominator=wrap_pdf_Int[pTV]=2465.37
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.191157 cHl3=1.72296 cHq3=0.00670752
RooFit::Detail::RooNormalizedPdf::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.0801556, denominator=wrap_pdf_Int[pTV]=5193.38
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0801556, denominator=wrap_pdf_Int[pTV]=5193.38
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0801556, denominator=wrap_pdf_Int[pTV]=5193.38
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0801556, denominator=wrap_pdf_Int[pTV]=5193.38
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.138238 cHl3=1.24892 cHq3=0.00485049
RooFit::Detail::RooNormalizedPdf::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=0.0688277, denominator=wrap_pdf_Int[pTV]=6538.46
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (39620.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.121667 cHl3=1.24299 cHq3=0.004269
RooFit::Detail::RooNormalizedPdf::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=0.0667251, denominator=wrap_pdf_Int[pTV]=6553.34
[#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 7388.38 for normalization of PDF projections
Date
January 2022
Author
Rahul Balasubramanian

Definition in file rf712_lagrangianmorphfit.py.