rf512_wsfactory_oper

'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #512

Illustration of operator expressions and expression-based basic p.d.f.s in the workspace factory syntax

Author: Clemens Lange, Wouter Verkerke (C version)
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Sunday, November 27, 2022 at 11:07 AM.

In [1]:
import ROOT


w = ROOT.RooWorkspace("w")
Welcome to JupyROOT 6.27/01

You can define typedefs for even shorter construction semantics

In [2]:
w.factory("$Typedef(Gaussian,Gaus)")
w.factory("$Typedef(Chebychev,Cheby)")
Out[2]:
<cppyy.gbl.RooAbsArg object at 0x(nil)>

Operator pdf examples

PDF addition is done with SUM (coef1*pdf1,pdf2)

In [3]:
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}) )")
Out[3]:
<cppyy.gbl.RooAddPdf object at 0x8f72530>

Extended PDF addition is done with SUM (yield1pdf1,yield2pdf2)

In [4]:
w.factory("SUM::extsummodel( Nsig[0,1000]*gx, Nbkg[0,1000]*ch )")
Out[4]:
<cppyy.gbl.RooAddPdf object at 0x914dc10>

PDF multiplication is done with PROD ( pdf1, pdf2 )

In [5]:
w.factory("PROD::gxz( gx, Gaussian::gz(z[-10,10],0,1) )")
Out[5]:
<cppyy.gbl.RooProdPdf object at 0x91a0ac0>

Conditional p.d.f multiplication is done with PROD ( pdf1|obs, pdf2 )

In [6]:
w.factory("Gaussian::gy( y[-10,10], x, 1.0 )")
w.factory("PROD::gxycond( gy|x, gx )")
Out[6]:
<cppyy.gbl.RooProdPdf object at 0x7b66760>

Convolution (numeric/ fft) is done with NCONV/FCONV (obs,pdf1,pdf2)

