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

Detailed Description

View in nbviewer Open in SWAN
Organization and simultaneous fits: basic use of the 'object factory' associated with a workspace to rapidly build pdfs functions and their parameter components

import ROOT
compact = False
w = ROOT.RooWorkspace("w")
# Creating and adding basic pdfs
# ----------------------------------------------------------------
# Remake example pdf of tutorial rs502_wspacewrite.C:
#
# Basic pdf construction: ClassName.ObjectName(constructor arguments)
# Variable construction : VarName[x,xlo,xhi], VarName[xlo,xhi], VarName[x]
# P.d.f. addition : SUM.ObjectName(coef1*pdf1,...coefM*pdfM,pdfN)
#
if not compact:
# Use object factory to build pdf of tutorial rs502_wspacewrite
w.factory("Gaussian::sig1(x[-10,10],mean[5,0,10],0.5)")
w.factory("Gaussian::sig2(x,mean,1)")
w.factory("Chebychev::bkg(x,{a0[0.5,0.,1],a1[-0.2,0.,1.]})")
w.factory("SUM::sig(sig1frac[0.8,0.,1.]*sig1,sig2)")
w.factory("SUM::model(bkgfrac[0.5,0.,1.]*bkg,sig)")
else:
# Use object factory to build pdf of tutorial rs502_wspacewrite but
# - Contracted to a single line recursive expression,
# - Omitting explicit names for components that are not referred to explicitly later
w.factory(
"SUM::model(bkgfrac[0.5,0.,1.]*Chebychev::bkg(x[-10,10],{a0[0.5,0.,1],a1[-0.2,0.,1.]}), "
"SUM(sig1frac[0.8,0.,1.]*Gaussian(x,mean[5,0,10],0.5), Gaussian(x,mean,1)))"
)
# Advanced pdf constructor arguments
# ----------------------------------------------------------------
#
# P.d.f. constructor arguments may by any type of ROOT.RooAbsArg, also
#
# Double_t -. converted to ROOT.RooConst(...)
# {a,b,c} -. converted to ROOT.RooArgSet() or ROOT.RooArgList() depending on required ctor arg
# dataset name -. converted to ROOT.RooAbsData reference for any dataset residing in the workspace
# enum -. Any enum label that belongs to an enum defined in the (base)
# class
# Make a dummy dataset pdf 'model' and import it in the workspace
data = w["model"].generate({w["x"]}, 1000)
# Cannot call 'import' directly because this is a python keyword:
w.Import(data, Rename="data")
# Construct a KEYS pdf passing a dataset name and an enum type defining the
# mirroring strategy
# w.factory("KeysPdf::k(x,data,NoMirror,0.2)")
# Workaround for pyROOT
x = w["x"]
k = ROOT.RooKeysPdf("k", "k", x, data, ROOT.RooKeysPdf.NoMirror, 0.2)
w.Import(k, RenameAllNodes="workspace")
# Print workspace contents
w.Print()
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing dataset modelData
[#1] INFO:ObjectHandling -- RooWorkSpace::import(w) changing name of dataset from modelData to data
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) Resolving name conflict in workspace by changing name of imported node k to k_workspace
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing RooKeysPdf::k_workspace
RooWorkspace(w) w contents
variables
---------
(a0,a1,bkgfrac,mean,sig1frac,x)
p.d.f.s
-------
RooChebychev::bkg[ x=x coefList=(a0,a1) ] = 1
RooKeysPdf::k_workspace[ x=x ] = 0.0139016
RooAddPdf::model[ bkgfrac * bkg + [%] * sig ] = 0.5/1
RooAddPdf::sig[ sig1frac * sig1 + [%] * sig2 ] = 7.45331e-07/1
RooGaussian::sig1[ x=x mean=mean sigma=0.5 ] = 1.92875e-22
RooGaussian::sig2[ x=x mean=mean sigma=1 ] = 3.72665e-06
datasets
--------
RooDataSet::data(x)
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf511_wsfactory_basic.py.