Logo ROOT   6.18/05
Reference Guide
rf801_mcstudy.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Validation and MC studies: toy Monte Carlo study that perform cycles of event generation and fitting
5///
6/// \macro_image
7/// \macro_output
8/// \macro_code
9/// \author
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 "RooMCStudy.h"
18#include "RooPlot.h"
19#include "TCanvas.h"
20#include "TAxis.h"
21#include "TH2.h"
22#include "RooFitResult.h"
23#include "TStyle.h"
24#include "TDirectory.h"
25
26using namespace RooFit;
27
28void rf801_mcstudy()
29{
30 // C r e a t e m o d e l
31 // -----------------------
32
33 // Declare observable x
34 RooRealVar x("x", "x", 0, 10);
35 x.setBins(40);
36
37 // Create two Gaussian PDFs g1(x,mean1,sigma) anf g2(x,mean2,sigma) and their parameters
38 RooRealVar mean("mean", "mean of gaussians", 5, 0, 10);
39 RooRealVar sigma1("sigma1", "width of gaussians", 0.5);
40 RooRealVar sigma2("sigma2", "width of gaussians", 1);
41
42 RooGaussian sig1("sig1", "Signal component 1", x, mean, sigma1);
43 RooGaussian sig2("sig2", "Signal component 2", x, mean, sigma2);
44
45 // Build Chebychev polynomial p.d.f.
46 RooRealVar a0("a0", "a0", 0.5, 0., 1.);
47 RooRealVar a1("a1", "a1", -0.2, -1, 1.);
48 RooChebychev bkg("bkg", "Background", x, RooArgSet(a0, a1));
49
50 // Sum the signal components into a composite signal p.d.f.
51 RooRealVar sig1frac("sig1frac", "fraction of component 1 in signal", 0.8, 0., 1.);
52 RooAddPdf sig("sig", "Signal", RooArgList(sig1, sig2), sig1frac);
53
54 // Sum the composite signal and background
55 RooRealVar nbkg("nbkg", "number of background events,", 150, 0, 1000);
56 RooRealVar nsig("nsig", "number of signal events", 150, 0, 1000);
57 RooAddPdf model("model", "g1+g2+a", RooArgList(bkg, sig), RooArgList(nbkg, nsig));
58
59 // C r e a t e m a n a g e r
60 // ---------------------------
61
62 // Instantiate RooMCStudy manager on model with x as observable and given choice of fit options
63 //
64 // The Silence() option kills all messages below the PROGRESS level, leaving only a single message
65 // per sample executed, and any error message that occur during fitting
66 //
67 // The Extended() option has two effects:
68 // 1) The extended ML term is included in the likelihood and
69 // 2) A poisson fluctuation is introduced on the number of generated events
70 //
71 // The FitOptions() given here are passed to the fitting stage of each toy experiment.
72 // If Save() is specified, the fit result of each experiment is saved by the manager
73 //
74 // A Binned() option is added in this example to bin the data between generation and fitting
75 // to speed up the study at the expense of some precision
76
77 RooMCStudy *mcstudy =
79
80 // G e n e r a t e a n d f i t e v e n t s
81 // ---------------------------------------------
82
83 // Generate and fit 1000 samples of Poisson(nExpected) events
84 mcstudy->generateAndFit(1000);
85
86 // E x p l o r e r e s u l t s o f s t u d y
87 // ------------------------------------------------
88
89 // Make plots of the distributions of mean, the error on mean and the pull of mean
90 RooPlot *frame1 = mcstudy->plotParam(mean, Bins(40));
91 RooPlot *frame2 = mcstudy->plotError(mean, Bins(40));
92 RooPlot *frame3 = mcstudy->plotPull(mean, Bins(40), FitGauss(kTRUE));
93
94 // Plot distribution of minimized likelihood
95 RooPlot *frame4 = mcstudy->plotNLL(Bins(40));
96
97 // Make some histograms from the parameter dataset
98 TH1 *hh_cor_a0_s1f = mcstudy->fitParDataSet().createHistogram("hh", a1, YVar(sig1frac));
99 TH1 *hh_cor_a0_a1 = mcstudy->fitParDataSet().createHistogram("hh", a0, YVar(a1));
100
101 // Access some of the saved fit results from individual toys
102 TH2 *corrHist000 = mcstudy->fitResult(0)->correlationHist("c000");
103 TH2 *corrHist127 = mcstudy->fitResult(127)->correlationHist("c127");
104 TH2 *corrHist953 = mcstudy->fitResult(953)->correlationHist("c953");
105
106 // Draw all plots on a canvas
107 gStyle->SetOptStat(0);
108 TCanvas *c = new TCanvas("rf801_mcstudy", "rf801_mcstudy", 900, 900);
109 c->Divide(3, 3);
110 c->cd(1);
111 gPad->SetLeftMargin(0.15);
112 frame1->GetYaxis()->SetTitleOffset(1.4);
113 frame1->Draw();
114 c->cd(2);
115 gPad->SetLeftMargin(0.15);
116 frame2->GetYaxis()->SetTitleOffset(1.4);
117 frame2->Draw();
118 c->cd(3);
119 gPad->SetLeftMargin(0.15);
120 frame3->GetYaxis()->SetTitleOffset(1.4);
121 frame3->Draw();
122 c->cd(4);
123 gPad->SetLeftMargin(0.15);
124 frame4->GetYaxis()->SetTitleOffset(1.4);
125 frame4->Draw();
126 c->cd(5);
127 gPad->SetLeftMargin(0.15);
128 hh_cor_a0_s1f->GetYaxis()->SetTitleOffset(1.4);
129 hh_cor_a0_s1f->Draw("box");
130 c->cd(6);
131 gPad->SetLeftMargin(0.15);
132 hh_cor_a0_a1->GetYaxis()->SetTitleOffset(1.4);
133 hh_cor_a0_a1->Draw("box");
134 c->cd(7);
135 gPad->SetLeftMargin(0.15);
136 corrHist000->GetYaxis()->SetTitleOffset(1.4);
137 corrHist000->Draw("colz");
138 c->cd(8);
139 gPad->SetLeftMargin(0.15);
140 corrHist127->GetYaxis()->SetTitleOffset(1.4);
141 corrHist127->Draw("colz");
142 c->cd(9);
143 gPad->SetLeftMargin(0.15);
144 corrHist953->GetYaxis()->SetTitleOffset(1.4);
145 corrHist953->Draw("colz");
146
147 // Make RooMCStudy object available on command line after
148 // macro finishes
149 gDirectory->Add(mcstudy);
150}
#define c(i)
Definition: RSha256.hxx:101
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define gDirectory
Definition: TDirectory.h:218
R__EXTERN TStyle * gStyle
Definition: TStyle.h:406
#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
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
TH2F * createHistogram(const RooAbsRealLValue &var1, const RooAbsRealLValue &var2, const char *cuts="", const char *name="hist") const
Create a TH2F histogram of the distribution of the specified variable using this dataset.
TH2 * correlationHist(const char *name="correlation_matrix") const
Return TH2D of correlation matrix.
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooMCStudy is a helper class to facilitate Monte Carlo studies such as 'goodness-of-fit' studies,...
Definition: RooMCStudy.h:32
RooPlot * plotPull(const RooRealVar &param, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Plot the distribution of pull values for the specified parameter on a newly created frame.
RooPlot * plotError(const RooRealVar &param, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Plot the distribution of the fit errors for the specified parameter on a newly created frame.
RooPlot * plotParam(const RooRealVar &param, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Plot the distribution of the fitted value of the given parameter on a newly created frame.
RooPlot * plotNLL(const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Plot the distribution of the -log(L) values on a newly created frame.
const RooFitResult * fitResult(Int_t sampleNum) const
Return the RooFitResult object of the fit to given sample.
const RooDataSet & fitParDataSet()
Return a RooDataSet the resulting fit parameters of each toy cycle.
Definition: RooMCStudy.cxx:966
Bool_t generateAndFit(Int_t nSamples, Int_t nEvtPerSample=0, Bool_t keepGenData=kFALSE, const char *asciiFilePat=0)
Generate and fit 'nSamples' samples of 'nEvtPerSample' events.
Definition: RooMCStudy.cxx:646
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
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
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetYaxis()
Definition: TH1.h:317
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1444
Double_t x[n]
Definition: legend1.C:17
Template specialisation used in RooAbsArg:
RooCmdArg Binned(Bool_t flag=kTRUE)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg Extended(Bool_t flag=kTRUE)
RooCmdArg FitGauss(Bool_t flag=kTRUE)
RooCmdArg PrintEvalErrors(Int_t numErrors)
RooCmdArg Silence(Bool_t flag=kTRUE)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg FitOptions(const char *opts)
RooCmdArg Bins(Int_t nbin)