Logo ROOT   6.16/01
Reference Guide
rf501_simultaneouspdf.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #501
5///
6/// Using simultaneous p.d.f.s to describe simultaneous fits to multiple
7/// datasets
8///
9/// \macro_image
10/// \macro_output
11/// \macro_code
12/// \author 07/2008 - Wouter Verkerke
13
14
15#include "RooRealVar.h"
16#include "RooDataSet.h"
17#include "RooGaussian.h"
18#include "RooConstVar.h"
19#include "RooChebychev.h"
20#include "RooAddPdf.h"
21#include "RooSimultaneous.h"
22#include "RooCategory.h"
23#include "TCanvas.h"
24#include "TAxis.h"
25#include "RooPlot.h"
26using namespace RooFit ;
27
28
30{
31 // C r e a t e m o d e l f o r p h y s i c s s a m p l e
32 // -------------------------------------------------------------
33
34 // Create observables
35 RooRealVar x("x","x",-8,8) ;
36
37 // Construct signal pdf
38 RooRealVar mean("mean","mean",0,-8,8) ;
39 RooRealVar sigma("sigma","sigma",0.3,0.1,10) ;
40 RooGaussian gx("gx","gx",x,mean,sigma) ;
41
42 // Construct background pdf
43 RooRealVar a0("a0","a0",-0.1,-1,1) ;
44 RooRealVar a1("a1","a1",0.004,-1,1) ;
45 RooChebychev px("px","px",x,RooArgSet(a0,a1)) ;
46
47 // Construct composite pdf
48 RooRealVar f("f","f",0.2,0.,1.) ;
49 RooAddPdf model("model","model",RooArgList(gx,px),f) ;
50
51
52
53 // C r e a t e m o d e l f o r c o n t r o l s a m p l e
54 // --------------------------------------------------------------
55
56 // Construct signal pdf.
57 // NOTE that sigma is shared with the signal sample model
58 RooRealVar mean_ctl("mean_ctl","mean_ctl",-3,-8,8) ;
59 RooGaussian gx_ctl("gx_ctl","gx_ctl",x,mean_ctl,sigma) ;
60
61 // Construct the background pdf
62 RooRealVar a0_ctl("a0_ctl","a0_ctl",-0.1,-1,1) ;
63 RooRealVar a1_ctl("a1_ctl","a1_ctl",0.5,-0.1,1) ;
64 RooChebychev px_ctl("px_ctl","px_ctl",x,RooArgSet(a0_ctl,a1_ctl)) ;
65
66 // Construct the composite model
67 RooRealVar f_ctl("f_ctl","f_ctl",0.5,0.,1.) ;
68 RooAddPdf model_ctl("model_ctl","model_ctl",RooArgList(gx_ctl,px_ctl),f_ctl) ;
69
70
71
72 // G e n e r a t e e v e n t s f o r b o t h s a m p l e s
73 // ---------------------------------------------------------------
74
75 // Generate 1000 events in x and y from model
76 RooDataSet *data = model.generate(RooArgSet(x),100) ;
77 RooDataSet *data_ctl = model_ctl.generate(RooArgSet(x),2000) ;
78
79
80
81 // C r e a t e i n d e x c a t e g o r y a n d j o i n s a m p l e s
82 // ---------------------------------------------------------------------------
83
84 // Define category to distinguish physics and control samples events
85 RooCategory sample("sample","sample") ;
86 sample.defineType("physics") ;
87 sample.defineType("control") ;
88
89 // Construct combined dataset in (x,sample)
90 RooDataSet combData("combData","combined data",x,Index(sample),Import("physics",*data),Import("control",*data_ctl)) ;
91
92
93
94 // C o n s t r u c t a s i m u l t a n e o u s p d f i n ( x , s a m p l e )
95 // -----------------------------------------------------------------------------------
96
97 // Construct a simultaneous pdf using category sample as index
98 RooSimultaneous simPdf("simPdf","simultaneous pdf",sample) ;
99
100 // Associate model with the physics state and model_ctl with the control state
101 simPdf.addPdf(model,"physics") ;
102 simPdf.addPdf(model_ctl,"control") ;
103
104
105
106 // P e r f o r m a s i m u l t a n e o u s f i t
107 // ---------------------------------------------------
108
109 // Perform simultaneous fit of model to data and model_ctl to data_ctl
110 simPdf.fitTo(combData) ;
111
112
113
114 // P l o t m o d e l s l i c e s o n d a t a s l i c e s
115 // ----------------------------------------------------------------
116
117 // Make a frame for the physics sample
118 RooPlot* frame1 = x.frame(Bins(30),Title("Physics sample")) ;
119
120 // Plot all data tagged as physics sample
121 combData.plotOn(frame1,Cut("sample==sample::physics")) ;
122
123 // Plot "physics" slice of simultaneous pdf.
124 // NBL You _must_ project the sample index category with data using ProjWData
125 // as a RooSimultaneous makes no prediction on the shape in the index category
126 // and can thus not be integrated
127 simPdf.plotOn(frame1,Slice(sample,"physics"),ProjWData(sample,combData)) ;
128 simPdf.plotOn(frame1,Slice(sample,"physics"),Components("px"),ProjWData(sample,combData),LineStyle(kDashed)) ;
129
130 // The same plot for the control sample slice
131 RooPlot* frame2 = x.frame(Bins(30),Title("Control sample")) ;
132 combData.plotOn(frame2,Cut("sample==sample::control")) ;
133 simPdf.plotOn(frame2,Slice(sample,"control"),ProjWData(sample,combData)) ;
134 simPdf.plotOn(frame2,Slice(sample,"control"),Components("px_ctl"),ProjWData(sample,combData),LineStyle(kDashed)) ;
135
136
137
138 TCanvas* c = new TCanvas("rf501_simultaneouspdf","rf403_simultaneouspdf",800,400) ;
139 c->Divide(2) ;
140 c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
141 c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
142
143
144}
#define f(i)
Definition: RSha256.hxx:104
#define c(i)
Definition: RSha256.hxx:101
@ kDashed
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
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
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
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
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
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
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Index(RooCategory &icat)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg Slice(const RooArgSet &sliceSet)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg Bins(Int_t nbin)
RooCmdArg LineStyle(Style_t style)
const char * Title
Definition: TXMLSetup.cxx:67