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

Detailed Description

View in nbviewer Open in SWAN
Organization and simultaneous fits: illustration use of ROOT.RooCustomizer and ROOT.RooSimWSTool interface in factory workspace tool in a complex standalone B physics example

import ROOT
w = ROOT.RooWorkspace("w")
# Build a complex example pdf
# -----------------------------------------------------------
# Make signal model for CPV: A bmixing decay function in t (convoluted with a triple Gaussian resolution model)
# times a Gaussian function the reconstructed mass
w.factory(
"PROD::sig( BMixDecay::sig_t( dt[-20,20], mixState[mixed=1,unmix=-1], tagFlav[B0=1,B0bar=-1], "
"tau[1.54], dm[0.472], w[0.05], dw[0], "
"AddModel::gm({GaussModel(dt,biasC[-10,10],sigmaC[0.1,3],dterr[0.01,0.2]), "
"GaussModel(dt,0,sigmaT[3,10]), "
"GaussModel(dt,0,20)},{fracC[0,1],fracT[0,1]}), "
"DoubleSided ), "
"Gaussian::sig_m( mes[5.20,5.30], mB0[5.20,5.30], sigmB0[0.01,0.05] ))"
)
# Make background component: A plain decay function in t times an Argus
# function in the reconstructed mass
w.factory("PROD::bkg( Decay::bkg_t( dt, tau, gm, DoubleSided), " "ArgusBG::bkg_m( mes, 5.291, k[-100,-10]))")
# Make composite model from the signal and background component
w.factory("SUM::model( Nsig[5000,0,10000]*sig, NBkg[500,0,10000]*bkg )")
# Example of RooSimWSTool interface
# ------------------------------------------------------------------
# Introduce a flavour tagging category tagCat as observable with 4 states corresponding
# to 4 flavour tagging techniques with different performance that require different
# parameterizations of the fit model
#
# ROOT.RooSimWSTool operation:
# - Make 4 clones of model (for each tagCat) state, will gain an individual
# copy of parameters w, and biasC. The other parameters remain common
# - Make a simultaneous pdf of the 4 clones assigning each to the appropriate
# state of the tagCat index category
# ROOT.RooSimWSTool is interfaced as meta-type SIMCLONE in the factory. The $SplitParam()
# argument maps to the SplitParam() named argument in the
# ROOT.RooSimWSTool constructor
w.factory("SIMCLONE::model_sim( model, $SplitParam({w,dw,biasC},tagCat[Lep,Kao,NT1,NT2]))")
# Example of RooCustomizer interface
# -------------------------------------------------------------------
#
# Class ROOT.RooCustomizer makes clones of existing pdfs with certain prescribed
# modifications (branch of leaf node replacements)
#
# Here we take our model (the original before ROOT.RooSimWSTool modifications)
# and request that the parameter w (the mistag rate) is replaced with
# an expression-based function that calculates w in terms of the Dilution
# parameter D that is defined D = 1-2*w
# Make a clone model_D of original 'model' replacing 'w' with
# 'expr('0.5-D/2',D[0,1])'
w.factory("EDIT::model_D(model, w=expr('0.5-D/2',D[0,1]) )")
# Print workspace contents
w.Print()
# Make workspace visible on command line
ROOT.gDirectory.Add(w)
RooWorkspace(w) w contents
variables
---------
(D,NBkg,Nsig,biasC,biasC_Kao,biasC_Lep,biasC_NT1,biasC_NT2,dm,dt,dterr,dw,dw_Kao,dw_Lep,dw_NT1,dw_NT2,fracC,fracT,k,mB0,mes,mixState,sigmB0,sigmaC,sigmaT,tagCat,tagFlav,tau,w,w_Kao,w_Lep,w_NT1,w_NT2)
p.d.f.s
-------
RooProdPdf::bkg[ bkg_t * bkg_m ] = 0.307193
RooProdPdf::bkg_Kao[ bkg_t_Kao * bkg_m ] = 0.307193
RooProdPdf::bkg_Lep[ bkg_t_Lep * bkg_m ] = 0.307193
RooProdPdf::bkg_NT1[ bkg_t_NT1 * bkg_m ] = 0.307193
RooProdPdf::bkg_NT2[ bkg_t_NT2 * bkg_m ] = 0.307193
RooArgusBG::bkg_m[ m=mes m0=5.291 c=k p=0.5 ] = 0.279062
RooDecay::bkg_t[ t=dt tau=tau ] = 1.10081
RooDecay::bkg_t_Kao[ t=dt tau=tau ] = 1.10081
RooDecay::bkg_t_Lep[ t=dt tau=tau ] = 1.10081
RooDecay::bkg_t_NT1[ t=dt tau=tau ] = 1.10081
RooDecay::bkg_t_NT2[ t=dt tau=tau ] = 1.10081
RooAddPdf::model[ Nsig * sig + NBkg * bkg ] = 1.88229/1
RooAddPdf::model_D[ Nsig * sig_model_D + NBkg * bkg ] = 1.5029/1
RooAddPdf::model_Kao[ Nsig * sig_Kao + NBkg * bkg_Kao ] = 1.88229/1
RooAddPdf::model_Lep[ Nsig * sig_Lep + NBkg * bkg_Lep ] = 1.88229/1
RooAddPdf::model_NT1[ Nsig * sig_NT1 + NBkg * bkg_NT1 ] = 1.88229/1
RooAddPdf::model_NT2[ Nsig * sig_NT2 + NBkg * bkg_NT2 ] = 1.88229/1
RooSimultaneous::model_sim[ indexCat=tagCat Kao=model_Kao Lep=model_Lep NT1=model_NT1 NT2=model_NT2 ] = 0.470573
RooProdPdf::sig[ sig_t * sig_m ] = 2.0398
RooProdPdf::sig_Kao[ sig_t_Kao * sig_m ] = 2.0398
RooProdPdf::sig_Lep[ sig_t_Lep * sig_m ] = 2.0398
RooProdPdf::sig_NT1[ sig_t_NT1 * sig_m ] = 2.0398
RooProdPdf::sig_NT2[ sig_t_NT2 * sig_m ] = 2.0398
RooGaussian::sig_m[ x=mes mean=mB0 sigma=sigmB0 ] = 1
RooProdPdf::sig_model_D[ sig_t_model_D * sig_m ] = 1.62247
RooBMixDecay::sig_t[ mistag=w delMistag=dw mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 2.0398
RooBMixDecay::sig_t_Kao[ mistag=w_Kao delMistag=dw_Kao mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 2.0398
RooBMixDecay::sig_t_Lep[ mistag=w_Lep delMistag=dw_Lep mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 2.0398
RooBMixDecay::sig_t_NT1[ mistag=w_NT1 delMistag=dw_NT1 mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 2.0398
RooBMixDecay::sig_t_NT2[ mistag=w_NT2 delMistag=dw_NT2 mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 2.0398
RooBMixDecay::sig_t_model_D[ mistag=model_D_2 delMistag=dw mixState=mixState tagFlav=tagFlav tau=tau dm=dm t=dt ] = 1.62247
analytical resolution models
----------------------------
RooAddModel::gm[ x=dt (fracC * gm_11 + fracT * gm_12 + [%] * gm_13) ] = 1.25632
RooGaussModel::gm_11[ x=dt mean=biasC sigma=sigmaC msf=dterr ssf=dterr ] = 2.45126
RooGaussModel::gm_11_Kao[ x=dt mean=biasC_Kao sigma=sigmaC msf=dterr ssf=dterr ] = 2.45126
RooGaussModel::gm_11_Lep[ x=dt mean=biasC_Lep sigma=sigmaC msf=dterr ssf=dterr ] = 2.45126
RooGaussModel::gm_11_NT1[ x=dt mean=biasC_NT1 sigma=sigmaC msf=dterr ssf=dterr ] = 2.45126
RooGaussModel::gm_11_NT2[ x=dt mean=biasC_NT2 sigma=sigmaC msf=dterr ssf=dterr ] = 2.45126
RooGaussModel::gm_12[ x=dt mean=0 sigma=sigmaT msf=1 ssf=1 ] = 0.0613757
RooGaussModel::gm_13[ x=dt mean=0 sigma=20 msf=1 ssf=1 ] = 0.0199471
RooAddModel::gm_Kao[ x=dt (fracC * gm_11_Kao + fracT * gm_12 + [%] * gm_13) ] = 1.25632
RooAddModel::gm_Lep[ x=dt (fracC * gm_11_Lep + fracT * gm_12 + [%] * gm_13) ] = 1.25632
RooAddModel::gm_NT1[ x=dt (fracC * gm_11_NT1 + fracT * gm_12 + [%] * gm_13) ] = 1.25632
RooAddModel::gm_NT2[ x=dt (fracC * gm_11_NT2 + fracT * gm_12 + [%] * gm_13) ] = 1.25632
functions
--------
RooFormulaVar::model_D_2[ actualVars=(D) formula="0.5-D/2" ] = 0.25
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf513_wsfactory_tools.py.