Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 pdfs with optional B physics effects (mixing and CP violation)
5///
6/// that can be analytically convolved with e.g. Gaussian resolution functions
7///
8/// ```
9/// pdf1 = decay(t,tau) (x) delta(t)
10/// pdf2 = decay(t,tau) (x) gauss(t,m,s)
11/// pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
12/// ```
13///
14/// \macro_image
15/// \macro_output
16/// \macro_code
17///
18/// \date July 2008
19/// \author Wouter Verkerke
20
21#include "RooRealVar.h"
22#include "RooDataSet.h"
23#include "RooGaussModel.h"
24#include "RooAddModel.h"
25#include "RooTruthModel.h"
26#include "RooDecay.h"
27#include "RooPlot.h"
28#include "TCanvas.h"
29#include "TAxis.h"
30#include "TH1.h"
31using namespace RooFit;
32
33void rf209_anaconv()
34{
35 // 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
36 // ---------------------------------------------------------------------
37
38 // Variables of decay pdf
39 RooRealVar dt("dt", "dt", -10, 10);
40 RooRealVar tau("tau", "tau", 1.548);
41
42 // Build a truth resolution model (delta function)
43 RooTruthModel tm1("tm", "truth model", dt);
44
45 // Construct decay(t) (x) delta(t)
46 RooDecay decay_tm("decay_tm", "decay", dt, tau, tm1, RooDecay::DoubleSided);
47
48 // Plot pdf (dashed)
49 RooPlot *frame = dt.frame(Title("Bdecay (x) resolution"));
50 decay_tm.plotOn(frame, LineStyle(kDashed));
51
52 // 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
53 // ----------------------------------------------------------------------------
54
55 // Build a gaussian resolution model
56 RooRealVar bias1("bias1", "bias1", 0);
57 RooRealVar sigma1("sigma1", "sigma1", 1);
58 RooGaussModel gm1("gm1", "gauss model 1", dt, bias1, sigma1);
59
60 // Construct decay(t) (x) gauss1(t)
61 RooDecay decay_gm1("decay_gm1", "decay", dt, tau, gm1, RooDecay::DoubleSided);
62
63 // Plot pdf
64 decay_gm1.plotOn(frame);
65
66 // 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
67 // ------------------------------------------------------------------------------------------
68
69 // Build another gaussian resolution model
70 RooRealVar bias2("bias2", "bias2", 0);
71 RooRealVar sigma2("sigma2", "sigma2", 5);
72 RooGaussModel gm2("gm2", "gauss model 2", dt, bias2, sigma2);
73
74 // Build a composite resolution model f*gm1+(1-f)*gm2
75 RooRealVar gm1frac("gm1frac", "fraction of gm1", 0.5);
76 RooAddModel gmsum("gmsum", "sum of gm1 and gm2", RooArgList(gm1, gm2), gm1frac);
77
78 // Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
79 RooDecay decay_gmsum("decay_gmsum", "decay", dt, tau, gmsum, RooDecay::DoubleSided);
80
81 // Plot pdf (red)
82 decay_gmsum.plotOn(frame, LineColor(kRed));
83
84 // Draw all frames on canvas
85 new TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600);
86 gPad->SetLeftMargin(0.15);
87 frame->GetYaxis()->SetTitleOffset(1.6);
88 frame->Draw();
89}
@ kRed
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
#define gPad
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:27
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.
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:1263
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:249
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:293
The Canvas class.
Definition TCanvas.h:23
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
const char * Title
Definition TXMLSetup.cxx:68