Logo ROOT  
Reference Guide
rf501_simultaneouspdf.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Organisation and simultaneous fits: using simultaneous p.d.f.s to describe simultaneous fits to multiple datasets
5///
6/// \macro_image
7/// \macro_output
8/// \macro_code
9/// \author 07/2008 - Wouter Verkerke
10
11#include "RooRealVar.h"
12#include "RooDataSet.h"
13#include "RooGaussian.h"
14#include "RooConstVar.h"
15#include "RooChebychev.h"
16#include "RooAddPdf.h"
17#include "RooSimultaneous.h"
18#include "RooCategory.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "RooPlot.h"
22using namespace RooFit;
23
25{
26 // 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
27 // -------------------------------------------------------------
28
29 // Create observables
30 RooRealVar x("x", "x", -8, 8);
31
32 // Construct signal pdf
33 RooRealVar mean("mean", "mean", 0, -8, 8);
34 RooRealVar sigma("sigma", "sigma", 0.3, 0.1, 10);
35 RooGaussian gx("gx", "gx", x, mean, sigma);
36
37 // Construct background pdf
38 RooRealVar a0("a0", "a0", -0.1, -1, 1);
39 RooRealVar a1("a1", "a1", 0.004, -1, 1);
40 RooChebychev px("px", "px", x, RooArgSet(a0, a1));
41
42 // Construct composite pdf
43 RooRealVar f("f", "f", 0.2, 0., 1.);
44 RooAddPdf model("model", "model", RooArgList(gx, px), f);
45
46 // 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
47 // --------------------------------------------------------------
48
49 // Construct signal pdf.
50 // NOTE that sigma is shared with the signal sample model
51 RooRealVar mean_ctl("mean_ctl", "mean_ctl", -3, -8, 8);
52 RooGaussian gx_ctl("gx_ctl", "gx_ctl", x, mean_ctl, sigma);
53
54 // Construct the background pdf
55 RooRealVar a0_ctl("a0_ctl", "a0_ctl", -0.1, -1, 1);
56 RooRealVar a1_ctl("a1_ctl", "a1_ctl", 0.5, -0.1, 1);
57 RooChebychev px_ctl("px_ctl", "px_ctl", x, RooArgSet(a0_ctl, a1_ctl));
58
59 // Construct the composite model
60 RooRealVar f_ctl("f_ctl", "f_ctl", 0.5, 0., 1.);
61 RooAddPdf model_ctl("model_ctl", "model_ctl", RooArgList(gx_ctl, px_ctl), f_ctl);
62
63 // 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
64 // ---------------------------------------------------------------
65
66 // Generate 1000 events in x and y from model
67 RooDataSet *data = model.generate(RooArgSet(x), 100);
68 RooDataSet *data_ctl = model_ctl.generate(RooArgSet(x), 2000);
69
70 // 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
71 // ---------------------------------------------------------------------------
72
73 // Define category to distinguish physics and control samples events
74 RooCategory sample("sample", "sample");
75 sample.defineType("physics");
76 sample.defineType("control");
77
78 // Construct combined dataset in (x,sample)
79 RooDataSet combData("combData", "combined data", x, Index(sample), Import("physics", *data),
80 Import("control", *data_ctl));
81
82 // 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 )
83 // -----------------------------------------------------------------------------------
84
85 // Construct a simultaneous pdf using category sample as index
86 RooSimultaneous simPdf("simPdf", "simultaneous pdf", sample);
87
88 // Associate model with the physics state and model_ctl with the control state
89 simPdf.addPdf(model, "physics");
90 simPdf.addPdf(model_ctl, "control");
91
92 // P e r f o r m a s i m u l t a n e o u s f i t
93 // ---------------------------------------------------
94
95 // Perform simultaneous fit of model to data and model_ctl to data_ctl
96 simPdf.fitTo(combData);
97
98 // 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
99 // ----------------------------------------------------------------
100
101 // Make a frame for the physics sample
102 RooPlot *frame1 = x.frame(Bins(30), Title("Physics sample"));
103
104 // Plot all data tagged as physics sample
105 combData.plotOn(frame1, Cut("sample==sample::physics"));
106
107 // Plot "physics" slice of simultaneous pdf.
108 // NBL You _must_ project the sample index category with data using ProjWData
109 // as a RooSimultaneous makes no prediction on the shape in the index category
110 // and can thus not be integrated
111 simPdf.plotOn(frame1, Slice(sample, "physics"), ProjWData(sample, combData));
112 simPdf.plotOn(frame1, Slice(sample, "physics"), Components("px"), ProjWData(sample, combData), LineStyle(kDashed));
113
114 // The same plot for the control sample slice
115 RooPlot *frame2 = x.frame(Bins(30), Title("Control sample"));
116 combData.plotOn(frame2, Cut("sample==sample::control"));
117 simPdf.plotOn(frame2, Slice(sample, "control"), ProjWData(sample, combData));
118 simPdf.plotOn(frame2, Slice(sample, "control"), Components("px_ctl"), ProjWData(sample, combData),
120
121 TCanvas *c = new TCanvas("rf501_simultaneouspdf", "rf403_simultaneouspdf", 800, 400);
122 c->Divide(2);
123 c->cd(1);
124 gPad->SetLeftMargin(0.15);
125 frame1->GetYaxis()->SetTitleOffset(1.4);
126 frame1->Draw();
127 c->cd(2);
128 gPad->SetLeftMargin(0.15);
129 frame2->GetYaxis()->SetTitleOffset(1.4);
130 frame2->Draw();
131}
#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
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
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:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
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.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
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