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