Logo ROOT   6.16/01
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' RooFit tutorial macro #209
5///
6/// Decay function p.d.fs with optional B physics
7/// effects (mixing and CP violation) that can be
8/// analytically convolved with e.g. Gaussian resolution
9/// functions
10///
11/// pdf1 = decay(t,tau) (x) delta(t)
12/// pdf2 = decay(t,tau) (x) gauss(t,m,s)
13/// pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
14///
15/// \macro_image
16/// \macro_output
17/// \macro_code
18/// \author 07/2008 - Wouter Verkerke
19
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
33
34
35void rf209_anaconv()
36{
37 // 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
38 // ---------------------------------------------------------------------
39
40 // Variables of decay p.d.f.
41 RooRealVar dt("dt","dt",-10,10) ;
42 RooRealVar tau("tau","tau",1.548) ;
43
44 // Build a truth resolution model (delta function)
45 RooTruthModel tm1("tm","truth model",dt) ;
46
47 // Construct decay(t) (x) delta(t)
48 RooDecay decay_tm("decay_tm","decay",dt,tau,tm1,RooDecay::DoubleSided) ;
49
50 // Plot p.d.f. (dashed)
51 RooPlot* frame = dt.frame(Title("Bdecay (x) resolution")) ;
52 decay_tm.plotOn(frame,LineStyle(kDashed)) ;
53
54
55 // 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
56 // ----------------------------------------------------------------------------
57
58 // Build a gaussian resolution model
59 RooRealVar bias1("bias1","bias1",0) ;
60 RooRealVar sigma1("sigma1","sigma1",1) ;
61 RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
62
63 // Construct decay(t) (x) gauss1(t)
64 RooDecay decay_gm1("decay_gm1","decay",dt,tau,gm1,RooDecay::DoubleSided) ;
65
66 // Plot p.d.f.
67 decay_gm1.plotOn(frame) ;
68
69
70 // 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
71 // ------------------------------------------------------------------------------------------
72
73 // Build another gaussian resolution model
74 RooRealVar bias2("bias2","bias2",0) ;
75 RooRealVar sigma2("sigma2","sigma2",5) ;
76 RooGaussModel gm2("gm2","gauss model 2",dt,bias2,sigma2) ;
77
78 // Build a composite resolution model f*gm1+(1-f)*gm2
79 RooRealVar gm1frac("gm1frac","fraction of gm1",0.5) ;
80 RooAddModel gmsum("gmsum","sum of gm1 and gm2",RooArgList(gm1,gm2),gm1frac) ;
81
82 // Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
83 RooDecay decay_gmsum("decay_gmsum","decay",dt,tau,gmsum,RooDecay::DoubleSided) ;
84
85 // Plot p.d.f. (red)
86 decay_gmsum.plotOn(frame,LineColor(kRed)) ;
87
88
89 // Draw all frames on canvas
90 new TCanvas("rf209_anaconv","rf209_anaconv",600, 600);
91 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.6) ; frame->Draw() ;
92
93}
@ kRed
Definition: Rtypes.h:63
@ kDashed
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
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:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
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 Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
const char * Title
Definition: TXMLSetup.cxx:67