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