00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
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
00058
00059 RooDataHist* asimov = w.pdf("p")->generateBinned(*w.var("x"),ExpectedData());
00060
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
00075
00076 RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00077
00078
00079
00080
00081
00082
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
00088 pi.plotOn(plot);
00089
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
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
00123
00124
00125 const RooArgSet* temp = w.set("poi");
00126 pi.getParameters(*temp)->Print();
00127
00128
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
00143
00144
00145
00146
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
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
00170
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
00176 pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
00177
00178 const RooArgSet* temp = w.set("poi");
00179 pi.getParameters(*temp)->Print();
00180
00181
00182
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
00198
00199
00200
00201
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
00207
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
00226
00227 w.defineSet("obs","x");
00228
00229 RooJeffreysPrior pi("jeffreys","jeffreys",*w.pdf("p"),*w.set("poi"),*w.set("obs"));
00230
00231 pi.specialIntegratorConfig(kTRUE)->getConfigSection("RooIntegrator1D").setRealValue("maxSteps",3);
00232
00233 const RooArgSet* temp = w.set("poi");
00234 pi.getParameters(*temp)->Print();
00235
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