Logo ROOT   6.16/01
Reference Guide
rf511_wsfactory_basic.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4##
5## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #511
6##
7## Basic use of the 'object factory' associated with a workspace
8## to rapidly build p.d.f.s functions and their parameter components
9##
10## \macro_code
11##
12## \date February 2018
13## \author Clemens Lange
14## \author Wouter Verkerke (C version)
15
16
17import ROOT
18
19
20compact=ROOT.kFALSE
21w = ROOT.RooWorkspace("w")
22
23# Creating and adding basic pdfs
24# ----------------------------------------------------------------
25
26# Remake example p.d.f. of tutorial rs502_wspacewrite.C:
27#
28# Basic p.d.f. construction: ClassName.ObjectName(constructor arguments)
29# Variable construction : VarName[x,xlo,xhi], VarName[xlo,xhi], VarName[x]
30# P.d.f. addition : SUM.ObjectName(coef1*pdf1,...coefM*pdfM,pdfN)
31#
32
33if not compact:
34 # Use object factory to build p.d.f. of tutorial rs502_wspacewrite
35 w.factory("Gaussian::sig1(x[-10,10],mean[5,0,10],0.5)")
36 w.factory("Gaussian::sig2(x,mean,1)")
37 w.factory("Chebychev::bkg(x,{a0[0.5,0.,1],a1[-0.2,0.,1.]})")
38 w.factory("SUM::sig(sig1frac[0.8,0.,1.]*sig1,sig2)")
39 w.factory("SUM::model(bkgfrac[0.5,0.,1.]*bkg,sig)")
40
41else:
42
43 # Use object factory to build p.d.f. of tutorial rs502_wspacewrite but
44 # - Contracted to a single line recursive expression,
45 # - Omitting explicit names for components that are not referred to explicitly later
46
47 w.factory("SUM::model(bkgfrac[0.5,0.,1.]*Chebychev::bkg(x[-10,10],{a0[0.5,0.,1],a1[-0.2,0.,1.]}), "
48 "SUM(sig1frac[0.8,0.,1.]*Gaussian(x,mean[5,0,10],0.5), Gaussian(x,mean,1)))")
49
50# Advanced pdf constructor arguments
51# ----------------------------------------------------------------
52#
53# P.d.f. constructor arguments may by any type of ROOT.RooAbsArg, also
54#
55# Double_t -. converted to ROOT.RooConst(...)
56# {a,b,c} -. converted to ROOT.RooArgSet() or ROOT.RooArgList() depending on required ctor arg
57# dataset name -. convered to ROOT.RooAbsData reference for any dataset residing in the workspace
58# enum -. Any enum label that belongs to an enum defined in the (base)
59# class
60
61# Make a dummy dataset p.d.f. 'model' and import it in the workspace
62data = w.pdf("model").generate(ROOT.RooArgSet(w.var("x")), 1000)
63getattr(w, 'import')(data, ROOT.RooFit.Rename("data"))
64
65# Construct a KEYS p.d.f. passing a dataset name and an enum type defining the
66# mirroring strategy
67# w.factory("KeysPdf::k(x,data,NoMirror,0.2)")
68# Workaround for pyROOT
69x = w.var("x")
70k = ROOT.RooKeysPdf("k", "k", x, data, ROOT.RooKeysPdf.NoMirror, 0.2)
71getattr(w, 'import')(k, ROOT.RooFit.Rename("k"))
72
73# Print workspace contents
74w.Print()
75
76# Make workspace visible on command line
77ROOT.gDirectory.Add(w)