ROOT
v6-30
Reference Guide
Loading...
Searching...
No Matches
rf512_wsfactory_oper.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
##
5
## 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #512
6
##
7
## Illustration of operator expressions and expression-based
8
## basic p.d.f.s in the workspace factory syntax
9
##
10
## \macro_code
11
## \macro_output
12
##
13
## \date February 2018
14
## \authors Clemens Lange, Wouter Verkerke (C version)
15
16
17
import
ROOT
18
19
20
w = ROOT.RooWorkspace(
"w"
)
21
22
# You can define typedefs for even shorter construction semantics
23
w.factory(
"$Typedef(Gaussian,Gaus)"
)
24
w.factory(
"$Typedef(Chebychev,Cheby)"
)
25
26
# Operator pdf examples
27
# ------------------------------------------------
28
29
# PDF addition is done with SUM (coef1*pdf1,pdf2)
30
w.factory(
"SUM::summodel( f[0,1]*Gaussian::gx(x[-10,10],m[0],1.0), Chebychev::ch(x,{0.1,0.2,-0.3}) )"
)
31
32
# Extended PDF addition is done with SUM (yield1*pdf1,yield2*pdf2)
33
w.factory(
"SUM::extsummodel( Nsig[0,1000]*gx, Nbkg[0,1000]*ch )"
)
34
35
# PDF multiplication is done with PROD ( pdf1, pdf2 )
36
w.factory(
"PROD::gxz( gx, Gaussian::gz(z[-10,10],0,1) )"
)
37
38
# Conditional p.d.f multiplication is done with PROD ( pdf1|obs, pdf2 )
39
w.factory(
"Gaussian::gy( y[-10,10], x, 1.0 )"
)
40
w.factory(
"PROD::gxycond( gy|x, gx )"
)
41
42
# Convolution (numeric/ fft) is done with NCONV/FCONV (obs,pdf1,pdf2)
43
w.factory(
"FCONV::lxg( x, Gaussian::g(x,mg[0],1), Landau::lc(x,0,1) )"
)
44
45
# Simultaneous p.d.f.s are constructed with SIMUL( index, state1=pdf1,
46
# state2=pdf2,...)
47
w.factory(
"SIMUL::smodel( c[A=0,B=1], A=Gaussian::gs(x,m,s[1.0, 0.01, 10.0]), B=Landau::ls(x,0,1) )"
)
48
49
# Operator function examples
50
# ---------------------------------------------------
51
52
# Function multiplication is done with prod (func1, func2,...)
53
w.factory(
"prod::uv(u[10],v[10])"
)
54
55
# Function addition is done with sum(func1,func2)
56
w.factory(
"sum::uv2(u,v)"
)
57
58
# Lagrangian morphing function for the example shown in rf711_lagrangianmorph
59
infilename = ROOT.gROOT.GetTutorialDir().Data() +
"/roofit/input_histos_rf_lagrangianmorph.root"
60
w.factory(
61
"lagrangianmorph::morph($observableName('pTV'),$fileName('"
62
+ infilename
63
+
"'),$couplings({cHq3[0,1],SM[1]}),$NewPhysics(cHq3=1,SM=0),$folders({'SM_NPsq0','cHq3_NPsq1','cHq3_NPsq2'}))"
64
)
65
66
# Taylor expansion is done with taylorexpand(func,{var1,var2,...},val,order)
67
w.factory(
"taylorexpand::te(expr::poly('x^4+5*x^3+2*x^2+x+1',x),{x},0.0,2)"
)
68
69
70
# Interpreted and compiled expression based pdfs
71
# ---------------------------------------------------------------------------------------------------
72
73
# Create a ROOT.RooGenericPdf interpreted p.d.f. You can use single quotes
74
# to pass the expression string argument
75
w.factory(
"EXPR::G('x*x+1',x)"
)
76
77
# Create a custom compiled p.d.f similar to the above interpreted p.d.f.
78
# The code required to make self p.d.f. is automatically embedded in
79
# the workspace
80
w.factory(
"CEXPR::GC('x*x+a',{x,a[1]})"
)
81
82
# Compiled and interpreted functions (rather than p.d.f.s) can be made with the lower case
83
# 'expr' and 'cexpr' types
84
85
# Print workspace contents
86
w.Print()
87
88
# Make workspace visible on command line
89
ROOT.gDirectory.Add(w)
tutorials
roofit
rf512_wsfactory_oper.py
ROOT v6-30 - Reference Guide Generated on Sun Dec 1 2024 05:18:21 (GVA Time) using Doxygen 1.9.8