import ROOT
x = ROOT.RooRealVar("x", "x", -20, 20)
g1mean = ROOT.RooRealVar("g1mean", "g1mean", -10)
g1 = ROOT.RooGaussian("g1", "g1", x, g1mean, 2.0)
g2 = ROOT.RooPolynomial("g2", "g2", x, [-0.03, -0.001])
alpha = ROOT.RooRealVar("alpha", "alpha", 0, 1.0)
x.setBins(1000, "cache")
alpha.setBins(50, "cache")
lmorph = ROOT.RooIntegralMorph("lmorph", "lmorph", g1, g2, x, alpha)
frame1 = x.frame()
g1.plotOn(frame1)
g2.plotOn(frame1)
alpha.setVal(0.125)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.25)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.375)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.50)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.625)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.75)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.875)
lmorph.plotOn(frame1, LineColor="r")
alpha.setVal(0.95)
lmorph.plotOn(frame1, LineColor="r")
hh = lmorph.createHistogram("hh", x, Binning=40, YVar=dict(var=alpha, Binning=40))
hh.SetLineColor(ROOT.kBlue)
alpha.setVal(0.8)
data = lmorph.generate({x}, 1000)
lmorph.setCacheAlpha(True)
lmorph.fitTo(data, Verbose=True, PrintLevel=-1)
frame2 = x.frame(Bins=100)
data.plotOn(frame2)
lmorph.plotOn(frame2)
frame3 = alpha.frame(Bins=100, Range=(0.1, 0.9))
nll = lmorph.createNLL(data)
nll.plotOn(frame3, ShiftToZero=True)
lmorph.setCacheAlpha(False)
c = ROOT.TCanvas("rf705_linearmorph", "rf705_linearmorph", 800, 800)
c.Divide(2, 2)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.6)
frame1.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.20)
hh.GetZaxis().SetTitleOffset(2.5)
hh.Draw("surf")
c.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
frame3.GetYaxis().SetTitleOffset(1.4)
frame3.Draw()
c.cd(4)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.4)
frame2.Draw()
c.SaveAs("rf705_linearmorph.png")
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8bddf40 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x8c083e0 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x_alpha for nset (x,alpha) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x91bacc0 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x90b1580 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0 from preexisting content.
[#0] PROGRESS:Eval -- RooIntegralMorph::fillCacheObject(lmorph) filling multi-dimensional cache..................................................
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x937f4b0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0
[#1] INFO:Fitting -- RooAbsPdf::fitTo(lmorph_over_lmorph_Int[x]) 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_lmorph_over_lmorph_Int[x]_lmorphData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for alpha: using 0.1
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x937e260 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x] for nset () with code 1 from preexisting content.
prevFCN = 9770.913047 alpha=0.807,
prevFCN = 9770.181345 alpha=0.7929,
prevFCN = 9772.015076 alpha=0.8008,
prevFCN = 9770.800603 alpha=0.7992,
prevFCN = 9770.785966 alpha=0.8001,
prevFCN = 9770.896896 alpha=0.7999,
prevFCN = 9770.682007
prevFCN = 9770.682007 alpha=0.7994,
prevFCN = 9770.748743 alpha=0.7997,
prevFCN = 9770.712263 alpha=0.8004,
prevFCN = 9770.849353 alpha=0.8002,
prevFCN = 9770.888877 alpha=0.8,
prevFCN = 9770.908995 alpha=0.8,
prevFCN = 9770.671799 alpha=0.8,
prevFCN = 9770.909805 alpha=0.8,
prevFCN = 9770.667124 alpha=0.8,
prevFCN = 9770.910615 alpha=0.8,
prevFCN = 9770.91254 alpha=0.8,
prevFCN = 9770.912239 alpha=0.8,
prevFCN = 9770.669355 alpha=0.8,
prevFCN = 9770.666818 alpha=0.8,
prevFCN = 9770.66743 alpha=0.8,
prevFCN = 9770.91141 alpha=0.8,
prevFCN = 9770.670187 alpha=0.8,
prevFCN = 9770.668617 alpha=0.8,
prevFCN = 9770.667852 alpha=0.8,
prevFCN = 9770.667479 alpha=0.8,
prevFCN = 9770.667297 alpha=0.8,
prevFCN = 9770.667208 alpha=0.8,
prevFCN = 9770.667165 alpha=0.8,
prevFCN = 9770.667144 alpha=0.8,
prevFCN = 9770.667134 alpha=0.8,
prevFCN = 9770.667129 alpha=0.8,
prevFCN = 9770.667126 alpha=0.8,
prevFCN = 9770.667124 alpha=0.8,
prevFCN = 9770.91141 alpha=0.8,
prevFCN = 9770.670187 alpha=0.8,
prevFCN = 9770.666827 alpha=0.8,
prevFCN = 9770.667421 alpha=0.8,
prevFCN = 9770.667064 alpha=0.8,
prevFCN = 9770.667183 alpha=0.8,
prevFCN = 9770.667124 alpha=0.8,
prevFCN = 9770.667064 alpha=0.8,
prevFCN = 9770.667183 alpha=0.8,
prevFCN = 9770.66653 alpha=0.8,
prevFCN = 9770.667718 alpha=0.8,
prevFCN = 9770.667005 alpha=0.8,
prevFCN = 9770.667243 alpha=0.8, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x2744ab0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x95fa7c0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Fitting -- RooAbsPdf::fitTo(lmorph_over_lmorph_Int[x]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x9805940 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x] for nset () with code 1 from preexisting content.