Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf201_composite.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Addition and convolution: composite p.d.f with signal and background component
5 ///
6 /// pdf = f_bkg * bkg(x,a0,a1) + (1-fbkg) * (f_sig1 * sig1(x,m,s1 + (1-f_sig1) * sig2(x,m,s2)))
7 ///
8 /// \macro_image
9 /// \macro_output
10 /// \macro_code
11 /// \author 07/2008 - Wouter Verkerke
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "RooChebychev.h"
17 #include "RooAddPdf.h"
18 #include "TCanvas.h"
19 #include "TAxis.h"
20 #include "RooPlot.h"
21 using namespace RooFit;
22 
23 void rf201_composite()
24 {
25  // S e t u p c o m p o n e n t p d f s
26  // ---------------------------------------
27 
28  // Declare observable x
29  RooRealVar x("x", "x", 0, 10);
30 
31  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
32  RooRealVar mean("mean", "mean of gaussians", 5);
33  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
34  RooRealVar sigma2("sigma2", "width of gaussians", 1);
35 
36  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
37  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
38 
39  // Build Chebychev polynomial p.d.f.
40  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
41  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
42  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
43 
44  // ---------------------------------------------
45  // M E T H O D 1 - T w o R o o A d d P d f s
46  // =============================================
47 
48  // A d d s i g n a l c o m p o n e n t s
49  // ------------------------------------------
50 
51  // Sum the signal components into a composite signal p.d.f.
52  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
53  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
54 
55  // A d d s i g n a l a n d b a c k g r o u n d
56  // ------------------------------------------------
57 
58  // Sum the composite signal and background
59  RooRealVar bkgfrac("bkgfrac", "fraction of background", 0.5, 0., 1.);
60  RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), bkgfrac);
61 
62  // S a m p l e , f i t a n d p l o t m o d e l
63  // ---------------------------------------------------
64 
65  // Generate a data sample of 1000 events in x from model
66  RooDataSet *data = model.generate(x, 1000);
67 
68  // Fit model to data
69  model.fitTo(*data);
70 
71  // Plot data and PDF overlaid
72  RooPlot *xframe = x.frame(Title("Example of composite pdf=(sig1+sig2)+bkg"));
73  data->plotOn(xframe);
74  model.plotOn(xframe);
75 
76  // Overlay the background component of model with a dashed line
77  model.plotOn(xframe, Components(bkg), LineStyle(kDashed));
78 
79  // Overlay the background+sig2 components of model with a dotted line
80  model.plotOn(xframe, Components(RooArgSet(bkg, sig2)), LineStyle(kDotted));
81 
82  // Print structure of composite p.d.f.
83  model.Print("t");
84 
85  // ---------------------------------------------------------------------------------------------
86  // M E T H O D 2 - O n e R o o A d d P d f w i t h r e c u r s i v e f r a c t i o n s
87  // =============================================================================================
88 
89  // Construct sum of models on one go using recursive fraction interpretations
90  //
91  // model2 = bkg + (sig1 + sig2)
92  //
93  RooAddPdf model2("model", "g1+g2+a", RooArgList(bkg, sig1, sig2), RooArgList(bkgfrac, sig1frac), kTRUE);
94 
95  // NB: Each coefficient is interpreted as the fraction of the
96  // left-hand component of the i-th recursive sum, i.e.
97  //
98  // sum4 = A + ( B + ( C + D) with fraction fA, fB and fC expands to
99  //
100  // sum4 = fA*A + (1-fA)*(fB*B + (1-fB)*(fC*C + (1-fC)*D))
101 
102  // P l o t r e c u r s i v e a d d i t i o n m o d e l
103  // ---------------------------------------------------------
104  model2.plotOn(xframe, LineColor(kRed), LineStyle(kDashed));
105  model2.plotOn(xframe, Components(RooArgSet(bkg, sig2)), LineColor(kRed), LineStyle(kDashed));
106  model2.Print("t");
107 
108  // Draw the frame on the canvas
109  new TCanvas("rf201_composite", "rf201_composite", 600, 600);
110  gPad->SetLeftMargin(0.15);
111  xframe->GetYaxis()->SetTitleOffset(1.4);
112  xframe->Draw();
113 }