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
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118 #include "TStopwatch.h"
00119 #include "TCanvas.h"
00120 #include "TROOT.h"
00121 #include "RooPlot.h"
00122 #include "RooAbsPdf.h"
00123 #include "RooWorkspace.h"
00124 #include "RooDataSet.h"
00125 #include "RooGlobalFunc.h"
00126 #include "RooFitResult.h"
00127 #include "RooRandom.h"
00128 #include "RooStats/ProfileLikelihoodCalculator.h"
00129 #include "RooStats/LikelihoodInterval.h"
00130 #include "RooStats/LikelihoodIntervalPlot.h"
00131 #include "RooStats/BayesianCalculator.h"
00132 #include "RooStats/MCMCCalculator.h"
00133 #include "RooStats/MCMCInterval.h"
00134 #include "RooStats/MCMCIntervalPlot.h"
00135 #include "RooStats/ProposalHelper.h"
00136 #include "RooStats/SimpleInterval.h"
00137 #include "RooStats/FeldmanCousins.h"
00138 #include "RooStats/PointSetInterval.h"
00139
00140 using namespace RooFit;
00141 using namespace RooStats;
00142
00143 void FourBinInstructional(bool doBayesian=false, bool doFeldmanCousins=false, bool doMCMC=false){
00144
00145
00146 TStopwatch t;
00147 t.Start();
00148
00149
00150 RooRandom::randomGenerator()->SetSeed(4357);
00151
00152
00153 RooWorkspace* wspace = new RooWorkspace("wspace");
00154 wspace->factory("Poisson::on(non[0,1000], sum::splusb(s[40,0,100],b[100,0,300]))");
00155 wspace->factory("Poisson::off(noff[0,5000], prod::taub(b,tau[5,3,7],rho[1,0,2]))");
00156 wspace->factory("Poisson::onbar(nonbar[0,10000], bbar[1000,500,2000])");
00157 wspace->factory("Poisson::offbar(noffbar[0,1000000], prod::lambdaoffbar(bbar, tau))");
00158 wspace->factory("Gaussian::mcCons(rhonom[1.,0,2], rho, sigma[.2])");
00159 wspace->factory("PROD::model(on,off,onbar,offbar,mcCons)");
00160 wspace->defineSet("obs","non,noff,nonbar,noffbar,rhonom");
00161
00162 wspace->factory("Uniform::prior_poi({s})");
00163 wspace->factory("Uniform::prior_nuis({b,bbar,tau, rho})");
00164 wspace->factory("PROD::prior(prior_poi,prior_nuis)");
00165
00166
00167
00168
00169
00170 wspace->defineSet("poi","s");
00171 wspace->defineSet("nuis","b,tau,rho,bbar");
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189 RooDataSet* data = wspace->pdf("model")->generate(*wspace->set("obs"),1);
00190
00191 wspace->import(*data);
00192
00193
00194
00195
00196 ModelConfig* modelConfig = new ModelConfig("FourBins");
00197 modelConfig->SetWorkspace(*wspace);
00198 modelConfig->SetPdf(*wspace->pdf("model"));
00199 modelConfig->SetPriorPdf(*wspace->pdf("prior"));
00200 modelConfig->SetParametersOfInterest(*wspace->set("poi"));
00201 modelConfig->SetNuisanceParameters(*wspace->set("nuis"));
00202 wspace->import(*modelConfig);
00203 wspace->writeToFile("FourBin.root");
00204
00205
00206
00207
00208
00209
00210 ProfileLikelihoodCalculator plc(*data, *modelConfig);
00211 plc.SetConfidenceLevel(0.95);
00212 LikelihoodInterval* plInt = plc.GetInterval();
00213 RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
00214 RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
00215 plInt->LowerLimit( *wspace->var("s") );
00216 RooMsgService::instance().setGlobalKillBelow(msglevel);
00217
00218
00219 FeldmanCousins fc(*data, *modelConfig);
00220 fc.SetConfidenceLevel(0.95);
00221
00222 fc.FluctuateNumDataEntries(false);
00223 fc.UseAdaptiveSampling(true);
00224 fc.SetNBins(40);
00225 PointSetInterval* fcInt = NULL;
00226 if(doFeldmanCousins){
00227 fcInt = (PointSetInterval*) fc.GetInterval();
00228 }
00229
00230
00231
00232 BayesianCalculator bc(*data, *modelConfig);
00233 bc.SetConfidenceLevel(0.95);
00234 SimpleInterval* bInt = NULL;
00235 if(doBayesian && wspace->set("poi")->getSize() == 1) {
00236 bInt = bc.GetInterval();
00237 } else{
00238 cout << "Bayesian Calc. only supports on parameter of interest" << endl;
00239 }
00240
00241
00242
00243
00244
00245 RooFitResult* fit = wspace->pdf("model")->fitTo(*data,Save());
00246 ProposalHelper ph;
00247 ph.SetVariables((RooArgSet&)fit->floatParsFinal());
00248 ph.SetCovMatrix(fit->covarianceMatrix());
00249 ph.SetUpdateProposalParameters(kTRUE);
00250 ph.SetCacheSize(100);
00251 ProposalFunction* pf = ph.GetProposalFunction();
00252
00253 MCMCCalculator mc(*data, *modelConfig);
00254 mc.SetConfidenceLevel(0.95);
00255 mc.SetProposalFunction(*pf);
00256 mc.SetNumBurnInSteps(500);
00257 mc.SetNumIters(50000);
00258 mc.SetLeftSideTailFraction(0.5);
00259 MCMCInterval* mcInt = NULL;
00260 if(doMCMC)
00261 mcInt = mc.GetInterval();
00262
00263
00264
00265 TCanvas* c1 = (TCanvas*) gROOT->Get("c1");
00266 if(!c1)
00267 c1 = new TCanvas("c1");
00268
00269 if(doBayesian && doMCMC){
00270 c1->Divide(3);
00271 c1->cd(1);
00272 }
00273 else if (doBayesian || doMCMC){
00274 c1->Divide(2);
00275 c1->cd(1);
00276 }
00277
00278 LikelihoodIntervalPlot* lrplot = new LikelihoodIntervalPlot(plInt);
00279 lrplot->Draw();
00280
00281 if(doBayesian && wspace->set("poi")->getSize() == 1) {
00282 c1->cd(2);
00283
00284
00285 bc.SetScanOfPosterior(20);
00286 RooPlot* bplot = bc.GetPosteriorPlot();
00287 bplot->Draw();
00288 }
00289
00290 if(doMCMC){
00291 if(doBayesian && wspace->set("poi")->getSize() == 1)
00292 c1->cd(3);
00293 else
00294 c1->cd(2);
00295 MCMCIntervalPlot mcPlot(*mcInt);
00296 mcPlot.Draw();
00297 }
00298
00299
00300
00301 cout << "Profile Likelihood interval on s = [" <<
00302 plInt->LowerLimit( *wspace->var("s") ) << ", " <<
00303 plInt->UpperLimit( *wspace->var("s") ) << "]" << endl;
00304
00305
00306
00307 if(doBayesian && wspace->set("poi")->getSize() == 1) {
00308 cout << "Bayesian interval on s = [" <<
00309 bInt->LowerLimit( ) << ", " <<
00310 bInt->UpperLimit( ) << "]" << endl;
00311 }
00312
00313 if(doFeldmanCousins){
00314 cout << "Feldman Cousins interval on s = [" <<
00315 fcInt->LowerLimit( *wspace->var("s") ) << ", " <<
00316 fcInt->UpperLimit( *wspace->var("s") ) << "]" << endl;
00317
00318 }
00319
00320 if(doMCMC){
00321 cout << "MCMC interval on s = [" <<
00322 mcInt->LowerLimit(*wspace->var("s") ) << ", " <<
00323 mcInt->UpperLimit(*wspace->var("s") ) << "]" << endl;
00324
00325
00326 }
00327
00328 t.Print();
00329
00330
00331 }