52 using namespace RooFit;
53 using namespace RooStats;
55 double StandardFrequentistDiscovery(
const char *infile =
"",
const char *workspaceName =
"channel1",
56 const char *modelConfigNameSB =
"ModelConfig",
const char *dataName =
"obsData",
57 int toys = 1000,
double poiValueForBackground = 0.0,
double poiValueForSignal = 1.0)
71 const char *filename =
"";
72 if (!strcmp(infile,
"")) {
73 filename =
"results/example_channel1_GammaExample_model.root";
74 bool fileExist = !gSystem->AccessPathName(filename);
78 cout <<
"HistFactory file cannot be generated on Windows - exit" << endl;
82 cout <<
"will run standard hist2workspace example" << endl;
83 gROOT->ProcessLine(
".! prepareHistFactory .");
84 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
85 cout <<
"\n\n---------------------" << endl;
86 cout <<
"Done creating example input" << endl;
87 cout <<
"---------------------\n\n" << endl;
94 TFile *file = TFile::Open(filename);
98 cout <<
"StandardRooStatsDemoMacro: Input file " << filename <<
" is not found" << endl;
106 TStopwatch *mn_t =
new TStopwatch;
110 RooWorkspace *w = (RooWorkspace *)file->Get(workspaceName);
112 cout <<
"workspace not found" << endl;
117 ModelConfig *mc = (ModelConfig *)w->obj(modelConfigNameSB);
120 RooAbsData *data = w->data(dataName);
125 cout <<
"data or ModelConfig was not found" << endl;
129 RooRealVar *firstPOI = (RooRealVar *)mc->GetParametersOfInterest()->first();
130 firstPOI->setVal(poiValueForSignal);
131 mc->SetSnapshot(*mc->GetParametersOfInterest());
133 ModelConfig *mcNull = mc->Clone(
"ModelConfigNull");
134 firstPOI->setVal(poiValueForBackground);
135 mcNull->SetSnapshot(*(RooArgSet *)mcNull->GetParametersOfInterest()->snapshot());
140 ProfileLikelihoodTestStat *plts =
new ProfileLikelihoodTestStat(*mc->GetPdf());
141 plts->SetOneSidedDiscovery(
true);
142 plts->SetVarName(
"q_{0}/2");
146 ToyMCSampler toymcs(*plts, 50);
155 if (!mc->GetPdf()->canBeExtended()) {
156 if (data->numEntries() == 1) {
157 toymcs.SetNEventsPerToy(1);
159 cout <<
"Not sure what to do about this model" << endl;
164 ProofConfig pc(*w, 2,
"",
false);
168 FrequentistCalculator freqCalc(*data, *mc, *mcNull, &toymcs);
169 freqCalc.SetToys(toys, toys);
172 HypoTestResult *freqCalcResult = freqCalc.GetHypoTest();
173 freqCalcResult->GetNullDistribution()->SetTitle(
"b only");
174 freqCalcResult->GetAltDistribution()->SetTitle(
"s+b");
175 freqCalcResult->Print();
176 double pvalue = freqCalcResult->NullPValue();
180 cout <<
"total CPU time: " << mn_t->CpuTime() << endl;
181 cout <<
"total real time: " << mn_t->RealTime() << endl;
184 TCanvas *c1 =
new TCanvas();
185 HypoTestPlot *plot =
new HypoTestPlot(*freqCalcResult, 100, -0.49, 9.51);
186 plot->SetLogYaxis(
true);
190 TF1 *f =
new TF1(
"f", TString::Format(
"1*ROOT::Math::chisquared_pdf(2*x,%d,0)", nPOI), 0, 20);
191 f->SetLineColor(kBlack);
193 plot->AddTF1(f, TString::Format(
"#chi^{2}(2x,%d)", nPOI));
196 c1->SaveAs(
"standard_discovery_output.pdf");