In [7]:
w.factory("FCONV::lxg( x, Gaussian::g(x,mg[0],1), Landau::lc(x,0,1) )")
Out[7]:
<cppyy.gbl.RooFFTConvPdf object at 0x924b300>
[#1] INFO:Caching -- Changing internal binning of variable 'x' in FFT 'lxg' from 100 to 930 to improve the precision of the numerical FFT. This can be done manually by setting an additional binning named 'cache'.

Simultaneous p.d.f.s are constructed with SIMUL( index, state1=pdf1, state2=pdf2,...)

In [8]:
w.factory("SIMUL::smodel( c[A=0,B=1], A=Gaussian::gs(x,m,s[1]), B=Landau::ls(x,0,1) )")
Out[8]:
<cppyy.gbl.RooSimultaneous object at 0x5a44c40>
[#0] WARNING:InputArguments -- The parameter 's' with range [-1e+30, 1e+30] of the RooGaussian 'gs' exceeds the safe range of (0, inf). Advise to limit its range.

Operator function examples

Function multiplication is done with prod (func1, func2,...)

In [9]:
w.factory("prod::uv(u[10],v[10])")
Out[9]:
<cppyy.gbl.RooProduct object at 0x93a1600>

Function addition is done with sum(func1,func2)

In [10]:
w.factory("sum::uv2(u,v)")
Out[10]:
<cppyy.gbl.RooAddition object at 0x93b89a0>

Lagrangian morphing function for the example shown in rf711_lagrangianmorph

In [11]:
infilename = ROOT.gROOT.GetTutorialDir().Data() + "/roofit/input_histos_rf_lagrangianmorph.root"
w.factory(
    "lagrangianmorph::morph($observableName('pTV'),$fileName('"
    + infilename
    + "'),$couplings({cHq3[0,1],SM[1]}),$NewPhysics(cHq3=1,SM=0),$folders({'SM_NPsq0','cHq3_NPsq1','cHq3_NPsq2'}))"
)
Out[11]:
<cppyy.gbl.RooLagrangianMorphFunc object at 0x96a57d0>
[#0] PROGRESS:InputArguments -- initializing physics inputs from file /home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/src/master.build/tutorials/roofit/input_histos_rf_lagrangianmorph.root with object name(s) 'pTV'
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing dataset dh_SM_NPsq0_morph
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing dataset dh_cHq3_NPsq1_morph
[#1] INFO:ObjectHandling -- RooWorkspace::import(w) importing dataset dh_cHq3_NPsq2_morph

Taylor expansion is done with taylorexpand(func,{var1,var2,...},val,order)

In [12]:
w.factory("taylorexpand::te(expr::poly('x^4+5*x^3+2*x^2+x+1',x),{x},0.0,2)")
Out[12]:
<cppyy.gbl.RooPolyFunc object at 0x982b930>

Interpreted and compiled expression based pdfs

Create a ROOT.RooGenericPdf interpreted p.d.f. You can use single quotes to pass the expression string argument

In [13]:
w.factory("EXPR::G('x*x+1',x)")
Out[13]:
<cppyy.gbl.RooGenericPdf object at 0x98475c0>

Create a custom compiled p.d.f similar to the above interpreted p.d.f. The code required to make self p.d.f. is automatically embedded in the workspace

In [14]:
w.factory("CEXPR::GC('x*x+a',{x,a[1]})")
Out[14]:
<cppyy.gbl.RooCFAuto000Pdf object at 0x9d7b3c0>
[#1] INFO:ObjectHandling -- RooWorkspace::autoImportClass(w) importing code of class RooCFAuto000Pdf from /home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/master_TMP/notebooks/RooCFAuto000Pdf.cxx and RooCFAuto000Pdf.h
Info in <TUnixSystem::ACLiC>: creating shared library /home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/master_TMP/notebooks/RooCFAuto000Pdf_cxx.so
In file included from input_line_9:6:
In file included from ./RooCFAuto000Pdf.cxx:11:
./RooCFAuto000Pdf.h:23:20: warning: 'clone' overrides a member function but is not marked 'override' [-Winconsistent-missing-override]
  virtual TObject* clone(const char* newname) const { return new RooCFAuto000Pdf(*this,newname); }
                   ^
/home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/src/master.build/include/RooAbsArg.h:82:20: note: overridden virtual function is here
  virtual TObject* clone(const char* newname=nullptr) const = 0 ;
                   ^
In file included from input_line_9:6:
In file included from ./RooCFAuto000Pdf.cxx:11:
./RooCFAuto000Pdf.h:31:10: warning: 'evaluate' overrides a member function but is not marked 'override' [-Winconsistent-missing-override]
  double evaluate() const ;
         ^
/home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/src/master.build/include/RooAbsReal.h:436:18: note: overridden virtual function is here
  virtual double evaluate() const = 0;
                 ^

Compiled and interpreted functions (rather than p.d.f.s) can be made with the lower case 'expr' and 'cexpr' types

Print workspace contents

In [15]:
w.Print()
RooWorkspace(w) w contents

variables
---------
(Nbkg,Nsig,SM,a,binWidth_pTV,c,cHq3,f,m,mg,nNP0,nNP1,nNP2,nNP3,nNP4,pTV,s,te_c0,te_c1,te_c2,te_x^0,te_x^1,te_x^2,u,v,x,y,z)

p.d.f.s
-------
RooGenericPdf::G[ actualVars=(x) formula="x[0]*x[0]+1" ] = 1
RooCFAuto000Pdf::GC[ x=x a=a ] = 1
RooChebychev::ch[ x=x coefList=(0.1,0.2,-0.3) ] = 0.8
RooAddPdf::extsummodel[ Nsig * gx + Nbkg * ch ] = 0.9/1
RooGaussian::g[ x=x mean=mg sigma=1 ] = 1
RooGaussian::gs[ x=x mean=m sigma=s ] = 1
RooGaussian::gx[ x=x mean=m sigma=1 ] = 1
RooProdPdf::gxycond[ gx * gy|x ] = 1
RooProdPdf::gxz[ gx * gz ] = 1
RooGaussian::gy[ x=y mean=x sigma=1 ] = 1
RooGaussian::gz[ x=z mean=0 sigma=1 ] = 1
RooLandau::lc[ x=x mean=0 sigma=1 ] = 0.178854
RooLandau::ls[ x=x mean=0 sigma=1 ] = 0.178854
RooFFTConvPdf::lxg[ g(x) (*) lc(x) ] = 0.171919
RooSimultaneous::smodel[ indexCat=c A=gs B=ls ] = 1
RooAddPdf::summodel[ f * gx + [%] * ch ] = 0.9/1

functions
--------
RooLagrangianMorphFunc::morph[ physics=(phys_SM_NPsq0_morph,phys_cHq3_NPsq1_morph,phys_cHq3_NPsq2_morph) operators=(cHq3,SM) observables=(pTV) binWidths=(binWidth_pTV) flags=(nNP0,nNP1,nNP2,nNP3,nNP4) binWidth_pTV * SM_NPsq0_morph + binWidth_pTV * cHq3_NPsq1_morph + binWidth_pTV * cHq3_NPsq2_morph ] = 35.9835
RooHistFunc::phys_SM_NPsq0_morph[ depList=(pTV) ] = 65.1309
RooHistFunc::phys_cHq3_NPsq1_morph[ depList=(pTV) ] = 544.966
RooHistFunc::phys_cHq3_NPsq2_morph[ depList=(pTV) ] = 1528.22
RooFormulaVar::poly[ actualVars=(x) formula="x^4+5*x^3+2*x^2+x+1" ] = 1
RooPolyFunc::te[ vars=(x) te_t0=(te_x^0,te_c0) te_t1=(te_x^1,te_c1) te_t2=(te_x^2,te_c2) ] = 1
RooProduct::uv[ u * v ] = 100
RooAddition::uv2[ u + v ] = 20

embedded datasets (in pdfs and functions)
-----------------------------------------
RooDataHist::dh_SM_NPsq0_morph(pTV)
RooDataHist::dh_cHq3_NPsq1_morph(pTV)
RooDataHist::dh_cHq3_NPsq2_morph(pTV)

embedded class code
-------------------
RooCFAuto000Pdf

embedded precalculated expensive components
-------------------------------------------
uid = 0 key=lxg_g_CONV_lc_CACHEHIST_Obs[x]_BufFrac0.1_BufStrat0 value=RooDataHist::lxg_g_CONV_lc_CACHEHIST_Obs[x]_BufFrac0.1_BufStrat0 parameters=( mg=0 )

Make workspace visible on command line

In [16]:
ROOT.gDirectory.Add(w)

Draw all canvases

In [17]:
from ROOT import gROOT 
gROOT.GetListOfCanvases().Draw()