rf209_anaconv

Addition and convolution: decay function pdfs with optional B physics effects (mixing and CP violation)

that can be analytically convolved with e.g. Gaussian resolution functions

 pdf1 = decay(t,tau) (x) delta(t)
 pdf2 = decay(t,tau) (x) gauss(t,m,s)
 pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))

Author: Wouter Verkerke
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Sunday, November 27, 2022 at 11:06 AM.

In [1]:
%%cpp -d
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussModel.h"
#include "RooAddModel.h"
#include "RooTruthModel.h"
#include "RooDecay.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit;

B-physics pdf with truth resolution

Variables of decay pdf

In [2]:
RooRealVar dt("dt", "dt", -10, 10);
RooRealVar tau("tau", "tau", 1.548);

Build a truth resolution model (delta function)

In [3]:
RooTruthModel tm1("tm", "truth model", dt);

Construct decay(t) (x) delta(t)

In [4]:
RooDecay decay_tm("decay_tm", "decay", dt, tau, tm1, RooDecay::DoubleSided);

Plot pdf (dashed)

In [5]:
RooPlot *frame = dt.frame(Title("Bdecay (x) resolution"));
decay_tm.plotOn(frame, LineStyle(kDashed));

B-physics pdf with Gaussian resolution

Build a gaussian resolution model

In [6]:
RooRealVar bias1("bias1", "bias1", 0);
RooRealVar sigma1("sigma1", "sigma1", 1);
RooGaussModel gm1("gm1", "gauss model 1", dt, bias1, sigma1);

Construct decay(t) (x) gauss1(t)

In [7]:
RooDecay decay_gm1("decay_gm1", "decay", dt, tau, gm1, RooDecay::DoubleSided);

Plot pdf

In [8]:
decay_gm1.plotOn(frame);

B-physics pdf with double Gaussian resolution

Build another gaussian resolution model

In [9]:
RooRealVar bias2("bias2", "bias2", 0);
RooRealVar sigma2("sigma2", "sigma2", 5);
RooGaussModel gm2("gm2", "gauss model 2", dt, bias2, sigma2);

Build a composite resolution model fgm1+(1-f)gm2

In [10]:
RooRealVar gm1frac("gm1frac", "fraction of gm1", 0.5);
RooAddModel gmsum("gmsum", "sum of gm1 and gm2", RooArgList(gm1, gm2), gm1frac);

Construct decay(t) (x) (fgm1 + (1-f)gm2)

In [11]:
RooDecay decay_gmsum("decay_gmsum", "decay", dt, tau, gmsum, RooDecay::DoubleSided);

Plot pdf (red)

In [12]:
decay_gmsum.plotOn(frame, LineColor(kRed));

Draw all frames on canvas

In [13]:
new TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();

Draw all canvases

In [14]:
%jsroot on
gROOT->GetListOfCanvases()->Draw()