Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf204_extrangefit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #204
5 ///
6 /// Extended maximum likelihood fit with alternate range definition
7 /// for observed number of events.
8 /// If multiple ranges are used, or only a part of the data is fitted,
9 /// it is advisable to use a RooAddPdf to extend the model. See tutorial
10 /// 204a.
11 ///
12 /// \macro_output
13 /// \macro_code
14 /// \author 07/2008 - Wouter Verkerke
15 
16 #include "RooRealVar.h"
17 #include "RooDataSet.h"
18 #include "RooGaussian.h"
19 #include "RooChebychev.h"
20 #include "RooAddPdf.h"
21 #include "RooExtendPdf.h"
22 #include "RooFitResult.h"
23 #include "TCanvas.h"
24 #include "TAxis.h"
25 #include "RooPlot.h"
26 using namespace RooFit;
27 
28 void rf204_extrangefit()
29 {
30 
31  // S e t u p c o m p o n e n t p d f s
32  // ---------------------------------------
33 
34  // Declare observable x
35  RooRealVar x("x", "x", 0, 10);
36 
37  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
38  RooRealVar mean("mean", "mean of gaussians", 5);
39  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
40  RooRealVar sigma2("sigma2", "width of gaussians", 1);
41 
42  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
43  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
44 
45  // Build Chebychev polynomial p.d.f.
46  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
47  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
48  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
49 
50  // Sum the signal components into a composite signal p.d.f.
51  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
52  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
53 
54  // C o n s t r u c t e x t e n d e d c o m p s wi t h r a n g e s p e c
55  // ------------------------------------------------------------------------------
56 
57  // Define signal range in which events counts are to be defined
58  x.setRange("signalRange", 4, 6);
59 
60  // Associated nsig/nbkg as expected number of events with sig/bkg _in_the_range_ "signalRange"
61  RooRealVar nsig("nsig", "number of signal events in signalRange", 500, 0., 10000) ;
62  RooRealVar nbkg("nbkg", "number of background events in signalRange", 500, 0, 10000) ;
63 
64  // Use AddPdf to extend the model:
65  RooAddPdf model("model","(g1+g2)+a", RooArgList(bkg,sig), RooArgList(nbkg,nsig)) ;
66 
67  // Clone these models here because the interpretation of normalisation coefficients changes
68  // when different ranges are used:
69  RooAddPdf model2(model);
70  RooAddPdf model3(model);
71 
72 
73  // S a m p l e d a t a , f i t m o d e l
74  // -------------------------------------------
75 
76  // Generate 1000 events from model so that nsig,nbkg come out to numbers <<500 in fit
77  RooDataSet *data = model.generate(x, 1000);
78 
79 
80  auto canv = new TCanvas("Canvas", "Canvas", 1500, 600);
81  canv->Divide(3,1);
82 
83  // Fit full range
84  // -------------------------------------------
85 
86  canv->cd(1);
87 
88  // Perform unbinned ML fit to data, full range
89  RooFitResult* r = model.fitTo(*data,Save()) ;
90  r->Print() ;
91 
92 }