Organization and simultaneous fits: easy interactive access to workspace contents - CINT to CLING code migration 
 
  
import ROOT
 
 
def fillWorkspace(w):
    
    
 
    
    x = ROOT.RooRealVar("x", "x", 0, 10)
 
    
    
    mean = ROOT.RooRealVar("mean", "mean of gaussians", 5, 0, 10)
    sigma1 = ROOT.RooRealVar("sigma1", "width of gaussians", 0.5)
    sigma2 = ROOT.RooRealVar("sigma2", "width of gaussians", 1)
 
    sig1 = ROOT.RooGaussian("sig1", "Signal component 1", x, mean, sigma1)
    sig2 = ROOT.RooGaussian("sig2", "Signal component 2", x, mean, sigma2)
 
    
    a0 = ROOT.RooRealVar("a0", "a0", 0.5, 0.0, 1.0)
    a1 = ROOT.RooRealVar("a1", "a1", -0.2, 0.0, 1.0)
    bkg = ROOT.RooChebychev("bkg", "Background", x, [a0, a1])
 
    
    sig1frac = ROOT.RooRealVar("sig1frac", "fraction of component 1 in signal", 0.8, 0.0, 1.0)
    sig = ROOT.RooAddPdf("sig", "Signal", [sig1, sig2], [sig1frac])
 
    
    bkgfrac = ROOT.RooRealVar("bkgfrac", "fraction of background", 0.5, 0.0, 1.0)
    model = ROOT.RooAddPdf("model", "g1+g2+a", [bkg, sig], [bkgfrac])
 
    w.Import(model)
 
 
 
 
w = ROOT.RooWorkspace("w", True)
 
fillWorkspace(w)
 
w.Print()
 
 
 
 
model = w["model"]
x = w["x"]
 
d = model.generate({x}, 1000)
r = model.fitTo(d, PrintLevel=-1)
 
 
frame = x.frame()
d.plotOn(frame)
 
 
bkg = w["bkg"]
model.plotOn(frame)
model.plotOn(frame, Components=bkg, LineStyle="--")
 
c = ROOT.TCanvas("rf509_wsinteractive", "rf509_wsinteractive", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()
 
c.SaveAs("rf509_wsinteractive.png")
  [#0] WARNING:InputArguments -- The parameter 'sigma1' with range [-1e+30, 1e+30] of the RooGaussian 'sig1' exceeds the safe range of (0, inf). Advise to limit its range.
[#0] WARNING:InputArguments -- The parameter 'sigma2' with range [-1e+30, 1e+30] of the RooGaussian 'sig2' exceeds the safe range of (0, inf). Advise to limit its range.
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooChebychev::bkg
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::x
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a0
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::a1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::bkgfrac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooAddPdf::sig
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::mean
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma1
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sig1frac
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooGaussian::sig2
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooRealVar::sigma2
 
RooWorkspace(w) w contents
 
variables
---------
(a0,a1,bkgfrac,mean,sig1frac,sigma1,sigma2,x)
 
p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 1
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 1/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 1/1
RooGaussian::sig1[ x=x mean=mean sigma=sigma1 ] = 1
RooGaussian::sig2[ x=x mean=mean sigma=sigma2 ] = 1
 
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --  The following expressions will be evaluated in cache-and-track mode: (bkg,sig1,sig2)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (bkg)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
- Date
 - February 2018 
 
- Authors
 - Clemens Lange, Wouter Verkerke (C++ version) 
 
Definition in file rf509_wsinteractive.py.