Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 pdfs to describe simultaneous
5/// fits to multiple datasets
6///
7/// \macro_image
8/// \macro_output
9/// \macro_code
10///
11/// \date July 2008
12/// \author Wouter Verkerke
13
14#include "RooRealVar.h"
15#include "RooDataSet.h"
16#include "RooGaussian.h"
17#include "RooConstVar.h"
18#include "RooChebychev.h"
19#include "RooAddPdf.h"
20#include "RooSimultaneous.h"
21#include "RooCategory.h"
22#include "TCanvas.h"
23#include "TAxis.h"
24#include "RooPlot.h"
25using namespace RooFit;
26
28{
29 // 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
30 // -------------------------------------------------------------
31
32 // Create observables
33 RooRealVar x("x", "x", -8, 8);
34
35 // Construct signal pdf
36 RooRealVar mean("mean", "mean", 0, -8, 8);
37 RooRealVar sigma("sigma", "sigma", 0.3, 0.1, 10);
38 RooGaussian gx("gx", "gx", x, mean, sigma);
39
40 // Construct background pdf
41 RooRealVar a0("a0", "a0", -0.1, -1, 1);
42 RooRealVar a1("a1", "a1", 0.004, -1, 1);
43 RooChebychev px("px", "px", x, RooArgSet(a0, a1));
44
45 // Construct composite pdf
46 RooRealVar f("f", "f", 0.2, 0., 1.);
47 RooAddPdf model("model", "model", RooArgList(gx, px), f);
48
49 // 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
50 // --------------------------------------------------------------
51
52 // Construct signal pdf.
53 // NOTE that sigma is shared with the signal sample model
54 RooRealVar mean_ctl("mean_ctl", "mean_ctl", -3, -8, 8);
55 RooGaussian gx_ctl("gx_ctl", "gx_ctl", x, mean_ctl, sigma);
56
57 // Construct the background pdf
58 RooRealVar a0_ctl("a0_ctl", "a0_ctl", -0.1, -1, 1);
59 RooRealVar a1_ctl("a1_ctl", "a1_ctl", 0.5, -0.1, 1);
60 RooChebychev px_ctl("px_ctl", "px_ctl", x, RooArgSet(a0_ctl, a1_ctl));
61
62 // Construct the composite model
63 RooRealVar f_ctl("f_ctl", "f_ctl", 0.5, 0., 1.);
64 RooAddPdf model_ctl("model_ctl", "model_ctl", RooArgList(gx_ctl, px_ctl), f_ctl);
65
66 // 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
67 // ---------------------------------------------------------------
68
69 // Generate 1000 events in x and y from model
70 RooDataSet *data = model.generate(RooArgSet(x), 100);
71 RooDataSet *data_ctl = model_ctl.generate(RooArgSet(x), 2000);
72
73 // 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
74 // ---------------------------------------------------------------------------
75
76 // Define category to distinguish physics and control samples events
77 RooCategory sample("sample", "sample");
78 sample.defineType("physics");
79 sample.defineType("control");
80
81 // Construct combined dataset in (x,sample)
82 RooDataSet combData("combData", "combined data", x, Index(sample), Import("physics", *data),
83 Import("control", *data_ctl));
84
85 // 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 )
86 // -----------------------------------------------------------------------------------
87
88 // Construct a simultaneous pdf using category sample as index
89 RooSimultaneous simPdf("simPdf", "simultaneous pdf", sample);
90
91 // Associate model with the physics state and model_ctl with the control state
92 simPdf.addPdf(model, "physics");
93 simPdf.addPdf(model_ctl, "control");
94
95 // P e r f o r m a s i m u l t a n e o u s f i t
96 // ---------------------------------------------------
97
98 // Perform simultaneous fit of model to data and model_ctl to data_ctl
99 simPdf.fitTo(combData);
100
101 // 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
102 // ----------------------------------------------------------------
103
104 // Make a frame for the physics sample
105 RooPlot *frame1 = x.frame(Bins(30), Title("Physics sample"));
106
107 // Plot all data tagged as physics sample
108 combData.plotOn(frame1, Cut("sample==sample::physics"));
109
110 // Plot "physics" slice of simultaneous pdf.
111 // NBL You _must_ project the sample index category with data using ProjWData
112 // as a RooSimultaneous makes no prediction on the shape in the index category
113 // and can thus not be integrated.
114 // In other words: Since the PDF doesn't know the number of events in the different
115 // category states, it doesn't know how much of each component it has to project out.
116 // This information is read from the data.
117 simPdf.plotOn(frame1, Slice(sample, "physics"), ProjWData(sample, combData));
118 simPdf.plotOn(frame1, Slice(sample, "physics"), Components("px"), ProjWData(sample, combData), LineStyle(kDashed));
119
120 // The same plot for the control sample slice
121 RooPlot *frame2 = x.frame(Bins(30), Title("Control sample"));
122 combData.plotOn(frame2, Cut("sample==sample::control"));
123 simPdf.plotOn(frame2, Slice(sample, "control"), ProjWData(sample, combData));
124 simPdf.plotOn(frame2, Slice(sample, "control"), Components("px_ctl"), ProjWData(sample, combData),
126
127 TCanvas *c = new TCanvas("rf501_simultaneouspdf", "rf403_simultaneouspdf", 800, 400);
128 c->Divide(2);
129 c->cd(1);
130 gPad->SetLeftMargin(0.15);
131 frame1->GetYaxis()->SetTitleOffset(1.4);
132 frame1->Draw();
133 c->cd(2);
134 gPad->SetLeftMargin(0.15);
135 frame2->GetYaxis()->SetTitleOffset(1.4);
136 frame2->Draw();
137}
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
@ kDashed
Definition TAttLine.h:48
#define gPad
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:32
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:29
RooCategory is an object to represent discrete states.
Definition RooCategory.h:27
Chebychev polynomial p.d.f.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:33
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
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:1263
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:249
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:691
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
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:293
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg Index(RooCategory &icat)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg Slice(const RooArgSet &sliceSet)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg LineStyle(Style_t style)
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...
const char * Title
Definition TXMLSetup.cxx:68