103 using namespace RooFit;
104 using namespace RooStats;
112 class BinCountTestStat :
public TestStatistic {
114 BinCountTestStat(
void) : fColumnName(
"tmp") {}
115 BinCountTestStat(
string columnName) : fColumnName(columnName) {}
117 virtual Double_t Evaluate(RooAbsData &data, RooArgSet & )
120 Double_t value = 0.0;
121 for (
int i = 0; i < data.numEntries(); i++) {
122 value += data.get(i)->getRealValue(fColumnName.c_str());
126 virtual const TString GetVarName()
const {
return fColumnName; }
132 ClassDef(BinCountTestStat, 1)
135 ClassImp(BinCountTestStat)
141 void HybridStandardForm()
162 TCanvas *c =
new TCanvas;
171 RooWorkspace *w =
new RooWorkspace(
"w");
176 w->factory(
"Uniform::f(m[0,1])");
177 w->factory(
"ExtendPdf::px(f,sum::splusb(s[0,0,100],b[100,0,300]))");
178 w->factory(
"Poisson::py(y[100,0,500],prod::taub(tau[1.],b))");
179 w->factory(
"PROD::model(px,py)");
180 w->factory(
"Uniform::prior_b(b)");
185 RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
188 ProofConfig *pc = NULL;
190 pc =
new ProofConfig(*w, 4,
"workers=4", kFALSE);
200 double p_Bi = NumberCountingUtils::BinomialWithTauObsP(150, 100, 1);
201 double Z_Bi = NumberCountingUtils::BinomialWithTauObsZ(150, 100, 1);
202 cout <<
"-----------------------------------------" << endl;
203 cout <<
"Part 3" << endl;
204 std::cout <<
"Z_Bi p-value (analytic): " << p_Bi << std::endl;
205 std::cout <<
"Z_Bi significance (analytic): " << Z_Bi << std::endl;
229 w->defineSet(
"obs",
"m");
230 w->defineSet(
"poi",
"s");
235 RooDataSet *data = w->pdf(
"px")->generate(*w->set(
"obs"), 150);
240 ModelConfig b_model(
"B_model", w);
241 b_model.SetPdf(*w->pdf(
"px"));
242 b_model.SetObservables(*w->set(
"obs"));
243 b_model.SetParametersOfInterest(*w->set(
"poi"));
244 w->var(
"s")->setVal(0.0);
245 b_model.SetSnapshot(*w->set(
"poi"));
248 ModelConfig sb_model(
"S+B_model", w);
249 sb_model.SetPdf(*w->pdf(
"px"));
250 sb_model.SetObservables(*w->set(
"obs"));
251 sb_model.SetParametersOfInterest(*w->set(
"poi"));
252 w->var(
"s")->setVal(50.0);
253 sb_model.SetSnapshot(*w->set(
"poi"));
263 NumEventsTestStat eventCount(*w->pdf(
"px"));
286 w->factory(
"Gaussian::gauss_prior(b,y, expr::sqrty('sqrt(y)',y))");
290 w->factory(
"Lognormal::lognorm_prior(b,y, expr::kappa('1+1./sqrt(y)',y))");
308 HybridCalculator hc1(*data, sb_model, b_model);
309 ToyMCSampler *toymcs1 = (ToyMCSampler *)hc1.GetTestStatSampler();
311 toymcs1->SetTestStatistic(&eventCount);
313 hc1.SetToys(30000, 1000);
314 hc1.ForcePriorNuisanceAlt(*w->pdf(
"py"));
315 hc1.ForcePriorNuisanceNull(*w->pdf(
"py"));
332 RooMsgService::instance().setGlobalKillBelow(RooFit::ERROR);
334 HypoTestResult *r1 = hc1.GetHypoTest();
335 RooMsgService::instance().setGlobalKillBelow(msglevel);
336 cout <<
"-----------------------------------------" << endl;
337 cout <<
"Part 4" << endl;
345 HypoTestPlot *p1 =
new HypoTestPlot(*r1, 30);
357 SimpleLikelihoodRatioTestStat slrts(*b_model.GetPdf(), *sb_model.GetPdf());
358 slrts.SetNullParameters(*b_model.GetSnapshot());
359 slrts.SetAltParameters(*sb_model.GetSnapshot());
362 HybridCalculator hc2(*data, sb_model, b_model);
363 ToyMCSampler *toymcs2 = (ToyMCSampler *)hc2.GetTestStatSampler();
365 toymcs2->SetTestStatistic(&slrts);
367 hc2.SetToys(20000, 1000);
368 hc2.ForcePriorNuisanceAlt(*w->pdf(
"py"));
369 hc2.ForcePriorNuisanceNull(*w->pdf(
"py"));
383 toymcs2->SetProofConfig(pc);
386 RooMsgService::instance().setGlobalKillBelow(RooFit::ERROR);
388 HypoTestResult *r2 = hc2.GetHypoTest();
389 cout <<
"-----------------------------------------" << endl;
390 cout <<
"Part 5" << endl;
396 RooMsgService::instance().setGlobalKillBelow(msglevel);
399 HypoTestPlot *p2 =
new HypoTestPlot(*r2, 30);