Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf202_extendedmlfit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// Addition and convolution: setting up an extended maximum likelihood fit
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 /// \author 07/2008 - Wouter Verkerke
10 
11 #include "RooRealVar.h"
12 #include "RooDataSet.h"
13 #include "RooGaussian.h"
14 #include "RooChebychev.h"
15 #include "RooAddPdf.h"
16 #include "RooExtendPdf.h"
17 #include "TCanvas.h"
18 #include "TAxis.h"
19 #include "RooPlot.h"
20 using namespace RooFit;
21 
22 void rf202_extendedmlfit()
23 {
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  // Sum the signal components into a composite signal p.d.f.
45  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
46  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
47 
48  //----------------
49  // M E T H O D 1
50  //================
51 
52  // C o n s t r u c t e x t e n d e d c o m p o s i t e m o d e l
53  // -------------------------------------------------------------------
54 
55  // Sum the composite signal and background into an extended pdf nsig*sig+nbkg*bkg
56  RooRealVar nsig("nsig", "number of signal events", 500, 0., 10000);
57  RooRealVar nbkg("nbkg", "number of background events", 500, 0, 10000);
58  RooAddPdf model("model", "(g1+g2)+a", RooArgList(bkg, sig), RooArgList(nbkg, nsig));
59 
60  // S a m p l e , f i t a n d p l o t e x t e n d e d m o d e l
61  // ---------------------------------------------------------------------
62 
63  // Generate a data sample of expected number events in x from model
64  // = model.expectedEvents() = nsig+nbkg
65  RooDataSet *data = model.generate(x);
66 
67  // Fit model to data, extended ML term automatically included
68  model.fitTo(*data);
69 
70  // Plot data and PDF overlaid, use expected number of events for p.d.f projection normalization
71  // rather than observed number of events (==data->numEntries())
72  RooPlot *xframe = x.frame(Title("extended ML fit example"));
73  data->plotOn(xframe);
74  model.plotOn(xframe, Normalization(1.0, RooAbsReal::RelativeExpected));
75 
76  // Overlay the background component of model with a dashed line
77  model.plotOn(xframe, Components(bkg), LineStyle(kDashed), Normalization(1.0, RooAbsReal::RelativeExpected));
78 
79  // Overlay the background+sig2 components of model with a dotted line
80  model.plotOn(xframe, Components(RooArgSet(bkg, sig2)), LineStyle(kDotted),
81  Normalization(1.0, RooAbsReal::RelativeExpected));
82 
83  // Print structure of composite p.d.f.
84  model.Print("t");
85 
86  //----------------
87  // M E T H O D 2
88  //================
89 
90  // C o n s t r u c t e x t e n d e d c o m p o n e n t s f i r s t
91  // ---------------------------------------------------------------------
92 
93  // Associated nsig/nbkg as expected number of events with sig/bkg
94  RooExtendPdf esig("esig", "extended signal p.d.f", sig, nsig);
95  RooExtendPdf ebkg("ebkg", "extended background p.d.f", bkg, nbkg);
96 
97  // S u m e x t e n d e d c o m p o n e n t s w i t h o u t c o e f s
98  // -------------------------------------------------------------------------
99 
100  // Construct sum of two extended p.d.f. (no coefficients required)
101  RooAddPdf model2("model2", "(g1+g2)+a", RooArgList(ebkg, esig));
102 
103  // Draw the frame on the canvas
104  new TCanvas("rf202_composite", "rf202_composite", 600, 600);
105  gPad->SetLeftMargin(0.15);
106  xframe->GetYaxis()->SetTitleOffset(1.4);
107  xframe->Draw();
108 }