ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
JeffreysPriorDemo.C
Go to the documentation of this file.
1 // tutorial demonstrating and validates the RooJeffreysPrior class
2 
3 /*
4 JeffreysPriorDemo.C
5 
6 author Kyle Cranmer
7 date Dec. 2010
8 
9 This tutorial demonstraites and validates the RooJeffreysPrior class
10 
11 Jeffreys's prior is an 'objective prior' based on formal rules.
12 It is calculated from the Fisher information matrix.
13 
14 Read more:
15 http://en.wikipedia.org/wiki/Jeffreys_prior
16 
17 The analytic form is not known for most PDFs, but it is for
18 simple cases like the Poisson mean, Gaussian mean, Gaussian sigma.
19 
20 This class uses numerical tricks to calculate the Fisher Information Matrix
21 efficiently. In particular, it takes advantage of a property of the
22 'Asimov data' as described in
23 Asymptotic formulae for likelihood-based tests of new physics
24 Glen Cowan, Kyle Cranmer, Eilam Gross, Ofer Vitells
25 http://arxiv.org/abs/arXiv:1007.1727
26 
27 This Demo has four parts:
28 TestJeffreysPriorDemo -- validates Poisson mean case 1/sqrt(mu)
29 TestJeffreysGaussMean -- validates Gaussian mean case
30 TestJeffreysGaussSigma -- validates Gaussian sigma case 1/sigma
31 TestJeffreysGaussMeanAndSigma -- demonstraites 2-d example
32 
33 */
34 
35 #include "RooJeffreysPrior.h"
36 
37 #include "RooWorkspace.h"
38 #include "RooDataHist.h"
39 #include "RooGenericPdf.h"
40 #include "TCanvas.h"
41 #include "RooPlot.h"
42 #include "RooFitResult.h"
43 #include "TMatrixDSym.h"
44 #include "RooRealVar.h"
45 #include "RooAbsPdf.h"
46 #include "RooNumIntConfig.h"
47 #include "TH1F.h"
48 
49 using namespace RooFit;
50 
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  // w.factory("Poisson::pois(n[0,inf],mu)");
58 
59  RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
60  // RooDataHist* asimov2 = w.pdf("pois")->generateBinned(*w.var("n"),ExpectedData());
61 
62  RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
63 
64  asimov->Print();
65  res->Print();
66  TMatrixDSym cov = res->covarianceMatrix();
67  cout << "variance = " << (cov.Determinant()) << endl;
68  cout << "stdev = " << sqrt(cov.Determinant()) << endl;
69  cov.Invert();
70  cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
71 
72  w.defineSet("poi","mu");
73  w.defineSet("obs","x");
74  // w.defineSet("obs2","n");
75 
76  RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
77  // pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
78  // pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",10);
79 
80  // JeffreysPrior pi2("jeffreys2","jeffreys",*w.pdf("pois"),*w.set("poi"),*w.set("obs2"));
81 
82  // return;
83  RooGenericPdf* test = new RooGenericPdf("test","test","1./sqrt(mu)",*w.set("poi"));
84 
85  TCanvas* c1 = new TCanvas;
86  RooPlot* plot = w.var("mu")->frame();
87  // pi.plotOn(plot, Normalization(1,RooAbsReal::Raw),Precision(.1));
88  pi.plotOn(plot);
89  // pi2.plotOn(plot,LineColor(kGreen),LineStyle(kDotted));
90  test->plotOn(plot,LineColor(kRed));
91  plot->Draw();
92 
93 }
94 
95 
96 //_________________________________________________
98  RooWorkspace w("w");
99  w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,0,10])");
100  w.factory("n[10,.1,200]");
101  w.factory("ExtendPdf::p(g,n)");
102  w.var("sigma")->setConstant();
103  w.var("n")->setConstant();
104 
105  RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
106 
107  RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
108 
109  asimov->Print();
110  res->Print();
111  TMatrixDSym cov = res->covarianceMatrix();
112  cout << "variance = " << (cov.Determinant()) << endl;
113  cout << "stdev = " << sqrt(cov.Determinant()) << endl;
114  cov.Invert();
115  cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
116 
117  // w.defineSet("poi","mu,sigma");
118  w.defineSet("poi","mu");
119  w.defineSet("obs","x");
120 
121  RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
122  // pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
123  // pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
124 
125  const RooArgSet* temp = w.set("poi");
126  pi.getParameters(*temp)->Print();
127 
128  // return;
129  RooGenericPdf* test = new RooGenericPdf("test","test","1",*w.set("poi"));
130 
131  TCanvas* c1 = new TCanvas;
132  RooPlot* plot = w.var("mu")->frame();
133  pi.plotOn(plot);
134  test->plotOn(plot,LineColor(kRed),LineStyle(kDotted));
135  plot->Draw();
136 
137 
138 }
139 
140 //_________________________________________________
142  // this one is VERY sensitive
143  // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
144  // and you get really bizzare shapes
145  // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
146  // and the PDF falls off too fast at high sigma
147  RooWorkspace w("w");
148  w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
149  w.factory("n[100,.1,2000]");
150  w.factory("ExtendPdf::p(g,n)");
151  // w.var("sigma")->setConstant();
152  w.var("mu")->setConstant();
153  w.var("n")->setConstant();
154  w.var("x")->setBins(301);
155 
156  RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
157 
158  RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
159 
160  asimov->Print();
161  res->Print();
162  TMatrixDSym cov = res->covarianceMatrix();
163  cout << "variance = " << (cov.Determinant()) << endl;
164  cout << "stdev = " << sqrt(cov.Determinant()) << endl;
165  cov.Invert();
166  cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
167 
168 
169  // w.defineSet("poi","mu,sigma");
170  //w.defineSet("poi","mu,sigma,n");
171  w.defineSet("poi","sigma");
172  w.defineSet("obs","x");
173 
174  RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
175  // pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
176  pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
177 
178  const RooArgSet* temp = w.set("poi");
179  pi.getParameters(*temp)->Print();
180  // return;
181 
182  // return;
183  RooGenericPdf* test = new RooGenericPdf("test","test","sqrt(2.)/sigma",*w.set("poi"));
184 
185  TCanvas* c1 = new TCanvas;
186  RooPlot* plot = w.var("sigma")->frame();
187  pi.plotOn(plot);
188  test->plotOn(plot,LineColor(kRed),LineStyle(kDotted));
189  plot->Draw();
190 
191 
192 }
193 
194 
195 //_________________________________________________
197  // this one is VERY sensitive
198  // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
199  // and you get really bizzare shapes
200  // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
201  // and the PDF falls off too fast at high sigma
202  RooWorkspace w("w");
203  w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
204  w.factory("n[100,.1,2000]");
205  w.factory("ExtendPdf::p(g,n)");
206  // w.var("sigma")->setConstant();
207  // w.var("mu")->setConstant();
208  w.var("n")->setConstant();
209  w.var("x")->setBins(301);
210 
211  RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
212 
213  RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
214 
215  asimov->Print();
216  res->Print();
217  TMatrixDSym cov = res->covarianceMatrix();
218  cout << "variance = " << (cov.Determinant()) << endl;
219  cout << "stdev = " << sqrt(cov.Determinant()) << endl;
220  cov.Invert();
221  cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
222 
223 
224  w.defineSet("poi","mu,sigma");
225  //w.defineSet("poi","mu,sigma,n");
226  // w.defineSet("poi","sigma");
227  w.defineSet("obs","x");
228 
229  RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
230  // pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D") ;
231  pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
232 
233  const RooArgSet* temp = w.set("poi");
234  pi.getParameters(*temp)->Print();
235  // return;
236 
237  TCanvas* c1 = new TCanvas;
238  TH1* Jeff2d = pi.createHistogram("2dJeffreys",*w.var("mu"),Binning(10),YVar(*w.var("sigma"),Binning(10)));
239  Jeff2d->Draw("surf");
240 
241 }
242 
243 
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooFitResult.h:65
void TestJeffreysGaussMean()
void TestJeffreysGaussSigma()
RooCmdArg LineColor(Color_t color)
const double pi
virtual RooDataHist * generateBinned(const RooArgSet &whatVars, Double_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:2175
TCanvas * c1
Definition: legend1.C:2
Definition: Rtypes.h:61
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsData.h:157
void plot(TString fname="data.root", TString var0="var0", TString var1="var1")
Definition: createData.C:17
double sqrt(double)
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
virtual Double_t Determinant() const
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooCmdArg LineStyle(Style_t style)
void setBins(Int_t nBins, const char *name=0)
Definition: RooRealVar.h:78
RooCmdArg ExpectedData(Bool_t flag=kTRUE)
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual RooPlot * plotOn(RooPlot *frame, 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(), const RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Plot (project) PDF on specified frame.
Definition: RooAbsPdf.h:105
void TestJeffreysGaussMeanAndSigma()
void setConstant(Bool_t value=kTRUE)
Bool_t defineSet(const char *name, const RooArgSet &aset, Bool_t importMissing=kFALSE)
Define a named RooArgSet with given constituents.
const TMatrixDSym & covarianceMatrix() const
Return covariance matrix.
tuple w
Definition: qtexample.py:51
RooPlot * frame(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()) const
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found...
RooCmdArg SumW2Error(Bool_t flag)
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
The TH1 histogram class.
Definition: TH1.h:80
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
RooFactoryWSTool & factory()
Return instance to factory tool.
RooCmdArg Save(Bool_t flag=kTRUE)
RooJeffreysPrior.
void JeffreysPriorDemo()
const RooArgSet * set(const char *name)
Return pointer to previously defined named set with given nmame If no such set is found a null pointe...
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
virtual RooFitResult * fitTo(RooAbsData &data, 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())
Fit PDF to given dataset.
Definition: RooAbsPdf.cxx:1056
const Bool_t kTRUE
Definition: Rtypes.h:91
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:42
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559
RooCmdArg Binning(const RooAbsBinning &binning)