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 "RooChebychev.h"
18#include "RooAddPdf.h"
19#include "RooSimultaneous.h"
20#include "RooCategory.h"
21#include "TCanvas.h"
22#include "TAxis.h"
23#include "RooPlot.h"
24using namespace RooFit;
25
27{
28 // 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
29 // -------------------------------------------------------------
30
31 // Create observables
32 RooRealVar x("x", "x", -8, 8);
33
34 // Construct signal pdf
35 RooRealVar mean("mean", "mean", 0, -8, 8);
36 RooRealVar sigma("sigma", "sigma", 0.3, 0.1, 10);
37 RooGaussian gx("gx", "gx", x, mean, sigma);
38
39 // Construct background pdf
40 RooRealVar a0("a0", "a0", -0.1, -1, 1);
41 RooRealVar a1("a1", "a1", 0.004, -1, 1);
42 RooChebychev px("px", "px", x, RooArgSet(a0, a1));
43
44 // Construct composite pdf
45 RooRealVar f("f", "f", 0.2, 0., 1.);
46 RooAddPdf model("model", "model", RooArgList(gx, px), f);
47
48 // 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
49 // --------------------------------------------------------------
50
51 // Construct signal pdf.
52 // NOTE that sigma is shared with the signal sample model
53 RooRealVar mean_ctl("mean_ctl", "mean_ctl", -3, -8, 8);
54 RooGaussian gx_ctl("gx_ctl", "gx_ctl", x, mean_ctl, sigma);
55
56 // Construct the background pdf
57 RooRealVar a0_ctl("a0_ctl", "a0_ctl", -0.1, -1, 1);
58 RooRealVar a1_ctl("a1_ctl", "a1_ctl", 0.5, -0.1, 1);
59 RooChebychev px_ctl("px_ctl", "px_ctl", x, RooArgSet(a0_ctl, a1_ctl));
60
61 // Construct the composite model
62 RooRealVar f_ctl("f_ctl", "f_ctl", 0.5, 0., 1.);
63 RooAddPdf model_ctl("model_ctl", "model_ctl", RooArgList(gx_ctl, px_ctl), f_ctl);
64
65 // 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
66 // ---------------------------------------------------------------
67
68 // Generate 1000 events in x and y from model
69 RooDataSet *data = model.generate(RooArgSet(x), 100);
70 RooDataSet *data_ctl = model_ctl.generate(RooArgSet(x), 2000);
71
72 // 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
73 // ---------------------------------------------------------------------------
74
75 // Define category to distinguish physics and control samples events
76 RooCategory sample("sample", "sample");
77 sample.defineType("physics");
78 sample.defineType("control");
79
80 // Construct combined dataset in (x,sample)
81 RooDataSet combData("combData", "combined data", x, Index(sample),
82 Import({{"physics", data}, {"control", data_ctl}}));
83
84 // 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 )
85 // -----------------------------------------------------------------------------------
86
87 // Construct a simultaneous pdf using category sample as index
88 RooSimultaneous simPdf("simPdf", "simultaneous pdf", sample);
89
90 // Associate model with the physics state and model_ctl with the control state
91 simPdf.addPdf(model, "physics");
92 simPdf.addPdf(model_ctl, "control");
93
94 // P e r f o r m a s i m u l t a n e o u s f i t
95 // ---------------------------------------------------
96
97 // Perform simultaneous fit of model to data and model_ctl to data_ctl
98 simPdf.fitTo(combData);
99
100 // 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
101 // ----------------------------------------------------------------
102
103 // Make a frame for the physics sample
104 RooPlot *frame1 = x.frame(Bins(30), Title("Physics sample"));
105
106 // Plot all data tagged as physics sample
107 combData.plotOn(frame1, Cut("sample==sample::physics"));
108
109 // Plot "physics" slice of simultaneous pdf.
110 // NBL You _must_ project the sample index category with data using ProjWData
111 // as a RooSimultaneous makes no prediction on the shape in the index category
112 // and can thus not be integrated.
113 // In other words: Since the PDF doesn't know the number of events in the different
114 // category states, it doesn't know how much of each component it has to project out.
115 // This information is read from the data.
116 simPdf.plotOn(frame1, Slice(sample, "physics"), ProjWData(sample, combData));
117 simPdf.plotOn(frame1, Slice(sample, "physics"), Components("px"), ProjWData(sample, combData), LineStyle(kDashed));
118
119 // The same plot for the control sample slice
120 RooPlot *frame2 = x.frame(Bins(30), Title("Control sample"));
121 combData.plotOn(frame2, Cut("sample==sample::control"));
122 simPdf.plotOn(frame2, Slice(sample, "control"), ProjWData(sample, combData));
123 simPdf.plotOn(frame2, Slice(sample, "control"), Components("px_ctl"), ProjWData(sample, combData),
125
126 TCanvas *c = new TCanvas("rf501_simultaneouspdf", "rf403_simultaneouspdf", 800, 400);
127 c->Divide(2);
128 c->cd(1);
129 gPad->SetLeftMargin(0.15);
130 frame1->GetYaxis()->SetTitleOffset(1.4);
131 frame1->Draw();
132 c->cd(2);
133 gPad->SetLeftMargin(0.15);
134 frame2->GetYaxis()->SetTitleOffset(1.4);
135 frame2->Draw();
136}
#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:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
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:36
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:1278
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:661
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:302
The Canvas class.
Definition TCanvas.h:23
RooCmdArg Bins(Int_t nbin)
RooCmdArg Index(RooCategory &icat)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg Components(Args_t &&... argsOrArgSet)
RooCmdArg ProjWData(const RooAbsData &projData, Bool_t binData=kFALSE)
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...
Definition Common.h:18
const char * Title
Definition TXMLSetup.cxx:68