Logo ROOT  
Reference Guide
rf202_extendedmlfit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -js
4 ///
5 /// Addition and convolution: setting up an extended maximum likelihood fit
6 ///
7 /// \macro_image
8 /// \macro_output
9 /// \macro_code
10 ///
11 /// \date 07/2008
12 /// \author Wouter Verkerke
13 
14 #include "RooRealVar.h"
15 #include "RooDataSet.h"
16 #include "RooGaussian.h"
17 #include "RooChebychev.h"
18 #include "RooAddPdf.h"
19 #include "RooExtendPdf.h"
20 #include "TCanvas.h"
21 #include "TAxis.h"
22 #include "RooPlot.h"
23 using namespace RooFit;
24 
26 {
27 
28  // S e t u p c o m p o n e n t p d f s
29  // ---------------------------------------
30 
31  // Declare observable x
32  RooRealVar x("x", "x", 0, 10);
33 
34  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
35  RooRealVar mean("mean", "mean of gaussians", 5);
36  RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
37  RooRealVar sigma2("sigma2", "width of gaussians", 1);
38 
39  RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
40  RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
41 
42  // Build Chebychev polynomial p.d.f.
43  RooRealVar a0("a0", "a0", 0.5, 0., 1.);
44  RooRealVar a1("a1", "a1", 0.2, 0., 1.);
45  RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
46 
47  // Sum the signal components into a composite signal p.d.f.
48  RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
49  RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
50 
51  //----------------
52  // M E T H O D 1
53  //================
54 
55  // 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
56  // -------------------------------------------------------------------
57 
58  // Sum the composite signal and background into an extended pdf nsig*sig+nbkg*bkg
59  RooRealVar nsig("nsig", "number of signal events", 500, 0., 10000);
60  RooRealVar nbkg("nbkg", "number of background events", 500, 0, 10000);
61  RooAddPdf model("model", "(g1+g2)+a", RooArgList(bkg, sig), RooArgList(nbkg, nsig));
62 
63  // 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
64  // ---------------------------------------------------------------------
65 
66  // Generate a data sample of expected number events in x from model
67  // = model.expectedEvents() = nsig+nbkg
68  RooDataSet *data = model.generate(x);
69 
70  // Fit model to data, extended ML term automatically included
71  model.fitTo(*data);
72 
73  // Plot data and PDF overlaid, use expected number of events for p.d.f projection normalization
74  // rather than observed number of events (==data->numEntries())
75  RooPlot *xframe = x.frame(Title("extended ML fit example"));
76  data->plotOn(xframe);
77  model.plotOn(xframe, Normalization(1.0, RooAbsReal::RelativeExpected));
78 
79  // Overlay the background component of model with a dashed line
81 
82  // Overlay the background+sig2 components of model with a dotted line
83  model.plotOn(xframe, Components(RooArgSet(bkg, sig2)), LineStyle(kDotted),
85 
86  // Print structure of composite p.d.f.
87  model.Print("t");
88 
89  //----------------
90  // M E T H O D 2
91  //================
92 
93  // 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
94  // ---------------------------------------------------------------------
95 
96  // Associated nsig/nbkg as expected number of events with sig/bkg
97  RooExtendPdf esig("esig", "extended signal p.d.f", sig, nsig);
98  RooExtendPdf ebkg("ebkg", "extended background p.d.f", bkg, nbkg);
99 
100  // 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
101  // -------------------------------------------------------------------------
102 
103  // Construct sum of two extended p.d.f. (no coefficients required)
104  RooAddPdf model2("model2", "(g1+g2)+a", RooArgList(ebkg, esig));
105 
106  // Draw the frame on the canvas
107  new TCanvas("rf202_composite", "rf202_composite", 600, 600);
108  gPad->SetLeftMargin(0.15);
109  xframe->GetYaxis()->SetTitleOffset(1.4);
110  xframe->Draw();
111 }
RooPlot::Draw
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:691
RooChebychev.h
RooAddPdf
Definition: RooAddPdf.h:32
RooExtendPdf.h
rf202_extendedmlfit
Definition: rf202_extendedmlfit.py:1
RooChebychev
Definition: RooChebychev.h:25
RooArgList
Definition: RooArgList.h:21
RooGaussian.h
RooFit::Normalization
RooCmdArg Normalization(Double_t scaleFactor)
Definition: RooGlobalFunc.cxx:51
RooExtendPdf
Definition: RooExtendPdf.h:22
x
Double_t x[n]
Definition: legend1.C:17
RooGaussian
Definition: RooGaussian.h:25
RooAddPdf.h
TCanvas.h
RooDataSet.h
RooFit
Definition: RooCFunction1Binding.h:29
RooAbsData::plotOn
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Definition: RooAbsData.cxx:547
RooPlot.h
RooPlot::GetYaxis
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1256
RooPlot
Definition: RooPlot.h:44
RooRealVar.h
TCanvas
Definition: TCanvas.h:23
TAxis.h
kDashed
@ kDashed
Definition: TAttLine.h:48
gPad
#define gPad
Definition: TVirtualPad.h:287
RooDataSet
Definition: RooDataSet.h:33
make_cnn_model.model
model
Definition: make_cnn_model.py:6
RooFit::LineStyle
RooCmdArg LineStyle(Style_t style)
Definition: RooGlobalFunc.cxx:57
TAttAxis::SetTitleOffset
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:293
RooRealVar
Definition: RooRealVar.h:35
RooFit::Components
RooCmdArg Components(const RooArgSet &compSet)
Definition: RooGlobalFunc.cxx:74
RooFit::Title
RooCmdArg Title(const char *name)
Definition: RooGlobalFunc.cxx:173
kDotted
@ kDotted
Definition: TAttLine.h:48
RooArgSet
Definition: RooArgSet.h:28
RooAbsReal::RelativeExpected
@ RelativeExpected
Definition: RooAbsReal.h:253