-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathStandardHypoTestInvDemo.C
356 lines (277 loc) · 12.1 KB
/
StandardHypoTestInvDemo.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
// Standard tutorial macro for performing an inverted hypothesis test
//
// This macro will perform a scan of tehe p-values for computing the limit
//
#include "TFile.h"
#include "RooWorkspace.h"
#include "RooAbsPdf.h"
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooStats/ModelConfig.h"
#include "TGraphErrors.h"
#include "TGraphAsymmErrors.h"
#include "TCanvas.h"
#include "TLine.h"
#include "RooStats/HybridCalculator.h"
#include "RooStats/FrequentistCalculator.h"
#include "RooStats/ToyMCSampler.h"
#include "RooStats/HypoTestPlot.h"
#include "RooStats/NumEventsTestStat.h"
//#include "MyProfileLikelihoodTestStat.h"
#include "RooStats/ProfileLikelihoodTestStat.h"
#include "RooStats/SimpleLikelihoodRatioTestStat.h"
#include "RooStats/RatioOfProfiledLikelihoodsTestStat.h"
#include "RooStats/MaxLikelihoodEstimateTestStat.h"
#include "RooStats/HypoTestInverter.h"
#include "RooStats/HypoTestInverterResult.h"
#include "RooStats/HypoTestInverterPlot.h"
using namespace RooFit;
using namespace RooStats;
bool plotHypoTestResult = true;
bool useProof = true;
bool optimize = false;
bool writeResult = false;
int nworkers = 8;
// internal routine to run the inverter
HypoTestInverterResult * RunInverter(RooWorkspace * w, const char * modelSBName, const char * modelBName, const char * dataName,
int type, int testStatType, int npoints, double poimin, double poimax, int ntoys, bool useCls );
void StandardHypoTestInvDemo(const char * fileName =0,
const char * wsName = "combined",
const char * modelSBName = "ModelConfig",
const char * modelBName = "",
const char * dataName = "obsData",
int calculatorType = 0,
int testStatType = 3,
bool useCls = true ,
int npoints = 5,
double poimin = 0,
double poimax = 5,
int ntoys=1000 )
{
/*
Other Parameter to pass in tutorial
apart from standard for filename, ws, modelconfig and data
type = 0 Freq calculator
type = 1 Hybrid
testStatType = 0 LEP
= 1 Tevatron
= 2 Profile Likelihood
= 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
useCLs scan for CLs (otherwise for CLs+b)
npoints: number of points to scan , for autoscan set npoints = -1
poimin,poimax: min/max value to scan in case of fixed scans
(if min >= max, try to find automatically)
ntoys: number of toys to use
extra options are available as global paramters of the macro. They are:
plotHypoTestResult plot result of tests at each point (TS distributions)
useProof = true;
writeResult = true;
nworkers = 4;
*/
if (fileName==0) {
fileName = "results/example_combined_GaussExample_model.root";
std::cout << "Use standard file generated with HistFactory :" << fileName << std::endl;
}
TFile * file = new TFile(fileName);
RooWorkspace * w = dynamic_cast<RooWorkspace*>( file->Get(wsName) );
HypoTestInverterResult * r = 0;
std::cout << w << "\t" << fileName << std::endl;
if (w != NULL) {
r = RunInverter(w, modelSBName, modelBName, dataName, calculatorType, testStatType, npoints, poimin, poimax, ntoys, useCls );
if (!r) {
std::cerr << "Error running the HypoTestInverter - Exit " << std::endl;
return;
}
}
else
{
// case workspace is not present look for the inverter result
std::cout << "Reading an HypoTestInverterResult with name " << wsName << " from file " << fileName << std::endl;
r = dynamic_cast<HypoTestInverterResult*>( file->Get(wsName) ); //
if (!r) {
std::cerr << "File " << fileName << " does not contain a workspace or an HypoTestInverterResult - Exit "
<< std::endl;
file->ls();
return;
}
}
double upperLimit = r->UpperLimit();
double ulError = r->UpperLimitEstimatedError();
std::cout << "The computed upper limit is: " << upperLimit << " +/- " << ulError << std::endl;
const int nEntries = r->ArraySize();
const char * typeName = (calculatorType == 0) ? "Frequentist" : "Hybrid";
const char * resultName = (w) ? w->GetName() : r->GetName();
TString plotTitle = TString::Format("%s CL Scan for workspace %s",typeName,resultName);
HypoTestInverterPlot *plot = new HypoTestInverterPlot("HTI_Result_Plot",plotTitle,r);
plot->Draw("CLb 2CL"); // plot all and Clb
if (plotHypoTestResult) {
TCanvas * c2 = new TCanvas();
c2->Divide( 2, TMath::Ceil(nEntries/2));
for (int i=0; i<nEntries; i++) {
c2->cd(i+1);
SamplingDistPlot * pl = plot->MakeTestStatPlot(i);
pl->SetLogYaxis(true);
pl->Draw();
}
}
std::cout << " expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
std::cout << " expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
std::cout << " expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
if (w != NULL && writeResult) {
// write to a file the results
const char * calcType = (calculatorType == 0) ? "Freq" : "Hybr";
const char * limitType = (useCls) ? "CLs" : "Cls+b";
const char * scanType = (npoints < 0) ? "auto" : "grid";
TString resultFileName = TString::Format("%s_%s_%s_ts%d_",calcType,limitType,scanType,testStatType);
resultFileName += fileName;
TFile * fileOut = new TFile(resultFileName,"RECREATE");
r->Write();
fileOut->Close();
}
}
// internal routine to run the inverter
HypoTestInverterResult * RunInverter(RooWorkspace * w, const char * modelSBName, const char * modelBName,
const char * dataName, int type, int testStatType,
int npoints, double poimin, double poimax,
int ntoys, bool useCls )
{
std::cout << "Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
w->Print();
RooAbsData * data = w->data(dataName);
if (!data) {
Error("StandardHypoTestDemo","Not existing data %s",dataName);
return 0;
}
else
std::cout << "Using data set " << dataName << std::endl;
// get models from WS
// get the modelConfig out of the file
ModelConfig* bModel = (ModelConfig*) w->obj(modelBName);
ModelConfig* sbModel = (ModelConfig*) w->obj(modelSBName);
if (!sbModel) {
Error("StandardHypoTestDemo","Not existing ModelConfig %s",modelSBName);
return 0;
}
// check the model
if (!sbModel->GetPdf()) {
Error("StandardHypoTestDemo","Model %s has no pdf ",modelSBName);
return 0;
}
if (!sbModel->GetParametersOfInterest()) {
Error("StandardHypoTestDemo","Model %s has no poi ",modelSBName);
return 0;
}
if (!sbModel->GetParametersOfInterest()) {
Error("StandardHypoTestInvDemo","Model %s has no poi ",modelSBName);
return 0;
}
if (!sbModel->GetSnapshot() ) {
Info("StandardHypoTestInvDemo","Model %s has no snapshot - make one using model poi",modelSBName);
sbModel->SetSnapshot( *sbModel->GetParametersOfInterest() );
}
if (!bModel || bModel == sbModel) {
Info("StandardHypoTestInvDemo","The background model %s does not exist",modelBName);
Info("StandardHypoTestInvDemo","Copy it from ModelConfig %s and set POI to zero",modelSBName);
bModel = (ModelConfig*) sbModel->Clone();
bModel->SetName(TString(modelSBName)+TString("_with_poi_0"));
RooRealVar * var = dynamic_cast<RooRealVar*>(bModel->GetParametersOfInterest()->first());
if (!var) return 0;
double oldval = var->getVal();
var->setVal(0);
bModel->SetSnapshot( RooArgSet(*var) );
var->setVal(oldval);
}
else {
if (!bModel->GetSnapshot() ) {
Info("StandardHypoTestInvDemo","Model %s has no snapshot - make one using model poi and 0 values ",modelBName);
RooRealVar * var = dynamic_cast<RooRealVar*>(bModel->GetParametersOfInterest()->first());
if (var) {
double oldval = var->getVal();
var->setVal(0);
bModel->SetSnapshot( RooArgSet(*var) );
var->setVal(oldval);
}
else {
Error("StandardHypoTestInvDemo","Model %s has no valid poi",modelBName);
return 0;
}
}
}
SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(),*bModel->GetPdf());
if (sbModel->GetSnapshot()) slrts.SetNullParameters(*sbModel->GetSnapshot());
if (bModel->GetSnapshot()) slrts.SetAltParameters(*bModel->GetSnapshot());
// ratio of profile likelihood - need to pass snapshot for the alt
RatioOfProfiledLikelihoodsTestStat
ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
ropl.SetSubtractMLE(false);
//MyProfileLikelihoodTestStat profll(*sbModel->GetPdf());
ProfileLikelihoodTestStat profll(*sbModel->GetPdf());
if (testStatType == 3) profll.SetOneSided(1);
if (optimize) profll.SetReuseNLL(true);
TestStatistic * testStat = &slrts;
if (testStatType == 1) testStat = &ropl;
if (testStatType == 2 || testStatType == 3) testStat = &profll;
HypoTestCalculatorGeneric * hc = 0;
if (type == 0) hc = new FrequentistCalculator(*data, *bModel, *sbModel);
else hc = new HybridCalculator(*data, *bModel, *sbModel);
ToyMCSampler *toymcs = (ToyMCSampler*)hc->GetTestStatSampler();
toymcs->SetNEventsPerToy(1);
toymcs->SetTestStatistic(testStat);
if (optimize) toymcs->SetUseMultiGen(true);
if (type == 1) {
HybridCalculator *hhc = (HybridCalculator*) hc;
hhc->SetToys(ntoys,ntoys);
// check for nuisance prior pdf
if (bModel->GetPriorPdf() && sbModel->GetPriorPdf() ) {
hhc->ForcePriorNuisanceAlt(*bModel->GetPriorPdf());
hhc->ForcePriorNuisanceNull(*sbModel->GetPriorPdf());
}
else {
if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters() ) {
Error("StandardHypoTestInvDemo","Cannnot run Hybrid calculator because no prior on the nuisance parameter is specified");
return 0;
}
}
}
else
((FrequentistCalculator*) hc)->SetToys(ntoys,ntoys);
// Get the result
RooMsgService::instance().getStream(1).removeTopic(RooFit::NumIntegration);
TStopwatch tw; tw.Start();
const RooArgSet * poiSet = sbModel->GetParametersOfInterest();
RooRealVar *poi = (RooRealVar*)poiSet->first();
// fit the data first
sbModel->GetPdf()->fitTo(*data);
double poihat = poi->getVal();
HypoTestInverter calc(*hc);
calc.SetConfidenceLevel(0.95);
calc.UseCLs(useCls);
calc.SetVerbose(true);
// can speed up using proof-lite
if (useProof && nworkers > 1) {
ProofConfig pc(*w, nworkers, "", kFALSE);
toymcs->SetProofConfig(&pc); // enable proof
}
if (npoints > 0) {
if (poimin >= poimax) {
// if no min/max given scan between MLE and +4 sigma
poimin = int(poihat);
poimax = int(poihat + 4 * poi->getError());
}
std::cout << "Doing a fixed scan in interval : " << poimin << " , " << poimax << std::endl;
calc.SetFixedScan(npoints,poimin,poimax);
}
else {
//poi->setMax(10*int( (poihat+ 10 *poi->getError() )/10 ) );
std::cout << "Doing an automatic scan in interval : " << poi->getMin() << " , " << poi->getMax() << std::endl;
}
HypoTestInverterResult * r = calc.GetInterval();
return r;
}
void ReadResult(const char * fileName, const char * resultName="", bool useCLs=true) {
// read a previous stored result from a file given the result name
StandardHypoTestInvDemo(fileName, resultName,"","","",0,0,useCLs);
}
int main() {
StandardHypoTestInvDemo();
}