Logo ROOT   6.10/09
Reference Guide
rf202_extendedmlfit.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////
2 //
3 // 'ADDITION AND CONVOLUTION' RooFit tutorial macro #202
4 //
5 // Setting up an extended maximum likelihood fit
6 //
7 //
8 //
9 // 07/2008 - Wouter Verkerke
10 //
11 /////////////////////////////////////////////////////////////////////////
12 
13 #ifndef __CINT__
14 #include "RooGlobalFunc.h"
15 #endif
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 "TCanvas.h"
23 #include "RooPlot.h"
24 using namespace RooFit ;
25 
26 
27 class TestBasic202 : public RooFitTestUnit
28 {
29 public:
30  TestBasic202(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Extended ML fits to addition operator p.d.f.s",refFile,writeRef,verbose) {} ;
31  Bool_t testCode() {
32 
33  // S e t u p c o m p o n e n t p d f s
34  // ---------------------------------------
35 
36  // Declare observable x
37  RooRealVar x("x","x",0,10) ;
38 
39  // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
40  RooRealVar mean("mean","mean of gaussians",5) ;
41  RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
42  RooRealVar sigma2("sigma2","width of gaussians",1) ;
43 
44  RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
45  RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
46 
47  // Build Chebychev polynomial p.d.f.
48  RooRealVar a0("a0","a0",0.5,0.,1.) ;
49  RooRealVar a1("a1","a1",-0.2,0.,1.) ;
50  RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
51 
52  // Sum the signal components into a composite signal p.d.f.
53  RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
54  RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
55 
56  /////////////////////
57  // M E T H O D 1 //
58  /////////////////////
59 
60 
61  // 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
62  // -------------------------------------------------------------------
63 
64  // Sum the composite signal and background into an extended pdf nsig*sig+nbkg*bkg
65  RooRealVar nsig("nsig","number of signal events",500,0.,10000) ;
66  RooRealVar nbkg("nbkg","number of background events",500,0,10000) ;
67  RooAddPdf model("model","(g1+g2)+a",RooArgList(bkg,sig),RooArgList(nbkg,nsig)) ;
68 
69 
70 
71  // 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
72  // ---------------------------------------------------------------------
73 
74  // Generate a data sample of expected number events in x from model
75  // = model.expectedEvents() = nsig+nbkg
76  RooDataSet *data = model.generate(x) ;
77 
78  // Fit model to data, extended ML term automatically included
79  model.fitTo(*data) ;
80 
81  // Plot data and PDF overlaid, use expected number of events for p.d.f projection normalization
82  // rather than observed number of events (==data->numEntries())
83  RooPlot* xframe = x.frame(Title("extended ML fit example")) ;
84  data->plotOn(xframe) ;
86 
87  // Overlay the background component of model with a dashed line
89 
90  // Overlay the background+sig2 components of model with a dotted line
92 
93 
94  /////////////////////
95  // M E T H O D 2 //
96  /////////////////////
97 
98  // 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
99  // ---------------------------------------------------------------------
100 
101  // Associated nsig/nbkg as expected number of events with sig/bkg
102  RooExtendPdf esig("esig","extended signal p.d.f",sig,nsig) ;
103  RooExtendPdf ebkg("ebkg","extended background p.d.f",bkg,nbkg) ;
104 
105 
106  // 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
107  // -------------------------------------------------------------------------
108 
109  // Construct sum of two extended p.d.f. (no coefficients required)
110  RooAddPdf model2("model2","(g1+g2)+a",RooArgList(ebkg,esig)) ;
111 
112 
113  regPlot(xframe,"rf202_plot1") ;
114 
115  delete data ;
116  return kTRUE ;
117 
118  }
119 
120 } ;
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:552
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooCmdArg Title(const char *name)
Double_t x[n]
Definition: legend1.C:17
RooCmdArg LineStyle(Style_t style)
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
bool verbose
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg Normalization(Double_t scaleFactor)
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
const Bool_t kTRUE
Definition: RtypesCore.h:91