47 using namespace RooFit;
48 using namespace RooStats;
50 void rs101_limitexample()
62 RooWorkspace *wspace =
new RooWorkspace();
63 wspace->factory(
"Poisson::countingModel(obs[150,0,300], "
64 "sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))");
67 wspace->factory(
"Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)");
68 wspace->factory(
"Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)");
69 wspace->factory(
"PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)");
72 RooAbsPdf *modelWithConstraints = wspace->pdf(
"modelWithConstraints");
73 RooRealVar *obs = wspace->var(
"obs");
74 RooRealVar *s = wspace->var(
"s");
79 RooRealVar *ratioSigEff = wspace->var(
"ratioSigEff");
80 RooRealVar *ratioBkgEff = wspace->var(
"ratioBkgEff");
81 RooArgSet constrainedParams(*ratioSigEff,
84 RooRealVar *gSigEff = wspace->var(
"gSigEff");
85 RooRealVar *gSigBkg = wspace->var(
"gSigBkg");
86 gSigEff->setConstant();
87 gSigBkg->setConstant();
91 RooDataSet *data =
new RooDataSet(
"exampleData",
"exampleData", RooArgSet(*obs));
94 RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
97 modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));
100 RooArgSet paramOfInterest(*s);
102 ModelConfig modelConfig(wspace);
103 modelConfig.SetPdf(*modelWithConstraints);
104 modelConfig.SetParametersOfInterest(paramOfInterest);
105 modelConfig.SetNuisanceParameters(constrainedParams);
106 modelConfig.SetObservables(*obs);
107 modelConfig.SetGlobalObservables(RooArgSet(*gSigEff, *gSigBkg));
108 modelConfig.SetName(
"ModelConfig");
109 wspace->import(modelConfig);
110 wspace->import(*data);
111 wspace->SetName(
"w");
112 wspace->writeToFile(
"rs101_ws.root");
116 ProfileLikelihoodCalculator plc(*data, modelConfig);
117 plc.SetTestSize(.05);
118 ConfInterval *lrinterval = plc.GetInterval();
121 TCanvas *dataCanvas =
new TCanvas(
"dataCanvas");
122 dataCanvas->Divide(2, 1);
125 LikelihoodIntervalPlot plotInt((LikelihoodInterval *)lrinterval);
126 plotInt.SetTitle(
"Profile Likelihood Ratio and Posterior for S");
130 FeldmanCousins fc(*data, modelConfig);
131 fc.UseAdaptiveSampling(
true);
132 fc.FluctuateNumDataEntries(
false);
136 ConfInterval *fcint = NULL;
137 fcint = fc.GetInterval();
139 RooFitResult *fit = modelWithConstraints->fitTo(*data, Save(
true));
145 ph.SetVariables((RooArgSet &)fit->floatParsFinal());
146 ph.SetCovMatrix(fit->covarianceMatrix());
147 ph.SetUpdateProposalParameters(
true);
148 ph.SetCacheSize(100);
149 ProposalFunction *pdfProp = ph.GetProposalFunction();
151 MCMCCalculator mc(*data, modelConfig);
152 mc.SetNumIters(20000);
154 mc.SetNumBurnInSteps(40);
155 mc.SetProposalFunction(*pdfProp);
156 mc.SetLeftSideTailFraction(0.5);
157 MCMCInterval *mcInt = (MCMCInterval *)mc.GetInterval();
160 cout <<
"Profile lower limit on s = " << ((LikelihoodInterval *)lrinterval)->LowerLimit(*s) << endl;
161 cout <<
"Profile upper limit on s = " << ((LikelihoodInterval *)lrinterval)->UpperLimit(*s) << endl;
165 double fcul = ((PointSetInterval *)fcint)->UpperLimit(*s);
166 double fcll = ((PointSetInterval *)fcint)->LowerLimit(*s);
167 cout <<
"FC lower limit on s = " << fcll << endl;
168 cout <<
"FC upper limit on s = " << fcul << endl;
169 TLine *fcllLine =
new TLine(fcll, 0, fcll, 1);
170 TLine *fculLine =
new TLine(fcul, 0, fcul, 1);
171 fcllLine->SetLineColor(kRed);
172 fculLine->SetLineColor(kRed);
173 fcllLine->Draw(
"same");
174 fculLine->Draw(
"same");
175 dataCanvas->Update();
179 MCMCIntervalPlot mcPlot(*mcInt);
180 mcPlot.SetLineColor(kMagenta);
181 mcPlot.SetLineWidth(2);
184 double mcul = mcInt->UpperLimit(*s);
185 double mcll = mcInt->LowerLimit(*s);
186 cout <<
"MCMC lower limit on s = " << mcll << endl;
187 cout <<
"MCMC upper limit on s = " << mcul << endl;
188 cout <<
"MCMC Actual confidence level: " << mcInt->GetActualConfidenceLevel() << endl;
193 RooDataSet *chainData = mcInt->GetChainAsDataSet();
196 std::cout <<
"plotting the chain data - nentries = " << chainData->numEntries() << std::endl;
197 TTree *chain = RooStats::GetAsTTree(
"chainTreeData",
"chainTreeData", *chainData);
199 chain->SetMarkerStyle(6);
200 chain->SetMarkerColor(kRed);
202 chain->Draw(
"s:ratioSigEff:ratioBkgEff",
"nll_MarkovChain_local_",
"box");
205 RooDataSet *parScanData = (RooDataSet *)fc.GetPointsToScan();
207 std::cout <<
"plotting the scanned points used in the frequentist construction - npoints = "
208 << parScanData->numEntries() << std::endl;
213 TGraph2D *gr =
new TGraph2D(parScanData->numEntries());
214 for (
int ievt = 0; ievt < parScanData->numEntries(); ++ievt) {
215 const RooArgSet *evt = parScanData->get(ievt);
216 double x = evt->getRealValue(
"ratioBkgEff");
217 double y = evt->getRealValue(
"ratioSigEff");
218 double z = evt->getRealValue(
"s");
219 gr->SetPoint(ievt, x, y, z);
222 gr->SetMarkerStyle(24);