Logo ROOT   6.18/05
Reference Guide
JeffreysPriorDemo.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook -js
4/// tutorial demonstrating and validates the RooJeffreysPrior class
5///
6/// Jeffreys's prior is an 'objective prior' based on formal rules.
7/// It is calculated from the Fisher information matrix.
8///
9/// Read more:
10/// http://en.wikipedia.org/wiki/Jeffreys_prior
11///
12/// The analytic form is not known for most PDFs, but it is for
13/// simple cases like the Poisson mean, Gaussian mean, Gaussian sigma.
14///
15/// This class uses numerical tricks to calculate the Fisher Information Matrix
16/// efficiently. In particular, it takes advantage of a property of the
17/// 'Asimov data' as described in
18/// Asymptotic formulae for likelihood-based tests of new physics
19/// Glen Cowan, Kyle Cranmer, Eilam Gross, Ofer Vitells
20/// http://arxiv.org/abs/arXiv:1007.1727
21///
22/// This Demo has four parts:
23/// 1. TestJeffreysPriorDemo -- validates Poisson mean case 1/sqrt(mu)
24/// 2. TestJeffreysGaussMean -- validates Gaussian mean case
25/// 3. TestJeffreysGaussSigma -- validates Gaussian sigma case 1/sigma
26/// 4. TestJeffreysGaussMeanAndSigma -- demonstrates 2-d example
27///
28/// \macro_image
29/// \macro_output
30/// \macro_code
31///
32/// \author Kyle Cranmer
33
34#include "RooJeffreysPrior.h"
35
36#include "RooWorkspace.h"
37#include "RooDataHist.h"
38#include "RooGenericPdf.h"
39#include "TCanvas.h"
40#include "RooPlot.h"
41#include "RooFitResult.h"
42#include "TMatrixDSym.h"
43#include "RooRealVar.h"
44#include "RooAbsPdf.h"
45#include "RooNumIntConfig.h"
46#include "TH1F.h"
47
48using namespace RooFit;
49
50void JeffreysPriorDemo()
51{
52 RooWorkspace w("w");
53 w.factory("Uniform::u(x[0,1])");
54 w.factory("mu[100,1,200]");
55 w.factory("ExtendPdf::p(u,mu)");
56
57 RooDataHist *asimov = w.pdf("p")->generateBinned(*w.var("x"), ExpectedData());
58
59 RooFitResult *res = w.pdf("p")->fitTo(*asimov, Save(), SumW2Error(kTRUE));
60
61 asimov->Print();
62 res->Print();
63 TMatrixDSym cov = res->covarianceMatrix();
64 cout << "variance = " << (cov.Determinant()) << endl;
65 cout << "stdev = " << sqrt(cov.Determinant()) << endl;
66 cov.Invert();
67 cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
68
69 w.defineSet("poi", "mu");
70 w.defineSet("obs", "x");
71
72 RooJeffreysPrior pi("jeffreys", "jeffreys", *w.pdf("p"), *w.set("poi"), *w.set("obs"));
73
74 RooGenericPdf *test = new RooGenericPdf("test", "test", "1./sqrt(mu)", *w.set("poi"));
75
76 TCanvas *c1 = new TCanvas;
77 RooPlot *plot = w.var("mu")->frame();
78 pi.plotOn(plot);
79 test->plotOn(plot, LineColor(kRed));
80 plot->Draw();
81}
82
83//_________________________________________________
84void TestJeffreysGaussMean()
85{
86 RooWorkspace w("w");
87 w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,0,10])");
88 w.factory("n[10,.1,200]");
89 w.factory("ExtendPdf::p(g,n)");
90 w.var("sigma")->setConstant();
91 w.var("n")->setConstant();
92
93 RooDataHist *asimov = w.pdf("p")->generateBinned(*w.var("x"), ExpectedData());
94
95 RooFitResult *res = w.pdf("p")->fitTo(*asimov, Save(), SumW2Error(kTRUE));
96
97 asimov->Print();
98 res->Print();
99 TMatrixDSym cov = res->covarianceMatrix();
100 cout << "variance = " << (cov.Determinant()) << endl;
101 cout << "stdev = " << sqrt(cov.Determinant()) << endl;
102 cov.Invert();
103 cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
104
105 w.defineSet("poi", "mu");
106 w.defineSet("obs", "x");
107
108 RooJeffreysPrior pi("jeffreys", "jeffreys", *w.pdf("p"), *w.set("poi"), *w.set("obs"));
109
110 const RooArgSet *temp = w.set("poi");
111 pi.getParameters(*temp)->Print();
112
113 // return;
114 RooGenericPdf *test = new RooGenericPdf("test", "test", "1", *w.set("poi"));
115
116 TCanvas *c1 = new TCanvas;
117 RooPlot *plot = w.var("mu")->frame();
118 pi.plotOn(plot);
119 test->plotOn(plot, LineColor(kRed), LineStyle(kDotted));
120 plot->Draw();
121}
122
123//_________________________________________________
124void TestJeffreysGaussSigma()
125{
126 // this one is VERY sensitive
127 // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
128 // and you get really bizarre shapes
129 // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
130 // and the PDF falls off too fast at high sigma
131 RooWorkspace w("w");
132 w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
133 w.factory("n[100,.1,2000]");
134 w.factory("ExtendPdf::p(g,n)");
135 // w.var("sigma")->setConstant();
136 w.var("mu")->setConstant();
137 w.var("n")->setConstant();
138 w.var("x")->setBins(301);
139
140 RooDataHist *asimov = w.pdf("p")->generateBinned(*w.var("x"), ExpectedData());
141
142 RooFitResult *res = w.pdf("p")->fitTo(*asimov, Save(), SumW2Error(kTRUE));
143
144 asimov->Print();
145 res->Print();
146 TMatrixDSym cov = res->covarianceMatrix();
147 cout << "variance = " << (cov.Determinant()) << endl;
148 cout << "stdev = " << sqrt(cov.Determinant()) << endl;
149 cov.Invert();
150 cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
151
152 w.defineSet("poi", "sigma");
153 w.defineSet("obs", "x");
154
155 RooJeffreysPrior pi("jeffreys", "jeffreys", *w.pdf("p"), *w.set("poi"), *w.set("obs"));
156 pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 3);
157
158 const RooArgSet *temp = w.set("poi");
159 pi.getParameters(*temp)->Print();
160
161 RooGenericPdf *test = new RooGenericPdf("test", "test", "sqrt(2.)/sigma", *w.set("poi"));
162
163 TCanvas *c1 = new TCanvas;
164 RooPlot *plot = w.var("sigma")->frame();
165 pi.plotOn(plot);
166 test->plotOn(plot, LineColor(kRed), LineStyle(kDotted));
167 plot->Draw();
168}
169
170//_________________________________________________
171void TestJeffreysGaussMeanAndSigma()
172{
173 // this one is VERY sensitive
174 // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
175 // and you get really bizarre shapes
176 // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
177 // and the PDF falls off too fast at high sigma
178 RooWorkspace w("w");
179 w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
180 w.factory("n[100,.1,2000]");
181 w.factory("ExtendPdf::p(g,n)");
182
183 w.var("n")->setConstant();
184 w.var("x")->setBins(301);
185
186 RooDataHist *asimov = w.pdf("p")->generateBinned(*w.var("x"), ExpectedData());
187
188 RooFitResult *res = w.pdf("p")->fitTo(*asimov, Save(), SumW2Error(kTRUE));
189
190 asimov->Print();
191 res->Print();
192 TMatrixDSym cov = res->covarianceMatrix();
193 cout << "variance = " << (cov.Determinant()) << endl;
194 cout << "stdev = " << sqrt(cov.Determinant()) << endl;
195 cov.Invert();
196 cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
197
198 w.defineSet("poi", "mu,sigma");
199 w.defineSet("obs", "x");
200
201 RooJeffreysPrior pi("jeffreys", "jeffreys", *w.pdf("p"), *w.set("poi"), *w.set("obs"));
202 pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps", 3);
203
204 const RooArgSet *temp = w.set("poi");
205 pi.getParameters(*temp)->Print();
206 // return;
207
208 TCanvas *c1 = new TCanvas;
209 TH1 *Jeff2d = pi.createHistogram("2dJeffreys", *w.var("mu"), Binning(10), YVar(*w.var("sigma"), Binning(10)));
210 Jeff2d->Draw("surf");
211}
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:64
@ kDotted
Definition: TAttLine.h:48
double sqrt(double)
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsData.h:162
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Definition: RooFitResult.h:40
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:66
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
RooJeffreysPrior.
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
The Canvas class.
Definition: TCanvas.h:31
The TH1 histogram class.
Definition: TH1.h:56
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
virtual Double_t Determinant() const
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
return c1
Definition: legend1.C:41
Template specialisation used in RooAbsArg:
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg SumW2Error(Bool_t flag)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg ExpectedData(Bool_t flag=kTRUE)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
static constexpr double pi
Definition: test.py:1