Logo ROOT  
Reference Guide
rf209_anaconv.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Addition and convolution: decay function p.d.fs with optional B physics effects (mixing and CP violation)
5///
6/// that can be analytically convolved with e.g. Gaussian resolution functions
7///
8/// pdf1 = decay(t,tau) (x) delta(t)
9/// pdf2 = decay(t,tau) (x) gauss(t,m,s)
10/// pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
11///
12/// \macro_image
13/// \macro_output
14/// \macro_code
15/// \author 07/2008 - Wouter Verkerke
16
17#include "RooRealVar.h"
18#include "RooDataSet.h"
19#include "RooGaussModel.h"
20#include "RooAddModel.h"
21#include "RooTruthModel.h"
22#include "RooDecay.h"
23#include "RooPlot.h"
24#include "TCanvas.h"
25#include "TAxis.h"
26#include "TH1.h"
27using namespace RooFit;
28
29void rf209_anaconv()
30{
31 // B - p h y s i c s p d f w i t h t r u t h r e s o l u t i o n
32 // ---------------------------------------------------------------------
33
34 // Variables of decay p.d.f.
35 RooRealVar dt("dt", "dt", -10, 10);
36 RooRealVar tau("tau", "tau", 1.548);
37
38 // Build a truth resolution model (delta function)
39 RooTruthModel tm1("tm", "truth model", dt);
40
41 // Construct decay(t) (x) delta(t)
42 RooDecay decay_tm("decay_tm", "decay", dt, tau, tm1, RooDecay::DoubleSided);
43
44 // Plot p.d.f. (dashed)
45 RooPlot *frame = dt.frame(Title("Bdecay (x) resolution"));
46 decay_tm.plotOn(frame, LineStyle(kDashed));
47
48 // B - p h y s i c s p d f w i t h G a u s s i a n r e s o l u t i o n
49 // ----------------------------------------------------------------------------
50
51 // Build a gaussian resolution model
52 RooRealVar bias1("bias1", "bias1", 0);
53 RooRealVar sigma1("sigma1", "sigma1", 1);
54 RooGaussModel gm1("gm1", "gauss model 1", dt, bias1, sigma1);
55
56 // Construct decay(t) (x) gauss1(t)
57 RooDecay decay_gm1("decay_gm1", "decay", dt, tau, gm1, RooDecay::DoubleSided);
58
59 // Plot p.d.f.
60 decay_gm1.plotOn(frame);
61
62 // B - p h y s i c s p d f w i t h d o u b l e G a u s s i a n r e s o l u t i o n
63 // ------------------------------------------------------------------------------------------
64
65 // Build another gaussian resolution model
66 RooRealVar bias2("bias2", "bias2", 0);
67 RooRealVar sigma2("sigma2", "sigma2", 5);
68 RooGaussModel gm2("gm2", "gauss model 2", dt, bias2, sigma2);
69
70 // Build a composite resolution model f*gm1+(1-f)*gm2
71 RooRealVar gm1frac("gm1frac", "fraction of gm1", 0.5);
72 RooAddModel gmsum("gmsum", "sum of gm1 and gm2", RooArgList(gm1, gm2), gm1frac);
73
74 // Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
75 RooDecay decay_gmsum("decay_gmsum", "decay", dt, tau, gmsum, RooDecay::DoubleSided);
76
77 // Plot p.d.f. (red)
78 decay_gmsum.plotOn(frame, LineColor(kRed));
79
80 // Draw all frames on canvas
81 new TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600);
82 gPad->SetLeftMargin(0.15);
83 frame->GetYaxis()->SetTitleOffset(1.6);
84 frame->Draw();
85}
@ kRed
Definition: Rtypes.h:64
@ kDashed
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
@ DoubleSided
Definition: RooDecay.h:25
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition: RooPlot.cxx:277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
const char * Title
Definition: TXMLSetup.cxx:67