Logo ROOT   6.10/09
Reference Guide
rf209_anaconv.cxx
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////
2 //
3 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #209
4 //
5 // Decay function p.d.fs with optional B physics
6 // effects (mixing and CP violation) that can be
7 // analytically convolved with e.g. Gaussian resolution
8 // functions
9 //
10 // pdf1 = decay(t,tau) (x) delta(t)
11 // pdf2 = decay(t,tau) (x) gauss(t,m,s)
12 // pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
13 //
14 // 07/2008 - Wouter Verkerke
15 //
16 /////////////////////////////////////////////////////////////////////////
17 
18 #ifndef __CINT__
19 #include "RooGlobalFunc.h"
20 #endif
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 "TH1.h"
30 using namespace RooFit ;
31 
32 
33 
34 class TestBasic209 : public RooFitTestUnit
35 {
36 public:
37  TestBasic209(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Analytical convolution operator",refFile,writeRef,verbose) {} ;
38  Bool_t testCode() {
39 
40  // 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
41  // ---------------------------------------------------------------------
42 
43  // Variables of decay p.d.f.
44  RooRealVar dt("dt","dt",-10,10) ;
45  RooRealVar tau("tau","tau",1.548) ;
46 
47  // Build a truth resolution model (delta function)
48  RooTruthModel tm("tm","truth model",dt) ;
49 
50  // Construct decay(t) (x) delta(t)
51  RooDecay decay_tm("decay_tm","decay",dt,tau,tm,RooDecay::DoubleSided) ;
52 
53  // Plot p.d.f. (dashed)
54  RooPlot* frame = dt.frame(Title("Bdecay (x) resolution")) ;
55  decay_tm.plotOn(frame,LineStyle(kDashed)) ;
56 
57 
58  // 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
59  // ----------------------------------------------------------------------------
60 
61  // Build a gaussian resolution model
62  RooRealVar bias1("bias1","bias1",0) ;
63  RooRealVar sigma1("sigma1","sigma1",1) ;
64  RooGaussModel gm1("gm1","gauss model 1",dt,bias1,sigma1) ;
65 
66  // Construct decay(t) (x) gauss1(t)
67  RooDecay decay_gm1("decay_gm1","decay",dt,tau,gm1,RooDecay::DoubleSided) ;
68 
69  // Plot p.d.f.
70  decay_gm1.plotOn(frame) ;
71 
72 
73  // 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
74  // ------------------------------------------------------------------------------------------
75 
76  // Build another gaussian resolution model
77  RooRealVar bias2("bias2","bias2",0) ;
78  RooRealVar sigma2("sigma2","sigma2",5) ;
79  RooGaussModel gm2("gm2","gauss model 2",dt,bias2,sigma2) ;
80 
81  // Build a composite resolution model f*gm1+(1-f)*gm2
82  RooRealVar gm1frac("gm1frac","fraction of gm1",0.5) ;
83  RooAddModel gmsum("gmsum","sum of gm1 and gm2",RooArgList(gm1,gm2),gm1frac) ;
84 
85  // Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
86  RooDecay decay_gmsum("decay_gmsum","decay",dt,tau,gmsum,RooDecay::DoubleSided) ;
87 
88  // Plot p.d.f. (red)
89  decay_gmsum.plotOn(frame,LineColor(kRed)) ;
90 
91  regPlot(frame,"rf209_plot1") ;
92 
93  return kTRUE ;
94 
95  }
96 } ;
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
RooCmdArg LineColor(Color_t color)
Definition: Rtypes.h:56
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
RooCmdArg LineStyle(Style_t style)
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
bool verbose
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
const Bool_t kTRUE
Definition: RtypesCore.h:91
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26