Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
rf204a_extrangefit_RooAddPdf.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 /// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #204a
5 ///
6 /// Extended maximum likelihood fit in multiple ranges.
7 /// When an extended pdf and multiple ranges are used, the
8 /// RooExtendPdf cannot correctly interpret the coefficients
9 /// used for extension.
10 /// This can be solved by using a RooAddPdf for extending the model.
11 ///
12 /// \macro_output
13 /// \macro_code
14 /// \author 12/2018 - Stephan Hageboeck
15 
16 
17 #include "RooRealVar.h"
18 #include "RooDataSet.h"
19 #include "RooGaussian.h"
20 #include "RooChebychev.h"
21 #include "RooAddPdf.h"
22 #include "RooExtendPdf.h"
23 #include "RooFitResult.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "RooPlot.h"
27 using namespace RooFit ;
28 
29 
30 void rf204a_extrangefit_RooAddPdf()
31 {
32 
33 
34  // S e t u p c o m p o n e n t p d f s
35  // ---------------------------------------
36 
37  // Declare observable x
38  RooRealVar x("x","x",0,11) ;
39 
40  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
41  RooRealVar mean("mean","mean of gaussians",5) ;
42  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
43  RooRealVar sigma2("sigma2","width of gaussians",1) ;
44 
45  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
46  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
47 
48  // Build Chebychev polynomial p.d.f.
49  RooRealVar a0("a0","a0",0.5,0.,1.) ;
50  RooRealVar a1("a1","a1",0.2,0.,1.) ;
51  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
52 
53  // Sum the signal components into a composite signal p.d.f.
54  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
55  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
56 
57 
58  // E x t e n d t h e p d f s
59  // -----------------------------
60 
61 
62  // Define signal range in which events counts are to be defined
63  x.setRange("signalRange",4,6) ;
64 
65  // Associated nsig/nbkg as expected number of events with sig/bkg _in_the_range_ "signalRange"
66  RooRealVar nsig("nsig","number of signal events in signalRange",500,0.,10000) ;
67  RooRealVar nbkg("nbkg","number of background events in signalRange",500,0,10000) ;
68 
69  // Use AddPdf to extend the model. Giving as many coefficients as pdfs switches
70  // on extension.
71  RooAddPdf model("model","(g1+g2)+a", RooArgList(bkg,sig), RooArgList(nbkg,nsig)) ;
72 
73 
74  // S a m p l e d a t a , f i t m o d e l
75  // -------------------------------------------
76 
77  // Generate 1000 events from model so that nsig,nbkg come out to numbers <<500 in fit
78  RooDataSet *data = model.generate(x,1000) ;
79 
80 
81 
82  auto canv = new TCanvas("Canvas", "Canvas", 1500, 600);
83  canv->Divide(3,1);
84 
85  // Fit full range
86  // -------------------------------------------
87 
88  canv->cd(1);
89 
90  // Perform unbinned ML fit to data, full range
91 
92  // IMPORTANT:
93  // The model needs to be copied when fitting with different ranges because
94  // the interpretation of the coefficients is tied to the fit range
95  // that's used in the first fit
96  RooAddPdf model1(model);
97  RooFitResult* r = model1.fitTo(*data,Save()) ;
98  r->Print() ;
99 
100  RooPlot * frame = x.frame(Title("Full range fitted"));
101  data->plotOn(frame);
102  model1.plotOn(frame, VisualizeError(*r));
103  model1.plotOn(frame);
104  model1.paramOn(frame);
105  frame->Draw();
106 
107 
108  // Fit in two regions
109  // -------------------------------------------
110 
111  canv->cd(2);
112  x.setRange("left", 0., 4.);
113  x.setRange("right", 6., 10.);
114 
115  RooAddPdf model2(model);
116  RooFitResult* r2 = model2.fitTo(*data,
117  Range("left,right"),
118  Save()) ;
119  r2->Print();
120 
121 
122  RooPlot * frame2 = x.frame(Title("Fit in left/right sideband"));
123  data->plotOn(frame2);
124  model2.plotOn(frame2, VisualizeError(*r2));
125  model2.plotOn(frame2);
126  model2.paramOn(frame2);
127  frame2->Draw();
128 
129 
130  // Fit in one region
131  // -------------------------------------------
132  // Note how restricting the region to only the left tail increases
133  // the fit uncertainty
134 
135  canv->cd(3);
136  x.setRange("leftToMiddle", 0., 5.);
137 
138  RooAddPdf model3(model);
139  RooFitResult* r3 = model3.fitTo(*data,
140  Range("leftToMiddle"),
141  Save()) ;
142  r3->Print();
143 
144 
145  RooPlot * frame3 = x.frame(Title("Fit from left to middle"));
146  data->plotOn(frame3);
147  model3.plotOn(frame3, VisualizeError(*r3));
148  model3.plotOn(frame3);
149  model3.paramOn(frame3);
150  frame3->Draw();
151 
152  canv->Draw();
153 }