Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
StandardHypoTestInvDemo.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook
4/// Standard tutorial macro for performing an inverted hypothesis test for computing an interval
5///
6/// This macro will perform a scan of the p-values for computing the interval or limit
7///
8/// Usage:
9///
10/// ~~~{.cpp}
11/// root>.L StandardHypoTestInvDemo.C
12/// root> StandardHypoTestInvDemo("fileName","workspace name","S+B modelconfig name","B model name","data set
13/// name",calculator type, test statistic type, use CLS,
14/// number of points, xmin, xmax, number of toys, use number counting)
15///
16/// type = 0 Freq calculator
17/// type = 1 Hybrid calculator
18/// type = 2 Asymptotic calculator
19/// type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
20///
21/// testStatType = 0 LEP
22/// = 1 Tevatron
23/// = 2 Profile Likelihood two sided
24/// = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
25/// = 4 Profile Likelihood signed ( pll = -pll if mu < mu_hat)
26/// = 5 Max Likelihood Estimate as test statistic
27/// = 6 Number of observed event as test statistic
28/// ~~~
29///
30/// \macro_image
31/// \macro_output
32/// \macro_code
33///
34/// \author Lorenzo Moneta
35
36#include "TFile.h"
37#include "RooWorkspace.h"
38#include "RooAbsPdf.h"
39#include "RooRealVar.h"
40#include "RooDataSet.h"
42#include "RooRandom.h"
43#include "TGraphErrors.h"
44#include "TGraphAsymmErrors.h"
45#include "TCanvas.h"
46#include "TLine.h"
47#include "TROOT.h"
48#include "TSystem.h"
49
55
62
66
67#include <cassert>
68
69using namespace RooFit;
70using namespace RooStats;
71using std::cout, std::endl;
72
73// structure defining the options
74struct HypoTestInvOptions {
75
76 bool plotHypoTestResult = true; // plot test statistic result at each point
77 bool writeResult = true; // write HypoTestInverterResult in a file
78 TString resultFileName; // file with results (by default is built automatically using the workspace input file name)
79 bool optimize = true; // optimize evaluation of test statistic
80 bool useVectorStore = true; // convert data to use new roofit data store
81 bool generateBinned = false; // generate binned data sets
82 bool noSystematics = false; // force all systematics to be off (i.e. set all nuisance parameters as constat
83 // to their nominal values)
84 double nToysRatio = 2; // ratio Ntoys S+b/ntoysB
85 double maxPOI = -1; // max value used of POI (in case of auto scan)
87 false; // enable detailed output with all fit information for each toys (output will be written in result file)
88 bool rebuild = false; // re-do extra toys for computing expected limits and rebuild test stat
89 // distributions (N.B this requires much more CPU (factor is equivalent to nToyToRebuild)
90 int nToyToRebuild = 100; // number of toys used to rebuild
91 int rebuildParamValues = 0; // = 0 do a profile of all the parameters on the B (alt snapshot) before performing a
92 // rebuild operation (default)
93 // = 1 use initial workspace parameters with B snapshot values
94 // = 2 use all initial workspace parameters with B
95 // Otherwise the rebuild will be performed using
96 int initialFit = -1; // do a first fit to the model (-1 : default, 0 skip fit, 1 do always fit)
97 int randomSeed = -1; // random seed (if = -1: use default value, if = 0 always random )
98
99 int nAsimovBins = 0; // number of bins in observables used for Asimov data sets (0 is the default and it is given by
100 // workspace, typically is 100)
101
102 bool reuseAltToys = false; // reuse same toys for alternate hypothesis (if set one gets more stable bands)
103 double confLevel = 0.95; // confidence level value
104
105 std::string minimizerType =
106 ""; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
107 std::string massValue = ""; // extra string to tag output file of result
108 int printLevel = 0; // print level for debugging PL test statistics and calculators
109
110 bool useNLLOffset = false; // use NLL offset when fitting (this increase stability of fits)
111};
112
114
115// internal class to run the inverter and more
116
117namespace RooStats {
118
119class HypoTestInvTool {
120
121public:
124
126 const char *dataName, int type, int testStatType, bool useCLs, int npoints,
127 double poimin, double poimax, int ntoys, bool useNumberCounting = false,
128 const char *nuisPriorName = nullptr);
129
131 const char *fileNameBase = nullptr);
132
133 void SetParameter(const char *name, const char *value);
134 void SetParameter(const char *name, bool value);
135 void SetParameter(const char *name, int value);
136 void SetParameter(const char *name, double value);
137
138private:
140 bool mWriteResult;
141 bool mOptimize;
142 bool mUseVectorStore;
143 bool mGenerateBinned;
144 bool mRebuild;
145 bool mReuseAltToys;
146 bool mEnableDetOutput;
147 int mNToyToRebuild;
149 int mPrintLevel;
150 int mInitialFit;
151 int mRandomSeed;
152 double mNToysRatio;
153 double mMaxPoi;
154 int mAsimovBins;
155 std::string mMassValue;
156 std::string
157 mMinimizerType; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
159};
160
161} // end namespace RooStats
162
163RooStats::HypoTestInvTool::HypoTestInvTool()
168{
169}
170
171void RooStats::HypoTestInvTool::SetParameter(const char *name, bool value)
172{
173 //
174 // set boolean parameters
175 //
176
177 std::string s_name(name);
178
179 if (s_name.find("PlotHypoTestResult") != std::string::npos)
181 if (s_name.find("WriteResult") != std::string::npos)
183 if (s_name.find("Optimize") != std::string::npos)
185 if (s_name.find("UseVectorStore") != std::string::npos)
187 if (s_name.find("GenerateBinned") != std::string::npos)
189 if (s_name.find("EnableDetailedOutput") != std::string::npos)
191 if (s_name.find("Rebuild") != std::string::npos)
192 mRebuild = value;
193 if (s_name.find("ReuseAltToys") != std::string::npos)
195
196 return;
197}
198
199void RooStats::HypoTestInvTool::SetParameter(const char *name, int value)
200{
201 //
202 // set integer parameters
203 //
204
205 std::string s_name(name);
206
207 if (s_name.find("NToyToRebuild") != std::string::npos)
209 if (s_name.find("RebuildParamValues") != std::string::npos)
211 if (s_name.find("PrintLevel") != std::string::npos)
213 if (s_name.find("InitialFit") != std::string::npos)
215 if (s_name.find("RandomSeed") != std::string::npos)
217 if (s_name.find("AsimovBins") != std::string::npos)
219
220 return;
221}
222
223void RooStats::HypoTestInvTool::SetParameter(const char *name, double value)
224{
225 //
226 // set double precision parameters
227 //
228
229 std::string s_name(name);
230
231 if (s_name.find("NToysRatio") != std::string::npos)
233 if (s_name.find("MaxPOI") != std::string::npos)
234 mMaxPoi = value;
235
236 return;
237}
238
239void RooStats::HypoTestInvTool::SetParameter(const char *name, const char *value)
240{
241 //
242 // set string parameters
243 //
244
245 std::string s_name(name);
246
247 if (s_name.find("MassValue") != std::string::npos)
248 mMassValue.assign(value);
249 if (s_name.find("MinimizerType") != std::string::npos)
250 mMinimizerType.assign(value);
251 if (s_name.find("ResultFileName") != std::string::npos)
253
254 return;
255}
256
257void wsName = "combined",
258 const char *modelSBName = "ModelConfig", const char *modelBName = "",
259 const char *dataName = "obsData", int calculatorType = 0, int testStatType = 0,
260 bool useCLs = true, int npoints = 6, double poimin = 0, double poimax = 5,
261 int ntoys = 1000, bool useNumberCounting = false, const char *nuisPriorName = nullptr)
262{
263 /*
264
265 Other Parameter to pass in tutorial
266 apart from standard for filename, ws, modelconfig and data
267
268 type = 0 Freq calculator
269 type = 1 Hybrid calculator
270 type = 2 Asymptotic calculator
271 type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
272
273 testStatType = 0 LEP
274 = 1 Tevatron
275 = 2 Profile Likelihood
276 = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
277 = 4 Profiel Likelihood signed ( pll = -pll if mu < mu_hat)
278 = 5 Max Likelihood Estimate as test statistic
279 = 6 Number of observed event as test statistic
280
281 useCLs scan for CLs (otherwise for CLs+b)
282
283 npoints: number of points to scan , for autoscan set npoints = -1
284
285 poimin,poimax: min/max value to scan in case of fixed scans
286 (if min > max, try to find automatically)
287
288 ntoys: number of toys to use
289
290 useNumberCounting: set to true when using number counting events
291
292 nuisPriorName: name of prior for the nuisance. This is often expressed as constraint term in the global model
293 It is needed only when using the HybridCalculator (type=1)
294 If not given by default the prior pdf from ModelConfig is used.
295
296 extra options are available as global parameters of the macro. They major ones are:
297
298 plotHypoTestResult plot result of tests at each point (TS distributions) (default is true)
299 writeResult write result of scan (default is true)
300 rebuild rebuild scan for expected limits (require extra toys) (default is false)
301 generateBinned generate binned data sets for toys (default is false) - be careful not to activate with
302 a too large (>=3) number of observables
303 nToyRatio ratio of S+B/B toys (default is 2)
304
305
306 */
307
309 if (filename.IsNull()) {
310 filename = "results/example_combined_GaussExample_model.root";
311 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
312 // if file does not exists generate with histfactory
313 if (!fileExist) {
314 // Normally this would be run on the command line
315 cout << "will run standard hist2workspace example" << endl;
316 gROOT->ProcessLine(".! prepareHistFactory .");
317 gROOT->ProcessLine(".! hist2workspace config/example.xml");
318 cout << "\n\n---------------------" << endl;
319 cout << "Done creating example input" << endl;
320 cout << "---------------------\n\n" << endl;
321 }
322
323 } else
325
326 // Try to open the file
327 TFile *file = TFile::Open(filename);
328
329 // if input file was specified but not found, quit
330 if (!file) {
331 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
332 return;
333 }
334
336
337 // set parameters
338 calc.SetParameter("PlotHypoTestResult", optHTInv.plotHypoTestResult);
339 calc.SetParameter("WriteResult", optHTInv.writeResult);
340 calc.SetParameter("Optimize", optHTInv.optimize);
341 calc.SetParameter("UseVectorStore", optHTInv.useVectorStore);
342 calc.SetParameter("GenerateBinned", optHTInv.generateBinned);
343 calc.SetParameter("NToysRatio", optHTInv.nToysRatio);
344 calc.SetParameter("MaxPOI", optHTInv.maxPOI);
345 calc.SetParameter("EnableDetailedOutput", optHTInv.enableDetailedOutput);
346 calc.SetParameter("Rebuild", optHTInv.rebuild);
347 calc.SetParameter("ReuseAltToys", optHTInv.reuseAltToys);
348 calc.SetParameter("NToyToRebuild", optHTInv.nToyToRebuild);
349 calc.SetParameter("RebuildParamValues", optHTInv.rebuildParamValues);
350 calc.SetParameter("MassValue", optHTInv.massValue.c_str());
351 calc.SetParameter("MinimizerType", optHTInv.minimizerType.c_str());
352 calc.SetParameter("PrintLevel", optHTInv.printLevel);
353 calc.SetParameter("InitialFit", optHTInv.initialFit);
354 calc.SetParameter("ResultFileName", optHTInv.resultFileName);
355 calc.SetParameter("RandomSeed", optHTInv.randomSeed);
356 calc.SetParameter("AsimovBins", optHTInv.nAsimovBins);
357
358 // enable offset for all roostats
359 if (optHTInv.useNLLOffset)
361
362 RooWorkspace *w = dynamic_cast<RooWorkspace *>(file->Get(wsName));
363 HypoTestInverterResult *r = nullptr;
364 std::cout << w << "\t" << filename << std::endl;
365 if (w != nullptr) {
368 if (!r) {
369 std::cerr << "Error running the HypoTestInverter - Exit " << std::endl;
370 return;
371 }
372 } else {
373 // case workspace is not present look for the inverter result
374 std::cout << "Reading an HypoTestInverterResult with name " << wsName << " from file " << filename << std::endl;
375 r = dynamic_cast<HypoTestInverterResult *>(file->Get(wsName)); //
376 if (!r) {
377 std::cerr << "File " << filename << " does not contain a workspace or an HypoTestInverterResult - Exit "
378 << std::endl;
379 file->ls();
380 return;
381 }
382 }
383
385
386 return;
387}
388
389void RooStats::HypoTestInvTool::AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType,
390 bool useCLs, int npoints, const char *fileNameBase)
391{
392
393 // analyze result produced by the inverter, optionally save it in a file
394
395 double lowerLimit = 0;
396 double llError = 0;
397#if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
398 if (r->IsTwoSided()) {
399 lowerLimit = r->LowerLimit();
400 llError = r->LowerLimitEstimatedError();
401 }
402#else
403 lowerLimit = r->LowerLimit();
404 llError = r->LowerLimitEstimatedError();
405#endif
406
407 double upperLimit = r->UpperLimit();
408 double ulError = r->UpperLimitEstimatedError();
409
410 // std::cout << "DEBUG : [ " << lowerLimit << " , " << upperLimit << " ] " << std::endl;
411
412 if (lowerLimit < upperLimit * (1. - 1.E-4) && lowerLimit != 0)
413 std::cout << "The computed lower limit is: " << lowerLimit << " +/- " << llError << std::endl;
414 std::cout << "The computed upper limit is: " << upperLimit << " +/- " << ulError << std::endl;
415
416 // compute expected limit
417 std::cout << "Expected upper limits, using the B (alternate) model : " << std::endl;
418 std::cout << " expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
419 std::cout << " expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
420 std::cout << " expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
421 std::cout << " expected limit (-2 sig) " << r->GetExpectedUpperLimit(-2) << std::endl;
422 std::cout << " expected limit (+2 sig) " << r->GetExpectedUpperLimit(2) << std::endl;
423
424 // detailed output
425 if (mEnableDetOutput) {
426 mWriteResult = true;
427 Info("StandardHypoTestInvDemo", "detailed output will be written in output result file");
428 }
429
430 // write result in a file
431 if (r != nullptr && mWriteResult) {
432
433 // write to a file the results
434 const char *calcType = (calculatorType == 0) ? "Freq" : (calculatorType == 1) ? "Hybr" : "Asym";
435 const char *limitType = (useCLs) ? "CLs" : "Cls+b";
436 const char *scanType = (npoints < 0) ? "auto" : "grid";
437 if (mResultFileName.IsNull()) {
439 // strip the / from the filename
440 if (!mMassValue.empty()) {
441 mResultFileName += mMassValue.c_str();
442 mResultFileName += "_";
443 }
444
446 name.Replace(0, name.Last('/') + 1, "");
448 }
449
450 // get (if existing) rebuilt UL distribution
451 TString uldistFile = "RULDist.root";
452 TObject *ulDist = nullptr;
454 if (existULDist) {
456 if (fileULDist)
457 ulDist = fileULDist->Get("RULDist");
458 }
459
460 mResultFileName, "RECREATE");
461 r->Write();
462 if (ulDist)
463 ulDist->Write();
464 Info("StandardHypoTestInvDemo", "HypoTestInverterResult has been written in the file %s", mResultFileName.Data());
465
466 fileOut->Close();
467 }
468
469 // plot the result ( p values vs scan points)
470 std::string typeName = "";
471 if (calculatorType == 0)
472 typeName = "Frequentist";
473 if (calculatorType == 1)
474 typeName = "Hybrid";
475 else if (calculatorType == 2 || calculatorType == 3) {
476 typeName = "Asymptotic";
477 mPlotHypoTestResult = false;
478 }
479
480 const char *resultName = r->GetName();
481 TString plotTitle = TString::Format("%s CL Scan for workspace %s", typeName.c_str(), resultName);
482 HypoTestInverterPlot *plot = new HypoTestInverterPlot("HTI_Result_Plot", plotTitle, r);
483
484 // plot in a new canvas with style
485 TString c1Name = TString::Format("%s_Scan", typeName.c_str());
486 TCanvas *c1 = new TCanvas(c1Name);
487 c1->SetLogy(false);
488
489 plot->Draw("CLb 2CL"); // plot all and Clb
490
491 // if (useCLs)
492 // plot->Draw("CLb 2CL"); // plot all and Clb
493 // else
494 // plot->Draw(""); // plot all and Clb
495
496 const int nEntries = r->ArraySize();
497
498 // plot test statistics distributions for the two hypothesis
500 TCanvas *c2 = new TCanvas("c2");
501 if (nEntries > 1) {
503 int nx = TMath::CeilNint(double(nEntries) / ny);
504 c2->Divide(nx, ny);
505 }
506 for (int i = 0; i < nEntries; i++) {
507 if (nEntries > 1)
508 c2->cd(i + 1);
509 SamplingDistPlot *pl = plot->MakeTestStatPlot(i);
510 pl->SetLogYaxis(true);
511 pl->Draw();
512 }
513 }
514 gPad = c1;
515}
516
517// internal routine to run the inverter
518HypoTestInverterResult *RooStats::HypoTestInvTool::RunInverter(RooWorkspace *w, const char *modelSBName,
519 const char *modelBName, const char *dataName, int type,
520 int testStatType, bool useCLs, int npoints,
521 double poimin, double poimax, int ntoys,
522 bool useNumberCounting, const char *nuisPriorName)
523{
524
525 std::cout << "Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
526
527 w->Print();
528
529 RooAbsData *data = w->data(dataName);
530 if (!data) {
531 Error("StandardHypoTestDemo", "Not existing data %s", dataName);
532 return nullptr;
533 } else
534 std::cout << "Using data set " << dataName << std::endl;
535
536 if (mUseVectorStore) {
538 data->convertToVectorStore();
539 }
540
541 // get models from WS
542 // get the modelConfig out of the file
545
546 if (!sbModel) {
547 Error("StandardHypoTestDemo", "Not existing ModelConfig %s", modelSBName);
548 return nullptr;
549 }
550 // check the model
551 if (!sbModel->GetPdf()) {
552 Error("StandardHypoTestDemo", "Model %s has no pdf ", modelSBName);
553 return nullptr;
554 }
555 if (!sbModel->GetParametersOfInterest()) {
556 Error("StandardHypoTestDemo", "Model %s has no poi ", modelSBName);
557 return nullptr;
558 }
559 if (!sbModel->GetObservables()) {
560 Error("StandardHypoTestInvDemo", "Model %s has no observables ", modelSBName);
561 return nullptr;
562 }
563 if (!sbModel->GetSnapshot()) {
564 Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi", modelSBName);
565 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
566 }
567
568 // case of no systematics
569 // remove nuisance parameters from model
570 if (optHTInv.noSystematics) {
571 const RooArgSet *nuisPar = sbModel->GetNuisanceParameters();
572 if (nuisPar && nuisPar->getSize() > 0) {
573 std::cout << "StandardHypoTestInvDemo"
574 << " - Switch off all systematics by setting them constant to their initial values" << std::endl;
576 }
577 if (bModel) {
578 const RooArgSet *bnuisPar = bModel->GetNuisanceParameters();
579 if (bnuisPar)
581 }
582 }
583
584 if (!bModel || bModel == sbModel) {
585 Info("StandardHypoTestInvDemo", "The background model %s does not exist", modelBName);
586 Info("StandardHypoTestInvDemo", "Copy it from ModelConfig %s and set POI to zero", modelSBName);
587 bModel = (ModelConfig *)sbModel->Clone();
588 bModel->SetName(TString(modelSBName) + TString("_with_poi_0"));
589 RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
590 if (!var)
591 return nullptr;
592 double oldval = var->getVal();
593 var->setVal(0);
594 bModel->SetSnapshot(RooArgSet(*var));
595 var->setVal(oldval);
596 } else {
597 if (!bModel->GetSnapshot()) {
598 Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi and 0 values ",
599 modelBName);
600 RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
601 if (var) {
602 double oldval = var->getVal();
603 var->setVal(0);
604 bModel->SetSnapshot(RooArgSet(*var));
605 var->setVal(oldval);
606 } else {
607 Error("StandardHypoTestInvDemo", "Model %s has no valid poi", modelBName);
608 return nullptr;
609 }
610 }
611 }
612
613 // check model has global observables when there are nuisance pdf
614 // for the hybrid case the globals are not needed
615 if (type != 1) {
616 bool hasNuisParam = (sbModel->GetNuisanceParameters() && sbModel->GetNuisanceParameters()->getSize() > 0);
617 bool hasGlobalObs = (sbModel->GetGlobalObservables() && sbModel->GetGlobalObservables()->getSize() > 0);
618 if (hasNuisParam && !hasGlobalObs) {
619 // try to see if model has nuisance parameters first
620 RooAbsPdf *constrPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisanceConstraintPdf_sbmodel");
621 if (constrPdf) {
622 Warning("StandardHypoTestInvDemo", "Model %s has nuisance parameters but no global observables associated",
623 sbModel->GetName());
624 Warning("StandardHypoTestInvDemo",
625 "\tThe effect of the nuisance parameters will not be treated correctly ");
626 }
627 }
628 }
629
630 // save all initial parameters of the model including the global observables
632 std::unique_ptr<RooArgSet> allParams{sbModel->GetPdf()->getParameters(*data)};
633 allParams->snapshot(initialParameters);
634
635 // run first a data fit
636
637 const RooArgSet *poiSet = sbModel->GetParametersOfInterest();
638 RooRealVar *poi = (RooRealVar *)poiSet->first();
639
640 std::cout << "StandardHypoTestInvDemo : POI initial value: " << poi->GetName() << " = " << poi->getVal()
641 << std::endl;
642
643 // fit the data first (need to use constraint )
645
646 bool doFit = mInitialFit;
647 if (testStatType == 0 && mInitialFit == -1)
648 doFit = false; // case of LEP test statistic
649 if (type == 3 && mInitialFit == -1)
650 doFit = false; // case of Asymptoticcalculator with nominal Asimov
651 double poihat = 0;
652
653 if (mMinimizerType.empty())
655 else
657
658 Info("StandardHypoTestInvDemo", "Using %s as minimizer for computing the test statistic",
660
661 if (doFit) {
662
663 // do the fit : By doing a fit the POI snapshot (for S+B) is set to the fit value
664 // and the nuisance parameters nominal values will be set to the fit value.
665 // This is relevant when using LEP test statistics
666
667 Info("StandardHypoTestInvDemo", " Doing a first fit to the observed data ");
669 if (sbModel->GetNuisanceParameters())
670 constrainParams.add(*sbModel->GetNuisanceParameters());
672 tw.Start();
673 std::unique_ptr<RooFitResult> fitres{sbModel->GetPdf()->fitTo(
674 *data, InitialHesse(false), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0),
676 if (fitres->status() != 0) {
677 Warning("StandardHypoTestInvDemo",
678 "Fit to the model failed - try with strategy 1 and perform first an Hesse computation");
679 fitres = std::unique_ptr<RooFitResult>{sbModel->GetPdf()->fitTo(
680 *data, InitialHesse(true), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(1),
682 }
683 if (fitres->status() != 0)
684 Warning("StandardHypoTestInvDemo", " Fit still failed - continue anyway.....");
685
686 poihat = poi->getVal();
687 std::cout << "StandardHypoTestInvDemo - Best Fit value : " << poi->GetName() << " = " << poihat << " +/- "
688 << poi->getError() << std::endl;
689 std::cout << "Time for fitting : ";
690 tw.Print();
691
692 // save best fit value in the poi snapshot
693 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
694 std::cout << "StandardHypoTestInvo: snapshot of S+B Model " << sbModel->GetName()
695 << " is set to the best fit value" << std::endl;
696 }
697
698 // print a message in case of LEP test statistics because it affects result by doing or not doing a fit
699 if (testStatType == 0) {
700 if (!doFit)
701 Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit is not done and the TS will use "
702 "the nuisances at the model value");
703 else
704 Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit has been done and the TS will use "
705 "the nuisances at the best fit value");
706 }
707
708 // build test statistics and hypotest calculators for running the inverter
709
710 SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(), *bModel->GetPdf());
711
712 // null parameters must includes snapshot of poi plus the nuisance values
713 RooArgSet nullParams(*sbModel->GetSnapshot());
714 if (sbModel->GetNuisanceParameters())
715 nullParams.add(*sbModel->GetNuisanceParameters());
716 if (sbModel->GetSnapshot())
717 slrts.SetNullParameters(nullParams);
718 RooArgSet altParams(*bModel->GetSnapshot());
719 if (bModel->GetNuisanceParameters())
720 altParams.add(*bModel->GetNuisanceParameters());
721 if (bModel->GetSnapshot())
722 slrts.SetAltParameters(altParams);
724 slrts.EnableDetailedOutput();
725
726 // ratio of profile likelihood - need to pass snapshot for the alt
727 RatioOfProfiledLikelihoodsTestStat ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
728 ropl.SetSubtractMLE(false);
729 if (testStatType == 11)
730 ropl.SetSubtractMLE(true);
731 ropl.SetPrintLevel(mPrintLevel);
732 ropl.SetMinimizer(mMinimizerType.c_str());
734 ropl.EnableDetailedOutput();
735
737 if (testStatType == 3)
738 profll.SetOneSided(true);
739 if (testStatType == 4)
740 profll.SetSigned(true);
741 profll.SetMinimizer(mMinimizerType.c_str());
742 profll.SetPrintLevel(mPrintLevel);
744 profll.EnableDetailedOutput();
745
746 profll.SetReuseNLL(mOptimize);
747 slrts.SetReuseNLL(mOptimize);
748 ropl.SetReuseNLL(mOptimize);
749
750 if (mOptimize) {
751 profll.SetStrategy(0);
752 ropl.SetStrategy(0);
754 }
755
756 if (mMaxPoi > 0)
757 poi->setMax(mMaxPoi); // increase limit
758
761
762 AsymptoticCalculator::SetPrintLevel(mPrintLevel);
763
764 // create the HypoTest calculator class
765 HypoTestCalculatorGeneric *hc = nullptr;
766 if (type == 0)
768 else if (type == 1)
770 // else if (type == 2 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false, mAsimovBins);
771 // else if (type == 3 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, true, mAsimovBins); // for using
772 // Asimov data generated with nominal values
773 else if (type == 2)
774 hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false);
775 else if (type == 3)
777 true); // for using Asimov data generated with nominal values
778 else {
779 Error("StandardHypoTestInvDemo", "Invalid - calculator type = %d supported values are only :\n\t\t\t 0 "
780 "(Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",
781 type);
782 return nullptr;
783 }
784
785 // set the test statistic
786 TestStatistic *testStat = nullptr;
787 if (testStatType == 0)
788 testStat = &slrts;
789 if (testStatType == 1 || testStatType == 11)
790 testStat = &ropl;
791 if (testStatType == 2 || testStatType == 3 || testStatType == 4)
792 testStat = &profll;
793 if (testStatType == 5)
794 testStat = &maxll;
795 if (testStatType == 6)
796 testStat = &nevtts;
797
798 if (testStat == nullptr) {
799 Error("StandardHypoTestInvDemo", "Invalid - test statistic type = %d supported values are only :\n\t\t\t 0 (SLR) "
800 ", 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",
802 return nullptr;
803 }
804
805 ToyMCSampler *toymcs = (ToyMCSampler *)hc->GetTestStatSampler();
806 if (toymcs && (type == 0 || type == 1)) {
807 // look if pdf is number counting or extended
808 if (sbModel->GetPdf()->canBeExtended()) {
810 Warning("StandardHypoTestInvDemo", "Pdf is extended: but number counting flag is set: ignore it ");
811 } else {
812 // for not extended pdf
813 if (!useNumberCounting) {
814 int nEvents = data->numEntries();
815 Info("StandardHypoTestInvDemo",
816 "Pdf is not extended: number of events to generate taken from observed data set is %d", nEvents);
817 toymcs->SetNEventsPerToy(nEvents);
818 } else {
819 Info("StandardHypoTestInvDemo", "using a number counting pdf");
820 toymcs->SetNEventsPerToy(1);
821 }
822 }
823
824 toymcs->SetTestStatistic(testStat);
825
826 if (data->isWeighted() && !mGenerateBinned) {
827 Info("StandardHypoTestInvDemo", "Data set is weighted, nentries = %d and sum of weights = %8.1f but toy "
828 "generation is unbinned - it would be faster to set mGenerateBinned to true\n",
829 data->numEntries(), data->sumEntries());
830 }
831 toymcs->SetGenerateBinned(mGenerateBinned);
832
833 toymcs->SetUseMultiGen(mOptimize);
834
835 if (mGenerateBinned && sbModel->GetObservables()->getSize() > 2) {
836 Warning("StandardHypoTestInvDemo", "generate binned is activated but the number of observable is %d. Too much "
837 "memory could be needed for allocating all the bins",
838 sbModel->GetObservables()->getSize());
839 }
840
841 // set the random seed if needed
842 if (mRandomSeed >= 0)
844 }
845
846 // specify if need to re-use same toys
847 if (mReuseAltToys) {
848 hc->UseSameAltToys();
849 }
850
851 if (type == 1) {
852 HybridCalculator *hhc = dynamic_cast<HybridCalculator *>(hc);
853 assert(hhc);
854
855 hhc->SetToys(ntoys, ntoys / mNToysRatio); // can use less ntoys for b hypothesis
856
857 // remove global observables from ModelConfig (this is probably not needed anymore in 5.32)
858 bModel->SetGlobalObservables(RooArgSet());
859 sbModel->SetGlobalObservables(RooArgSet());
860
861 // check for nuisance prior pdf in case of nuisance parameters
862 if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters()) {
863
864 // fix for using multigen (does not work in this case)
865 toymcs->SetUseMultiGen(false);
866 ToyMCSampler::SetAlwaysUseMultiGen(false);
867
868 RooAbsPdf *nuisPdf = nullptr;
869 if (nuisPriorName)
870 nuisPdf = w->pdf(nuisPriorName);
871 // use prior defined first in bModel (then in SbModel)
872 if (!nuisPdf) {
873 Info("StandardHypoTestInvDemo",
874 "No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
875 if (bModel->GetPdf() && bModel->GetObservables())
876 nuisPdf = RooStats::MakeNuisancePdf(*bModel, "nuisancePdf_bmodel");
877 else
878 nuisPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisancePdf_sbmodel");
879 }
880 if (!nuisPdf) {
881 if (bModel->GetPriorPdf()) {
882 nuisPdf = bModel->GetPriorPdf();
883 Info("StandardHypoTestInvDemo",
884 "No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
885 nuisPdf->GetName());
886 } else {
887 Error("StandardHypoTestInvDemo", "Cannot run Hybrid calculator because no prior on the nuisance "
888 "parameter is specified or can be derived");
889 return nullptr;
890 }
891 }
893 Info("StandardHypoTestInvDemo", "Using as nuisance Pdf ... ");
894 nuisPdf->Print();
895
896 const RooArgSet *nuisParams =
897 (bModel->GetNuisanceParameters()) ? bModel->GetNuisanceParameters() : sbModel->GetNuisanceParameters();
898 std::unique_ptr<RooArgSet> np{nuisPdf->getObservables(*nuisParams)};
899 if (np->getSize() == 0) {
900 Warning("StandardHypoTestInvDemo",
901 "Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
902 }
903
904 hhc->ForcePriorNuisanceAlt(*nuisPdf);
905 hhc->ForcePriorNuisanceNull(*nuisPdf);
906 }
907 } else if (type == 2 || type == 3) {
908 if (testStatType == 3)
909 ((AsymptoticCalculator *)hc)->SetOneSided(true);
910 if (testStatType != 2 && testStatType != 3)
911 Warning("StandardHypoTestInvDemo",
912 "Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
913 } else if (type == 0) {
915 // store also the fit information for each poi point used by calculator based on toys
917 ((FrequentistCalculator *)hc)->StoreFitInfo(true);
918 } else if (type == 1) {
919 ((HybridCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
920 // store also the fit information for each poi point used by calculator based on toys
921 // if (mEnableDetOutput) ((HybridCalculator*) hc)->StoreFitInfo(true);
922 }
923
924 // Get the result
925 RooMsgService::instance().getStream(1).removeTopic(RooFit::NumIntegration);
926
928 calc.SetConfidenceLevel(optHTInv.confLevel);
929
930 calc.UseCLs(useCLs);
931 calc.SetVerbose(true);
932
933 if (npoints > 0) {
934 if (poimin > poimax) {
935 // if no min/max given scan between MLE and +4 sigma
936 poimin = int(poihat);
937 poimax = int(poihat + 4 * poi->getError());
938 }
939 std::cout << "Doing a fixed scan in interval : " << poimin << " , " << poimax << std::endl;
940 calc.SetFixedScan(npoints, poimin, poimax);
941 } else {
942 // poi->setMax(10*int( (poihat+ 10 *poi->getError() )/10 ) );
943 std::cout << "Doing an automatic scan in interval : " << poi->getMin() << " , " << poi->getMax() << std::endl;
944 }
945
946 tw.Start();
947 HypoTestInverterResult *r = calc.GetInterval();
948 std::cout << "Time to perform limit scan \n";
949 tw.Print();
950
951 if (mRebuild) {
952
953 std::cout << "\n***************************************************************\n";
954 std::cout << "Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute "
955 "for each of them a new upper limit\n\n";
956
957 allParams = std::unique_ptr<RooArgSet>{sbModel->GetPdf()->getParameters(*data)};
958
959 // define on which value of nuisance parameters to do the rebuild
960 // default is best fit value for bmodel snapshot
961
962 if (mRebuildParamValues != 0) {
963 // set all parameters to their initial workspace values
965 }
966 if (mRebuildParamValues == 0 || mRebuildParamValues == 1) {
968 if (sbModel->GetNuisanceParameters())
969 constrainParams.add(*sbModel->GetNuisanceParameters());
971
972 const RooArgSet *poiModel = sbModel->GetParametersOfInterest();
973 bModel->LoadSnapshot();
974
975 // do a profile using the B model snapshot
976 if (mRebuildParamValues == 0) {
977
979
980 sbModel->GetPdf()->fitTo(*data, InitialHesse(false), Hesse(false),
981 Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0), PrintLevel(mPrintLevel),
983
984 std::cout << "rebuild using fitted parameter value for B-model snapshot" << std::endl;
985 constrainParams.Print("v");
986
988 }
989 }
990 std::cout << "StandardHypoTestInvDemo: Initial parameters used for rebuilding: ";
992
993 tw.Start();
994 SamplingDistribution *limDist = calc.GetUpperLimitDistribution(true, mNToyToRebuild);
995 std::cout << "Time to rebuild distributions " << std::endl;
996 tw.Print();
997
998 if (limDist) {
999 std::cout << "Expected limits after rebuild distribution " << std::endl;
1000 std::cout << "expected upper limit (median of limit distribution) " << limDist->InverseCDF(0.5) << std::endl;
1001 std::cout << "expected -1 sig limit (0.16% quantile of limit dist) "
1002 << limDist->InverseCDF(ROOT::Math::normal_cdf(-1)) << std::endl;
1003 std::cout << "expected +1 sig limit (0.84% quantile of limit dist) "
1004 << limDist->InverseCDF(ROOT::Math::normal_cdf(1)) << std::endl;
1005 std::cout << "expected -2 sig limit (.025% quantile of limit dist) "
1006 << limDist->InverseCDF(ROOT::Math::normal_cdf(-2)) << std::endl;
1007 std::cout << "expected +2 sig limit (.975% quantile of limit dist) "
1008 << limDist->InverseCDF(ROOT::Math::normal_cdf(2)) << std::endl;
1009
1010 // Plot the upper limit distribution
1011 SamplingDistPlot limPlot((mNToyToRebuild < 200) ? 50 : 100);
1012 limPlot.AddSamplingDistribution(limDist);
1013 limPlot.GetTH1F()->SetStats(true); // display statistics
1014 limPlot.SetLineColor(kBlue);
1015 new TCanvas("limPlot", "Upper Limit Distribution");
1016 limPlot.Draw();
1017
1018 /// save result in a file
1019 limDist->SetName("RULDist");
1020 TFile("RULDist.root", "RECREATE");
1021 limDist->Write();
1022 fileOut->Close();
1023
1024 // update r to a new updated result object containing the rebuilt expected p-values distributions
1025 // (it will not recompute the expected limit)
1026 if (r)
1027 delete r; // need to delete previous object since GetInterval will return a cloned copy
1028 r = calc.GetInterval();
1029
1030 } else
1031 std::cout << "ERROR : failed to re-build distributions " << std::endl;
1032 }
1033
1034 return r;
1035}
1036
1037void ReadResult(const char *fileName, const char *resultName = "", bool useCLs = true)
1038{
1039 // read a previous stored result from a file given the result name
1040
1041 useCLs);
1042}
1043
1044#ifdef USE_AS_MAIN
1045int main()
1046{
1047 StandardHypoTestInvDemo();
1048}
1049#endif
int main()
Definition Prototype.cxx:12
@ kBlue
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
winID h TVirtualViewer3D TVirtualGLPainter char TVirtualGLPainter plot
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define gPad
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
static void SetDefaultStrategy(int strat)
Set the default strategy.
static const std::string & DefaultMinimizerType()
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
static void setDefaultStorageType(StorageType s)
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
static RooMsgService & instance()
Return reference to singleton instance.
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
double getError() const
Definition RooRealVar.h:58
void setMax(const char *name, double value)
Set maximum of name range to given value.
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio.
Does a frequentist hypothesis test.
Same purpose as HybridCalculatorOriginal, but different implementation.
Common base class for the Hypothesis Test Calculators.
Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
A class for performing a hypothesis test inversion by scanning the hypothesis test results of a HypoT...
MaxLikelihoodEstimateTestStat: TestStatistic that returns maximum likelihood estimate of a specified ...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
NumEventsTestStat is a simple implementation of the TestStatistic interface used for simple number co...
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
TestStatistic that returns the ratio of profiled likelihoods.
This class provides simple and straightforward utilities to plot SamplingDistribution objects.
This class simply holds a sampling distribution of some test statistic.
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
ToyMCSampler is an implementation of the TestStatSampler interface.
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:131
void ls(Option_t *option="") const override
List file contents.
Definition TFile.cxx:1462
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4130
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Mother of all ROOT objects.
Definition TObject.h:41
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:139
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition TSystem.cxx:1296
RooCmdArg InitialHesse(bool flag=true)
RooCmdArg Offset(std::string const &mode)
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg Minimizer(const char *type, const char *alg=nullptr)
RooCmdArg Hesse(bool flag=true)
RooCmdArg Strategy(Int_t code)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
std::ostream & Info()
Definition hadd.cxx:171
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
@ NumIntegration
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
std::string const & NLLOffsetMode()
Test what offsetting mode RooStats should use by default.
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
void SetNLLOffsetMode(std::string const &mode)
Function to set a global flag in RooStats to use NLL offset when performing nll computations.
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
useful function to print in one line the content of a set with their values
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:678