JeffreysPriorDemo.C

Go to the documentation of this file.
00001 // tutorial demonstrating and validates the RooJeffreysPrior class
00002 
00003 /*
00004 JeffreysPriorDemo.C
00005 
00006 author Kyle Cranmer
00007 date   Dec. 2010
00008 
00009 This tutorial demonstraites and validates the RooJeffreysPrior class
00010 
00011 Jeffreys's prior is an 'objective prior' based on formal rules.  
00012 It is calculated from the Fisher information matrix.  
00013 
00014 Read more:
00015 http://en.wikipedia.org/wiki/Jeffreys_prior
00016 
00017 The analytic form is not known for most PDFs, but it is for 
00018 simple cases like the Poisson mean, Gaussian mean, Gaussian sigma.
00019 
00020 This class uses numerical tricks to calculate the Fisher Information Matrix
00021 efficiently.  In particular, it takes advantage of a property of the
00022 'Asimov data' as described in 
00023 Asymptotic formulae for likelihood-based tests of new physics
00024 Glen Cowan, Kyle Cranmer, Eilam Gross, Ofer Vitells
00025 http://arxiv.org/abs/arXiv:1007.1727
00026 
00027 This Demo has four parts:
00028 TestJeffreysPriorDemo -- validates Poisson mean case 1/sqrt(mu)
00029 TestJeffreysGaussMean -- validates Gaussian mean case
00030 TestJeffreysGaussSigma -- validates Gaussian sigma case 1/sigma
00031 TestJeffreysGaussMeanAndSigma -- demonstraites 2-d example
00032 
00033 */
00034 
00035 #include "RooJeffreysPrior.h"
00036 
00037 #include "RooWorkspace.h"
00038 #include "RooDataHist.h"
00039 #include "RooGenericPdf.h"
00040 #include "TCanvas.h"
00041 #include "RooPlot.h"
00042 #include "RooFitResult.h"
00043 #include "TMatrixDSym.h"
00044 #include "RooRealVar.h"
00045 #include "RooAbsPdf.h"
00046 #include "RooNumIntConfig.h"
00047 #include "TH1F.h"
00048 
00049 using namespace RooFit;
00050 
00051 void JeffreysPriorDemo(){
00052   RooWorkspace w("w");
00053   w.factory("Uniform::u(x[0,1])");
00054   w.factory("mu[100,1,200]");
00055   w.factory("ExtendPdf::p(u,mu)");
00056 
00057   //  w.factory("Poisson::pois(n[0,inf],mu)");
00058 
00059   RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
00060   //  RooDataHist* asimov2 = w.pdf("pois")->generateBinned(*w.var("n"),ExpectedData());
00061 
00062   RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
00063 
00064   asimov->Print();
00065   res->Print();
00066   TMatrixDSym cov = res->covarianceMatrix();
00067   cout << "variance = " << (cov.Determinant()) << endl;
00068   cout << "stdev = " << sqrt(cov.Determinant()) << endl;
00069   cov.Invert();
00070   cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
00071 
00072   w.defineSet("poi","mu");
00073   w.defineSet("obs","x");
00074   //  w.defineSet("obs2","n");
00075 
00076   RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00077   //  pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")  ;
00078   //  pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",10);
00079 
00080   //  JeffreysPrior pi2("jeffreys2","jeffreys",*w.pdf("pois"),*w.set("poi"),*w.set("obs2"));
00081 
00082   //  return;
00083   RooGenericPdf* test = new RooGenericPdf("test","test","1./sqrt(mu)",*w.set("poi"));
00084   
00085   TCanvas* c1 = new TCanvas;
00086   RooPlot* plot = w.var("mu")->frame();
00087   //  pi.plotOn(plot, Normalization(1,RooAbsReal::Raw),Precision(.1));
00088   pi.plotOn(plot);
00089   //  pi2.plotOn(plot,LineColor(kGreen),LineStyle(kDotted));
00090   test->plotOn(plot,LineColor(kRed));
00091   plot->Draw();
00092  
00093 }
00094 
00095 
00096 //_________________________________________________
00097 void TestJeffreysGaussMean(){
00098   RooWorkspace w("w");
00099   w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,0,10])");
00100   w.factory("n[10,.1,200]");
00101   w.factory("ExtendPdf::p(g,n)");
00102   w.var("sigma")->setConstant();
00103   w.var("n")->setConstant();
00104 
00105   RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
00106 
00107   RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
00108 
00109   asimov->Print();
00110   res->Print();
00111   TMatrixDSym cov = res->covarianceMatrix();
00112   cout << "variance = " << (cov.Determinant()) << endl;
00113   cout << "stdev = " << sqrt(cov.Determinant()) << endl;
00114   cov.Invert();
00115   cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
00116 
00117   //  w.defineSet("poi","mu,sigma");
00118   w.defineSet("poi","mu");
00119   w.defineSet("obs","x");
00120 
00121   RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00122   //  pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")  ;
00123   //  pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
00124 
00125   const RooArgSet* temp = w.set("poi");
00126   pi.getParameters(*temp)->Print();
00127 
00128   //  return;
00129   RooGenericPdf* test = new RooGenericPdf("test","test","1",*w.set("poi"));
00130   
00131   TCanvas* c1 = new TCanvas;
00132   RooPlot* plot = w.var("mu")->frame();
00133   pi.plotOn(plot);
00134   test->plotOn(plot,LineColor(kRed),LineStyle(kDotted));
00135   plot->Draw();
00136 
00137   
00138 }
00139 
00140 //_________________________________________________
00141 void TestJeffreysGaussSigma(){
00142   // this one is VERY sensitive
00143   // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
00144   //   and you get really bizzare shapes
00145   // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
00146   //   and the PDF falls off too fast at high sigma
00147   RooWorkspace w("w");
00148   w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
00149   w.factory("n[100,.1,2000]");
00150   w.factory("ExtendPdf::p(g,n)");
00151   //  w.var("sigma")->setConstant();
00152   w.var("mu")->setConstant();
00153   w.var("n")->setConstant();
00154   w.var("x")->setBins(301);
00155 
00156   RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
00157 
00158   RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
00159 
00160   asimov->Print();
00161   res->Print();
00162   TMatrixDSym cov = res->covarianceMatrix();
00163   cout << "variance = " << (cov.Determinant()) << endl;
00164   cout << "stdev = " << sqrt(cov.Determinant()) << endl;
00165   cov.Invert();
00166   cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
00167 
00168 
00169   //  w.defineSet("poi","mu,sigma");
00170   //w.defineSet("poi","mu,sigma,n");
00171   w.defineSet("poi","sigma");
00172   w.defineSet("obs","x");
00173 
00174   RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00175   //  pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")  ;
00176   pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
00177 
00178   const RooArgSet* temp = w.set("poi");
00179   pi.getParameters(*temp)->Print();
00180   //  return;
00181 
00182   //  return;
00183   RooGenericPdf* test = new RooGenericPdf("test","test","sqrt(2.)/sigma",*w.set("poi"));
00184   
00185   TCanvas* c1 = new TCanvas;
00186   RooPlot* plot = w.var("sigma")->frame();
00187   pi.plotOn(plot);
00188   test->plotOn(plot,LineColor(kRed),LineStyle(kDotted));
00189   plot->Draw();
00190 
00191   
00192 }
00193 
00194 
00195 //_________________________________________________
00196 void TestJeffreysGaussMeanAndSigma(){
00197   // this one is VERY sensitive
00198   // if the Gaussian is narrow ~ range(x)/nbins(x) then the peak isn't resolved
00199   //   and you get really bizzare shapes
00200   // if the Gaussian is too wide range(x) ~ sigma then PDF gets renormalized
00201   //   and the PDF falls off too fast at high sigma
00202   RooWorkspace w("w");
00203   w.factory("Gaussian::g(x[0,-20,20],mu[0,-5,5],sigma[1,1,5])");
00204   w.factory("n[100,.1,2000]");
00205   w.factory("ExtendPdf::p(g,n)");
00206   //  w.var("sigma")->setConstant();
00207   //  w.var("mu")->setConstant();
00208   w.var("n")->setConstant();
00209   w.var("x")->setBins(301);
00210 
00211   RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
00212 
00213   RooFitResult* res = w.pdf("p")->fitTo(*asimov,Save(),SumW2Error(kTRUE));
00214 
00215   asimov->Print();
00216   res->Print();
00217   TMatrixDSym cov = res->covarianceMatrix();
00218   cout << "variance = " << (cov.Determinant()) << endl;
00219   cout << "stdev = " << sqrt(cov.Determinant()) << endl;
00220   cov.Invert();
00221   cout << "jeffreys = " << sqrt(cov.Determinant()) << endl;
00222 
00223 
00224   w.defineSet("poi","mu,sigma");
00225   //w.defineSet("poi","mu,sigma,n");
00226   //  w.defineSet("poi","sigma");
00227   w.defineSet("obs","x");
00228 
00229   RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00230   //  pi.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooAdaptiveGaussKronrodIntegrator1D")  ;
00231   pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
00232 
00233   const RooArgSet* temp = w.set("poi");
00234   pi.getParameters(*temp)->Print();
00235   //  return;
00236 
00237   TCanvas* c1 = new TCanvas;
00238   TH1* Jeff2d = pi.createHistogram("2dJeffreys",*w.var("mu"),Binning(10),YVar(*w.var("sigma"),Binning(10)));
00239   Jeff2d->Draw("surf");
00240   
00241 }
00242 
00243 

Generated on Tue Jul 5 15:45:09 2011 for ROOT_528-00b_version by  doxygen 1.5.1