40 using namespace RooFit;
41 using namespace RooStats;
44 void rs_numberCountingCombination_expected();
45 void rs_numberCountingCombination_observed();
46 void rs_numberCountingCombination_observedWithTau();
50 void rs_numberCountingCombination(
int flag = 1)
53 rs_numberCountingCombination_expected();
55 rs_numberCountingCombination_observed();
57 rs_numberCountingCombination_observedWithTau();
61 void rs_numberCountingCombination_expected()
74 Double_t s[2] = {20., 10.};
75 Double_t b[2] = {100., 100.};
76 Double_t db[2] = {.0100, .0100};
84 NumberCountingPdfFactory f;
85 RooWorkspace *wspace =
new RooWorkspace();
86 f.AddModel(s, 2, wspace,
"TopLevelPdf",
"masterSignal");
91 f.AddExpData(s, b, db, 2, wspace,
"ExpectedNumberCountingData");
101 RooRealVar *mu = wspace->var(
"masterSignal");
102 RooArgSet *poi =
new RooArgSet(*mu);
103 RooArgSet *nullParams =
new RooArgSet(
"nullParams");
104 nullParams->addClone(*mu);
106 nullParams->setRealValue(
"masterSignal", 0);
110 ProfileLikelihoodCalculator plc(*wspace->data(
"ExpectedNumberCountingData"), *wspace->pdf(
"TopLevelPdf"), *poi, 0.05,
114 HypoTestResult *htr = plc.GetHypoTest();
116 cout <<
"-------------------------------------------------" << endl;
117 cout <<
"The p-value for the null is " << htr->NullPValue() << endl;
118 cout <<
"Corresponding to a significance of " << htr->Significance() << endl;
119 cout <<
"-------------------------------------------------\n\n" << endl;
133 RooArgSet *paramsOfInterest = nullParams;
134 plc.SetParameters(*paramsOfInterest);
135 LikelihoodInterval *lrint = (LikelihoodInterval *)plc.GetInterval();
136 lrint->SetConfidenceLevel(0.95);
141 double lower = lrint->LowerLimit(*mu);
142 double upper = lrint->UpperLimit(*mu);
144 LikelihoodIntervalPlot lrPlot(lrint);
145 lrPlot.SetMaximum(3.);
149 cout <<
"lower limit on master signal = " << lower << endl;
150 cout <<
"upper limit on master signal = " << upper << endl;
157 paramsOfInterest->setRealValue(
"masterSignal", 0.);
158 cout <<
"-------------------------------------------------" << endl;
159 std::cout <<
"Consider this parameter point:" << std::endl;
160 paramsOfInterest->first()->Print();
161 if (lrint->IsInInterval(*paramsOfInterest))
162 std::cout <<
"It IS in the interval." << std::endl;
164 std::cout <<
"It is NOT in the interval." << std::endl;
165 cout <<
"-------------------------------------------------\n\n" << endl;
168 paramsOfInterest->setRealValue(
"masterSignal", 2.);
169 cout <<
"-------------------------------------------------" << endl;
170 std::cout <<
"Consider this parameter point:" << std::endl;
171 paramsOfInterest->first()->Print();
172 if (lrint->IsInInterval(*paramsOfInterest))
173 std::cout <<
"It IS in the interval." << std::endl;
175 std::cout <<
"It is NOT in the interval." << std::endl;
176 cout <<
"-------------------------------------------------\n\n" << endl;
219 void rs_numberCountingCombination_observed()
234 Double_t s[2] = {20., 10.};
242 NumberCountingPdfFactory f;
243 RooWorkspace *wspace =
new RooWorkspace();
244 f.AddModel(s, 2, wspace,
"TopLevelPdf",
"masterSignal");
248 Double_t mainMeas[2] = {123., 117.};
249 Double_t bkgMeas[2] = {111.23, 98.76};
250 Double_t dbMeas[2] = {.011, .0095};
251 f.AddData(mainMeas, bkgMeas, dbMeas, 2, wspace,
"ObservedNumberCountingData");
261 RooRealVar *mu = wspace->var(
"masterSignal");
262 RooArgSet *poi =
new RooArgSet(*mu);
263 RooArgSet *nullParams =
new RooArgSet(
"nullParams");
264 nullParams->addClone(*mu);
266 nullParams->setRealValue(
"masterSignal", 0);
270 ProfileLikelihoodCalculator plc(*wspace->data(
"ObservedNumberCountingData"), *wspace->pdf(
"TopLevelPdf"), *poi, 0.05,
273 wspace->var(
"tau_0")->Print();
274 wspace->var(
"tau_1")->Print();
277 HypoTestResult *htr = plc.GetHypoTest();
278 cout <<
"-------------------------------------------------" << endl;
279 cout <<
"The p-value for the null is " << htr->NullPValue() << endl;
280 cout <<
"Corresponding to a significance of " << htr->Significance() << endl;
281 cout <<
"-------------------------------------------------\n\n" << endl;
295 RooArgSet *paramsOfInterest = nullParams;
296 plc.SetParameters(*paramsOfInterest);
297 LikelihoodInterval *lrint = (LikelihoodInterval *)plc.GetInterval();
298 lrint->SetConfidenceLevel(0.95);
301 cout <<
"lower limit on master signal = " << lrint->LowerLimit(*mu) << endl;
302 cout <<
"upper limit on master signal = " << lrint->UpperLimit(*mu) << endl;
311 void rs_numberCountingCombination_observedWithTau()
326 Double_t s[2] = {20., 10.};
334 NumberCountingPdfFactory f;
335 RooWorkspace *wspace =
new RooWorkspace();
336 f.AddModel(s, 2, wspace,
"TopLevelPdf",
"masterSignal");
340 Double_t mainMeas[2] = {123., 117.};
341 Double_t sideband[2] = {11123., 9876.};
342 Double_t tau[2] = {100., 100.};
343 f.AddDataWithSideband(mainMeas, sideband, tau, 2, wspace,
"ObservedNumberCountingDataWithSideband");
353 RooRealVar *mu = wspace->var(
"masterSignal");
354 RooArgSet *poi =
new RooArgSet(*mu);
355 RooArgSet *nullParams =
new RooArgSet(
"nullParams");
356 nullParams->addClone(*mu);
358 nullParams->setRealValue(
"masterSignal", 0);
362 ProfileLikelihoodCalculator plc(*wspace->data(
"ObservedNumberCountingDataWithSideband"), *wspace->pdf(
"TopLevelPdf"),
363 *poi, 0.05, nullParams);
366 HypoTestResult *htr = plc.GetHypoTest();
367 cout <<
"-------------------------------------------------" << endl;
368 cout <<
"The p-value for the null is " << htr->NullPValue() << endl;
369 cout <<
"Corresponding to a significance of " << htr->Significance() << endl;
370 cout <<
"-------------------------------------------------\n\n" << endl;
384 RooArgSet *paramsOfInterest = nullParams;
385 plc.SetParameters(*paramsOfInterest);
386 LikelihoodInterval *lrint = (LikelihoodInterval *)plc.GetInterval();
387 lrint->SetConfidenceLevel(0.95);
390 cout <<
"lower limit on master signal = " << lrint->LowerLimit(*mu) << endl;
391 cout <<
"upper limit on master signal = " << lrint->UpperLimit(*mu) << endl;