Logo ROOT   6.16/01
Reference Guide
rf202_extendedmlfit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #202
5///
6/// Setting up an extended maximum likelihood fit
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 "RooExtendPdf.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit ;
23
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
56 // 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
57 // -------------------------------------------------------------------
58
59 // Sum the composite signal and background into an extended pdf nsig*sig+nbkg*bkg
60 RooRealVar nsig("nsig","number of signal events",500,0.,10000) ;
61 RooRealVar nbkg("nbkg","number of background events",500,0,10000) ;
62 RooAddPdf model("model","(g1+g2)+a",RooArgList(bkg,sig),RooArgList(nbkg,nsig)) ;
63
64
65
66 // 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
67 // ---------------------------------------------------------------------
68
69 // Generate a data sample of expected number events in x from model
70 // = model.expectedEvents() = nsig+nbkg
71 RooDataSet *data = model.generate(x) ;
72
73 // Fit model to data, extended ML term automatically included
74 model.fitTo(*data) ;
75
76 // Plot data and PDF overlaid, use expected number of events for p.d.f projection normalization
77 // rather than observed number of events (==data->numEntries())
78 RooPlot* xframe = x.frame(Title("extended ML fit example")) ;
79 data->plotOn(xframe) ;
81
82 // Overlay the background component of model with a dashed line
84
85 // Overlay the background+sig2 components of model with a dotted line
87
88 // Print structure of composite p.d.f.
89 model.Print("t") ;
90
91
92 //----------------
93 // M E T H O D 2
94 //================
95
96 // 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
97 // ---------------------------------------------------------------------
98
99 // Associated nsig/nbkg as expected number of events with sig/bkg
100 RooExtendPdf esig("esig","extended signal p.d.f",sig,nsig) ;
101 RooExtendPdf ebkg("ebkg","extended background p.d.f",bkg,nbkg) ;
102
103
104 // 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
105 // -------------------------------------------------------------------------
106
107 // Construct sum of two extended p.d.f. (no coefficients required)
108 RooAddPdf model2("model2","(g1+g2)+a",RooArgList(ebkg,esig)) ;
109
110
111
112 // Draw the frame on the canvas
113 new TCanvas("rf202_composite","rf202_composite",600,600) ;
114 gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;
115
116
117}
@ kDashed
Definition: TAttLine.h:48
@ kDotted
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
@ RelativeExpected
Definition: RooAbsReal.h:218
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Definition: RooExtendPdf.h:22
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1123
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
Double_t x[n]
Definition: legend1.C:17
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg Normalization(Double_t scaleFactor)
RooCmdArg LineStyle(Style_t style)
const char * Title
Definition: TXMLSetup.cxx:67