165 using namespace RooFit;
166 using namespace RooStats;
168 void FourBinInstructional(
bool doBayesian =
false,
bool doFeldmanCousins =
false,
bool doMCMC =
false)
176 RooRandom::randomGenerator()->SetSeed(4357);
179 RooWorkspace *wspace =
new RooWorkspace(
"wspace");
180 wspace->factory(
"Poisson::on(non[0,1000], sum::splusb(s[40,0,100],b[100,0,300]))");
181 wspace->factory(
"Poisson::off(noff[0,5000], prod::taub(b,tau[5,3,7],rho[1,0,2]))");
182 wspace->factory(
"Poisson::onbar(nonbar[0,10000], bbar[1000,500,2000])");
183 wspace->factory(
"Poisson::offbar(noffbar[0,1000000], prod::lambdaoffbar(bbar, tau))");
184 wspace->factory(
"Gaussian::mcCons(rhonom[1.,0,2], rho, sigma[.2])");
185 wspace->factory(
"PROD::model(on,off,onbar,offbar,mcCons)");
186 wspace->defineSet(
"obs",
"non,noff,nonbar,noffbar,rhonom");
188 wspace->factory(
"Uniform::prior_poi({s})");
189 wspace->factory(
"Uniform::prior_nuis({b,bbar,tau, rho})");
190 wspace->factory(
"PROD::prior(prior_poi,prior_nuis)");
196 wspace->defineSet(
"poi",
"s");
197 wspace->defineSet(
"nuis",
"b,tau,rho,bbar");
215 RooDataSet *data = wspace->pdf(
"model")->generate(*wspace->set(
"obs"), 1);
217 wspace->import(*data);
222 ModelConfig *modelConfig =
new ModelConfig(
"FourBins");
223 modelConfig->SetWorkspace(*wspace);
224 modelConfig->SetPdf(*wspace->pdf(
"model"));
225 modelConfig->SetPriorPdf(*wspace->pdf(
"prior"));
226 modelConfig->SetParametersOfInterest(*wspace->set(
"poi"));
227 modelConfig->SetNuisanceParameters(*wspace->set(
"nuis"));
228 wspace->import(*modelConfig);
229 wspace->writeToFile(
"FourBin.root");
236 ProfileLikelihoodCalculator plc(*data, *modelConfig);
237 plc.SetConfidenceLevel(0.95);
238 LikelihoodInterval *plInt = plc.GetInterval();
239 RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
240 RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
241 plInt->LowerLimit(*wspace->var(
"s"));
242 RooMsgService::instance().setGlobalKillBelow(msglevel);
245 FeldmanCousins fc(*data, *modelConfig);
246 fc.SetConfidenceLevel(0.95);
248 fc.FluctuateNumDataEntries(
false);
249 fc.UseAdaptiveSampling(
true);
251 PointSetInterval *fcInt = NULL;
252 if (doFeldmanCousins) {
253 fcInt = (PointSetInterval *)fc.GetInterval();
257 BayesianCalculator bc(*data, *modelConfig);
258 bc.SetConfidenceLevel(0.95);
259 SimpleInterval *bInt = NULL;
260 if (doBayesian && wspace->set(
"poi")->getSize() == 1) {
261 bInt = bc.GetInterval();
263 cout <<
"Bayesian Calc. only supports on parameter of interest" << endl;
269 RooFitResult *fit = wspace->pdf(
"model")->fitTo(*data, Save());
271 ph.SetVariables((RooArgSet &)fit->floatParsFinal());
272 ph.SetCovMatrix(fit->covarianceMatrix());
273 ph.SetUpdateProposalParameters(kTRUE);
274 ph.SetCacheSize(100);
275 ProposalFunction *pf = ph.GetProposalFunction();
277 MCMCCalculator mc(*data, *modelConfig);
278 mc.SetConfidenceLevel(0.95);
279 mc.SetProposalFunction(*pf);
280 mc.SetNumBurnInSteps(500);
281 mc.SetNumIters(50000);
282 mc.SetLeftSideTailFraction(0.5);
283 MCMCInterval *mcInt = NULL;
285 mcInt = mc.GetInterval();
289 TCanvas *c1 = (TCanvas *)gROOT->Get(
"c1");
291 c1 =
new TCanvas(
"c1");
293 if (doBayesian && doMCMC) {
296 }
else if (doBayesian || doMCMC) {
301 LikelihoodIntervalPlot *lrplot =
new LikelihoodIntervalPlot(plInt);
304 if (doBayesian && wspace->set(
"poi")->getSize() == 1) {
308 bc.SetScanOfPosterior(20);
309 RooPlot *bplot = bc.GetPosteriorPlot();
314 if (doBayesian && wspace->set(
"poi")->getSize() == 1)
318 MCMCIntervalPlot mcPlot(*mcInt);
324 cout <<
"Profile Likelihood interval on s = [" << plInt->LowerLimit(*wspace->var(
"s")) <<
", "
325 << plInt->UpperLimit(*wspace->var(
"s")) <<
"]" << endl;
328 if (doBayesian && wspace->set(
"poi")->getSize() == 1) {
329 cout <<
"Bayesian interval on s = [" << bInt->LowerLimit() <<
", " << bInt->UpperLimit() <<
"]" << endl;
332 if (doFeldmanCousins) {
333 cout <<
"Feldman Cousins interval on s = [" << fcInt->LowerLimit(*wspace->var(
"s")) <<
", "
334 << fcInt->UpperLimit(*wspace->var(
"s")) <<
"]" << endl;
339 cout <<
"MCMC interval on s = [" << mcInt->LowerLimit(*wspace->var(
"s")) <<
", "
340 << mcInt->UpperLimit(*wspace->var(
"s")) <<
"]" << endl;