import ROOT
ROOT.gStyle.SetOptStat(0)
ROOT.PyConfig.IgnoreCommandLineOptions = True
ROOT.gROOT.SetBatch(True)
observablename = "pTV"
obsvar = ROOT.RooRealVar(observablename, "observable of pTV", 10, 600)
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)
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",
]
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)
morphfunc = ROOT.RooLagrangianMorphFunc("morphfunc", "morphed dist. of pTV", config)
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)
morphfunc.setParameter("cHq3", 0.0)
morphfunc.setParameter("cHl3", 0.0)
morphfunc.setParameter("cHDD", 0.0)
cHq3.setError(0.1)
cHl3.setError(0.1)
cHDD.setError(0.1)
model = ROOT.RooWrapperPdf("wrap_pdf", "wrap_pdf", morphfunc)
fitres = model.fitTo(pseudo_dh, SumW2Error=True, Optimize=False, Save=True, PrintLevel=-1)
hcorr = fitres.correlationHist()
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")
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-v628/rootspi/rdoc/src/v6-28-00-patches.build/tutorials/roofit/input_histos_rf_lagrangianmorph.root with object name(s) 'pTV'
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x7d6dd60
[#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 0x8dfa7e0
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
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
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)
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
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
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)
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)
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
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)
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
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
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)
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
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
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)
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
[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf) Calculating sum-of-weights-squared correction matrix for covariance matrix
[#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