Logo ROOT  
Reference Guide
 
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)
86 bool useProof = false; // use Proof Lite when using toys (for freq or hybrid)
87 int nworkers = 0; // number of worker for ProofLite (default use all available cores)
88 bool enableDetailedOutput =
89 false; // enable detailed output with all fit information for each toys (output will be written in result file)
90 bool rebuild = false; // re-do extra toys for computing expected limits and rebuild test stat
91 // distributions (N.B this requires much more CPU (factor is equivalent to nToyToRebuild)
92 int nToyToRebuild = 100; // number of toys used to rebuild
93 int rebuildParamValues = 0; // = 0 do a profile of all the parameters on the B (alt snapshot) before performing a
94 // rebuild operation (default)
95 // = 1 use initial workspace parameters with B snapshot values
96 // = 2 use all initial workspace parameters with B
97 // Otherwise the rebuild will be performed using
98 int initialFit = -1; // do a first fit to the model (-1 : default, 0 skip fit, 1 do always fit)
99 int randomSeed = -1; // random seed (if = -1: use default value, if = 0 always random )
100 // NOTE: Proof uses automatically a random seed
101
102 int nAsimovBins = 0; // number of bins in observables used for Asimov data sets (0 is the default and it is given by
103 // workspace, typically is 100)
104
105 bool reuseAltToys = false; // reuse same toys for alternate hypothesis (if set one gets more stable bands)
106 double confLevel = 0.95; // confidence level value
107
108 std::string minimizerType =
109 ""; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
110 std::string massValue = ""; // extra string to tag output file of result
111 int printLevel = 0; // print level for debugging PL test statistics and calculators
112
113 bool useNLLOffset = false; // use NLL offset when fitting (this increase stability of fits)
114};
115
116HypoTestInvOptions optHTInv;
117
118// internal class to run the inverter and more
119
120namespace RooStats {
121
122class HypoTestInvTool {
123
124public:
125 HypoTestInvTool();
126 ~HypoTestInvTool(){};
127
128 HypoTestInverterResult *RunInverter(RooWorkspace *w, const char *modelSBName, const char *modelBName,
129 const char *dataName, int type, int testStatType, bool useCLs, int npoints,
130 double poimin, double poimax, int ntoys, bool useNumberCounting = false,
131 const char *nuisPriorName = nullptr);
132
133 void AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType, bool useCLs, int npoints,
134 const char *fileNameBase = nullptr);
135
136 void SetParameter(const char *name, const char *value);
137 void SetParameter(const char *name, bool value);
138 void SetParameter(const char *name, int value);
139 void SetParameter(const char *name, double value);
140
141private:
142 bool mPlotHypoTestResult;
143 bool mWriteResult;
144 bool mOptimize;
145 bool mUseVectorStore;
146 bool mGenerateBinned;
147 bool mUseProof;
148 bool mRebuild;
149 bool mReuseAltToys;
150 bool mEnableDetOutput;
151 int mNWorkers;
152 int mNToyToRebuild;
153 int mRebuildParamValues;
154 int mPrintLevel;
155 int mInitialFit;
156 int mRandomSeed;
157 double mNToysRatio;
158 double mMaxPoi;
159 int mAsimovBins;
160 std::string mMassValue;
161 std::string
162 mMinimizerType; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
163 TString mResultFileName;
164};
165
166} // end namespace RooStats
167
168RooStats::HypoTestInvTool::HypoTestInvTool()
169 : mPlotHypoTestResult(true), mWriteResult(false), mOptimize(true), mUseVectorStore(true), mGenerateBinned(false),
170 mUseProof(false), mEnableDetOutput(false), mRebuild(false), mReuseAltToys(false), mNWorkers(4),
171 mNToyToRebuild(100), mRebuildParamValues(0), mPrintLevel(0), mInitialFit(-1), mRandomSeed(-1), mNToysRatio(2),
172 mMaxPoi(-1), mAsimovBins(0), mMassValue(""), mMinimizerType(""), mResultFileName()
173{
174}
175
176void RooStats::HypoTestInvTool::SetParameter(const char *name, bool value)
177{
178 //
179 // set boolean parameters
180 //
181
182 std::string s_name(name);
183
184 if (s_name.find("PlotHypoTestResult") != std::string::npos)
185 mPlotHypoTestResult = value;
186 if (s_name.find("WriteResult") != std::string::npos)
187 mWriteResult = value;
188 if (s_name.find("Optimize") != std::string::npos)
189 mOptimize = value;
190 if (s_name.find("UseVectorStore") != std::string::npos)
191 mUseVectorStore = value;
192 if (s_name.find("GenerateBinned") != std::string::npos)
193 mGenerateBinned = value;
194 if (s_name.find("UseProof") != std::string::npos)
195 mUseProof = value;
196 if (s_name.find("EnableDetailedOutput") != std::string::npos)
197 mEnableDetOutput = value;
198 if (s_name.find("Rebuild") != std::string::npos)
199 mRebuild = value;
200 if (s_name.find("ReuseAltToys") != std::string::npos)
201 mReuseAltToys = value;
202
203 return;
204}
205
206void RooStats::HypoTestInvTool::SetParameter(const char *name, int value)
207{
208 //
209 // set integer parameters
210 //
211
212 std::string s_name(name);
213
214 if (s_name.find("NWorkers") != std::string::npos)
215 mNWorkers = value;
216 if (s_name.find("NToyToRebuild") != std::string::npos)
217 mNToyToRebuild = value;
218 if (s_name.find("RebuildParamValues") != std::string::npos)
219 mRebuildParamValues = value;
220 if (s_name.find("PrintLevel") != std::string::npos)
221 mPrintLevel = value;
222 if (s_name.find("InitialFit") != std::string::npos)
223 mInitialFit = value;
224 if (s_name.find("RandomSeed") != std::string::npos)
225 mRandomSeed = value;
226 if (s_name.find("AsimovBins") != std::string::npos)
227 mAsimovBins = value;
228
229 return;
230}
231
232void RooStats::HypoTestInvTool::SetParameter(const char *name, double value)
233{
234 //
235 // set double precision parameters
236 //
237
238 std::string s_name(name);
239
240 if (s_name.find("NToysRatio") != std::string::npos)
241 mNToysRatio = value;
242 if (s_name.find("MaxPOI") != std::string::npos)
243 mMaxPoi = value;
244
245 return;
246}
247
248void RooStats::HypoTestInvTool::SetParameter(const char *name, const char *value)
249{
250 //
251 // set string parameters
252 //
253
254 std::string s_name(name);
255
256 if (s_name.find("MassValue") != std::string::npos)
257 mMassValue.assign(value);
258 if (s_name.find("MinimizerType") != std::string::npos)
259 mMinimizerType.assign(value);
260 if (s_name.find("ResultFileName") != std::string::npos)
261 mResultFileName = value;
262
263 return;
264}
265
266void StandardHypoTestInvDemo(const char *infile = nullptr, const char *wsName = "combined",
267 const char *modelSBName = "ModelConfig", const char *modelBName = "",
268 const char *dataName = "obsData", int calculatorType = 0, int testStatType = 0,
269 bool useCLs = true, int npoints = 6, double poimin = 0, double poimax = 5,
270 int ntoys = 1000, bool useNumberCounting = false, const char *nuisPriorName = nullptr)
271{
272 /*
273
274 Other Parameter to pass in tutorial
275 apart from standard for filename, ws, modelconfig and data
276
277 type = 0 Freq calculator
278 type = 1 Hybrid calculator
279 type = 2 Asymptotic calculator
280 type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
281
282 testStatType = 0 LEP
283 = 1 Tevatron
284 = 2 Profile Likelihood
285 = 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
286 = 4 Profiel Likelihood signed ( pll = -pll if mu < mu_hat)
287 = 5 Max Likelihood Estimate as test statistic
288 = 6 Number of observed event as test statistic
289
290 useCLs scan for CLs (otherwise for CLs+b)
291
292 npoints: number of points to scan , for autoscan set npoints = -1
293
294 poimin,poimax: min/max value to scan in case of fixed scans
295 (if min > max, try to find automatically)
296
297 ntoys: number of toys to use
298
299 useNumberCounting: set to true when using number counting events
300
301 nuisPriorName: name of prior for the nuisance. This is often expressed as constraint term in the global model
302 It is needed only when using the HybridCalculator (type=1)
303 If not given by default the prior pdf from ModelConfig is used.
304
305 extra options are available as global parameters of the macro. They major ones are:
306
307 plotHypoTestResult plot result of tests at each point (TS distributions) (default is true)
308 useProof use Proof (default is true)
309 writeResult write result of scan (default is true)
310 rebuild rebuild scan for expected limits (require extra toys) (default is false)
311 generateBinned generate binned data sets for toys (default is false) - be careful not to activate with
312 a too large (>=3) number of observables
313 nToyRatio ratio of S+B/B toys (default is 2)
314
315
316 */
317
318 TString filename(infile);
319 if (filename.IsNull()) {
320 filename = "results/example_combined_GaussExample_model.root";
321 bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
322 // if file does not exists generate with histfactory
323 if (!fileExist) {
324 // Normally this would be run on the command line
325 cout << "will run standard hist2workspace example" << endl;
326 gROOT->ProcessLine(".! prepareHistFactory .");
327 gROOT->ProcessLine(".! hist2workspace config/example.xml");
328 cout << "\n\n---------------------" << endl;
329 cout << "Done creating example input" << endl;
330 cout << "---------------------\n\n" << endl;
331 }
332
333 } else
334 filename = infile;
335
336 // Try to open the file
337 TFile *file = TFile::Open(filename);
338
339 // if input file was specified byt not found, quit
340 if (!file) {
341 cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
342 return;
343 }
344
345 HypoTestInvTool calc;
346
347 // set parameters
348 calc.SetParameter("PlotHypoTestResult", optHTInv.plotHypoTestResult);
349 calc.SetParameter("WriteResult", optHTInv.writeResult);
350 calc.SetParameter("Optimize", optHTInv.optimize);
351 calc.SetParameter("UseVectorStore", optHTInv.useVectorStore);
352 calc.SetParameter("GenerateBinned", optHTInv.generateBinned);
353 calc.SetParameter("NToysRatio", optHTInv.nToysRatio);
354 calc.SetParameter("MaxPOI", optHTInv.maxPOI);
355 calc.SetParameter("UseProof", optHTInv.useProof);
356 calc.SetParameter("EnableDetailedOutput", optHTInv.enableDetailedOutput);
357 calc.SetParameter("NWorkers", optHTInv.nworkers);
358 calc.SetParameter("Rebuild", optHTInv.rebuild);
359 calc.SetParameter("ReuseAltToys", optHTInv.reuseAltToys);
360 calc.SetParameter("NToyToRebuild", optHTInv.nToyToRebuild);
361 calc.SetParameter("RebuildParamValues", optHTInv.rebuildParamValues);
362 calc.SetParameter("MassValue", optHTInv.massValue.c_str());
363 calc.SetParameter("MinimizerType", optHTInv.minimizerType.c_str());
364 calc.SetParameter("PrintLevel", optHTInv.printLevel);
365 calc.SetParameter("InitialFit", optHTInv.initialFit);
366 calc.SetParameter("ResultFileName", optHTInv.resultFileName);
367 calc.SetParameter("RandomSeed", optHTInv.randomSeed);
368 calc.SetParameter("AsimovBins", optHTInv.nAsimovBins);
369
370 // enable offset for all roostats
371 if (optHTInv.useNLLOffset)
373
374 RooWorkspace *w = dynamic_cast<RooWorkspace *>(file->Get(wsName));
375 HypoTestInverterResult *r = nullptr;
376 std::cout << w << "\t" << filename << std::endl;
377 if (w != nullptr) {
378 r = calc.RunInverter(w, modelSBName, modelBName, dataName, calculatorType, testStatType, useCLs, npoints, poimin,
379 poimax, ntoys, useNumberCounting, nuisPriorName);
380 if (!r) {
381 std::cerr << "Error running the HypoTestInverter - Exit " << std::endl;
382 return;
383 }
384 } else {
385 // case workspace is not present look for the inverter result
386 std::cout << "Reading an HypoTestInverterResult with name " << wsName << " from file " << filename << std::endl;
387 r = dynamic_cast<HypoTestInverterResult *>(file->Get(wsName)); //
388 if (!r) {
389 std::cerr << "File " << filename << " does not contain a workspace or an HypoTestInverterResult - Exit "
390 << std::endl;
391 file->ls();
392 return;
393 }
394 }
395
396 calc.AnalyzeResult(r, calculatorType, testStatType, useCLs, npoints, infile);
397
398 return;
399}
400
401void RooStats::HypoTestInvTool::AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType,
402 bool useCLs, int npoints, const char *fileNameBase)
403{
404
405 // analyze result produced by the inverter, optionally save it in a file
406
407 double lowerLimit = 0;
408 double llError = 0;
409#if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
410 if (r->IsTwoSided()) {
411 lowerLimit = r->LowerLimit();
412 llError = r->LowerLimitEstimatedError();
413 }
414#else
415 lowerLimit = r->LowerLimit();
416 llError = r->LowerLimitEstimatedError();
417#endif
418
419 double upperLimit = r->UpperLimit();
420 double ulError = r->UpperLimitEstimatedError();
421
422 // std::cout << "DEBUG : [ " << lowerLimit << " , " << upperLimit << " ] " << std::endl;
423
424 if (lowerLimit < upperLimit * (1. - 1.E-4) && lowerLimit != 0)
425 std::cout << "The computed lower limit is: " << lowerLimit << " +/- " << llError << std::endl;
426 std::cout << "The computed upper limit is: " << upperLimit << " +/- " << ulError << std::endl;
427
428 // compute expected limit
429 std::cout << "Expected upper limits, using the B (alternate) model : " << std::endl;
430 std::cout << " expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
431 std::cout << " expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
432 std::cout << " expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
433 std::cout << " expected limit (-2 sig) " << r->GetExpectedUpperLimit(-2) << std::endl;
434 std::cout << " expected limit (+2 sig) " << r->GetExpectedUpperLimit(2) << std::endl;
435
436 // detailed output
437 if (mEnableDetOutput) {
438 mWriteResult = true;
439 Info("StandardHypoTestInvDemo", "detailed output will be written in output result file");
440 }
441
442 // write result in a file
443 if (r != nullptr && mWriteResult) {
444
445 // write to a file the results
446 const char *calcType = (calculatorType == 0) ? "Freq" : (calculatorType == 1) ? "Hybr" : "Asym";
447 const char *limitType = (useCLs) ? "CLs" : "Cls+b";
448 const char *scanType = (npoints < 0) ? "auto" : "grid";
449 if (mResultFileName.IsNull()) {
450 mResultFileName = TString::Format("%s_%s_%s_ts%d_", calcType, limitType, scanType, testStatType);
451 // strip the / from the filename
452 if (!mMassValue.empty()) {
453 mResultFileName += mMassValue.c_str();
454 mResultFileName += "_";
455 }
456
457 TString name = fileNameBase;
458 name.Replace(0, name.Last('/') + 1, "");
459 mResultFileName += name;
460 }
461
462 // get (if existing) rebuilt UL distribution
463 TString uldistFile = "RULDist.root";
464 TObject *ulDist = nullptr;
465 bool existULDist = !gSystem->AccessPathName(uldistFile);
466 if (existULDist) {
467 TFile *fileULDist = TFile::Open(uldistFile);
468 if (fileULDist)
469 ulDist = fileULDist->Get("RULDist");
470 }
471
472 TFile *fileOut = new TFile(mResultFileName, "RECREATE");
473 r->Write();
474 if (ulDist)
475 ulDist->Write();
476 Info("StandardHypoTestInvDemo", "HypoTestInverterResult has been written in the file %s", mResultFileName.Data());
477
478 fileOut->Close();
479 }
480
481 // plot the result ( p values vs scan points)
482 std::string typeName = "";
483 if (calculatorType == 0)
484 typeName = "Frequentist";
485 if (calculatorType == 1)
486 typeName = "Hybrid";
487 else if (calculatorType == 2 || calculatorType == 3) {
488 typeName = "Asymptotic";
489 mPlotHypoTestResult = false;
490 }
491
492 const char *resultName = r->GetName();
493 TString plotTitle = TString::Format("%s CL Scan for workspace %s", typeName.c_str(), resultName);
494 HypoTestInverterPlot *plot = new HypoTestInverterPlot("HTI_Result_Plot", plotTitle, r);
495
496 // plot in a new canvas with style
497 TString c1Name = TString::Format("%s_Scan", typeName.c_str());
498 TCanvas *c1 = new TCanvas(c1Name);
499 c1->SetLogy(false);
500
501 plot->Draw("CLb 2CL"); // plot all and Clb
502
503 // if (useCLs)
504 // plot->Draw("CLb 2CL"); // plot all and Clb
505 // else
506 // plot->Draw(""); // plot all and Clb
507
508 const int nEntries = r->ArraySize();
509
510 // plot test statistics distributions for the two hypothesis
511 if (mPlotHypoTestResult) {
512 TCanvas *c2 = new TCanvas("c2");
513 if (nEntries > 1) {
514 int ny = TMath::CeilNint(TMath::Sqrt(nEntries));
515 int nx = TMath::CeilNint(double(nEntries) / ny);
516 c2->Divide(nx, ny);
517 }
518 for (int i = 0; i < nEntries; i++) {
519 if (nEntries > 1)
520 c2->cd(i + 1);
521 SamplingDistPlot *pl = plot->MakeTestStatPlot(i);
522 pl->SetLogYaxis(true);
523 pl->Draw();
524 }
525 }
526 gPad = c1;
527}
528
529// internal routine to run the inverter
530HypoTestInverterResult *RooStats::HypoTestInvTool::RunInverter(RooWorkspace *w, const char *modelSBName,
531 const char *modelBName, const char *dataName, int type,
532 int testStatType, bool useCLs, int npoints,
533 double poimin, double poimax, int ntoys,
534 bool useNumberCounting, const char *nuisPriorName)
535{
536
537 std::cout << "Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
538
539 w->Print();
540
541 RooAbsData *data = w->data(dataName);
542 if (!data) {
543 Error("StandardHypoTestDemo", "Not existing data %s", dataName);
544 return nullptr;
545 } else
546 std::cout << "Using data set " << dataName << std::endl;
547
548 if (mUseVectorStore) {
550 data->convertToVectorStore();
551 }
552
553 // get models from WS
554 // get the modelConfig out of the file
555 ModelConfig *bModel = (ModelConfig *)w->obj(modelBName);
556 ModelConfig *sbModel = (ModelConfig *)w->obj(modelSBName);
557
558 if (!sbModel) {
559 Error("StandardHypoTestDemo", "Not existing ModelConfig %s", modelSBName);
560 return nullptr;
561 }
562 // check the model
563 if (!sbModel->GetPdf()) {
564 Error("StandardHypoTestDemo", "Model %s has no pdf ", modelSBName);
565 return nullptr;
566 }
567 if (!sbModel->GetParametersOfInterest()) {
568 Error("StandardHypoTestDemo", "Model %s has no poi ", modelSBName);
569 return nullptr;
570 }
571 if (!sbModel->GetObservables()) {
572 Error("StandardHypoTestInvDemo", "Model %s has no observables ", modelSBName);
573 return nullptr;
574 }
575 if (!sbModel->GetSnapshot()) {
576 Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi", modelSBName);
577 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
578 }
579
580 // case of no systematics
581 // remove nuisance parameters from model
582 if (optHTInv.noSystematics) {
583 const RooArgSet *nuisPar = sbModel->GetNuisanceParameters();
584 if (nuisPar && nuisPar->getSize() > 0) {
585 std::cout << "StandardHypoTestInvDemo"
586 << " - Switch off all systematics by setting them constant to their initial values" << std::endl;
587 RooStats::SetAllConstant(*nuisPar);
588 }
589 if (bModel) {
590 const RooArgSet *bnuisPar = bModel->GetNuisanceParameters();
591 if (bnuisPar)
592 RooStats::SetAllConstant(*bnuisPar);
593 }
594 }
595
596 if (!bModel || bModel == sbModel) {
597 Info("StandardHypoTestInvDemo", "The background model %s does not exist", modelBName);
598 Info("StandardHypoTestInvDemo", "Copy it from ModelConfig %s and set POI to zero", modelSBName);
599 bModel = (ModelConfig *)sbModel->Clone();
600 bModel->SetName(TString(modelSBName) + TString("_with_poi_0"));
601 RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
602 if (!var)
603 return nullptr;
604 double oldval = var->getVal();
605 var->setVal(0);
606 bModel->SetSnapshot(RooArgSet(*var));
607 var->setVal(oldval);
608 } else {
609 if (!bModel->GetSnapshot()) {
610 Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi and 0 values ",
611 modelBName);
612 RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
613 if (var) {
614 double oldval = var->getVal();
615 var->setVal(0);
616 bModel->SetSnapshot(RooArgSet(*var));
617 var->setVal(oldval);
618 } else {
619 Error("StandardHypoTestInvDemo", "Model %s has no valid poi", modelBName);
620 return nullptr;
621 }
622 }
623 }
624
625 // check model has global observables when there are nuisance pdf
626 // for the hybrid case the globals are not needed
627 if (type != 1) {
628 bool hasNuisParam = (sbModel->GetNuisanceParameters() && sbModel->GetNuisanceParameters()->getSize() > 0);
629 bool hasGlobalObs = (sbModel->GetGlobalObservables() && sbModel->GetGlobalObservables()->getSize() > 0);
630 if (hasNuisParam && !hasGlobalObs) {
631 // try to see if model has nuisance parameters first
632 RooAbsPdf *constrPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisanceConstraintPdf_sbmodel");
633 if (constrPdf) {
634 Warning("StandardHypoTestInvDemo", "Model %s has nuisance parameters but no global observables associated",
635 sbModel->GetName());
636 Warning("StandardHypoTestInvDemo",
637 "\tThe effect of the nuisance parameters will not be treated correctly ");
638 }
639 }
640 }
641
642 // save all initial parameters of the model including the global observables
643 RooArgSet initialParameters;
644 std::unique_ptr<RooArgSet> allParams{sbModel->GetPdf()->getParameters(*data)};
645 allParams->snapshot(initialParameters);
646
647 // run first a data fit
648
649 const RooArgSet *poiSet = sbModel->GetParametersOfInterest();
650 RooRealVar *poi = (RooRealVar *)poiSet->first();
651
652 std::cout << "StandardHypoTestInvDemo : POI initial value: " << poi->GetName() << " = " << poi->getVal()
653 << std::endl;
654
655 // fit the data first (need to use constraint )
656 TStopwatch tw;
657
658 bool doFit = mInitialFit;
659 if (testStatType == 0 && mInitialFit == -1)
660 doFit = false; // case of LEP test statistic
661 if (type == 3 && mInitialFit == -1)
662 doFit = false; // case of Asymptoticcalculator with nominal Asimov
663 double poihat = 0;
664
665 if (mMinimizerType.empty())
667 else
669
670 Info("StandardHypoTestInvDemo", "Using %s as minimizer for computing the test statistic",
672
673 if (doFit) {
674
675 // do the fit : By doing a fit the POI snapshot (for S+B) is set to the fit value
676 // and the nuisance parameters nominal values will be set to the fit value.
677 // This is relevant when using LEP test statistics
678
679 Info("StandardHypoTestInvDemo", " Doing a first fit to the observed data ");
680 RooArgSet constrainParams;
681 if (sbModel->GetNuisanceParameters())
682 constrainParams.add(*sbModel->GetNuisanceParameters());
683 RooStats::RemoveConstantParameters(&constrainParams);
684 tw.Start();
685 std::unique_ptr<RooFitResult> fitres{sbModel->GetPdf()->fitTo(
686 *data, InitialHesse(false), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0),
687 PrintLevel(mPrintLevel), Constrain(constrainParams), Save(true), Offset(RooStats::IsNLLOffset()))};
688 if (fitres->status() != 0) {
689 Warning("StandardHypoTestInvDemo",
690 "Fit to the model failed - try with strategy 1 and perform first an Hesse computation");
691 fitres = std::unique_ptr<RooFitResult>{sbModel->GetPdf()->fitTo(
692 *data, InitialHesse(true), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(1),
693 PrintLevel(mPrintLevel + 1), Constrain(constrainParams), Save(true), Offset(RooStats::IsNLLOffset()))};
694 }
695 if (fitres->status() != 0)
696 Warning("StandardHypoTestInvDemo", " Fit still failed - continue anyway.....");
697
698 poihat = poi->getVal();
699 std::cout << "StandardHypoTestInvDemo - Best Fit value : " << poi->GetName() << " = " << poihat << " +/- "
700 << poi->getError() << std::endl;
701 std::cout << "Time for fitting : ";
702 tw.Print();
703
704 // save best fit value in the poi snapshot
705 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
706 std::cout << "StandardHypoTestInvo: snapshot of S+B Model " << sbModel->GetName()
707 << " is set to the best fit value" << std::endl;
708 }
709
710 // print a message in case of LEP test statistics because it affects result by doing or not doing a fit
711 if (testStatType == 0) {
712 if (!doFit)
713 Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit is not done and the TS will use "
714 "the nuisances at the model value");
715 else
716 Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit has been done and the TS will use "
717 "the nuisances at the best fit value");
718 }
719
720 // build test statistics and hypotest calculators for running the inverter
721
722 SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(), *bModel->GetPdf());
723
724 // null parameters must includes snapshot of poi plus the nuisance values
725 RooArgSet nullParams(*sbModel->GetSnapshot());
726 if (sbModel->GetNuisanceParameters())
727 nullParams.add(*sbModel->GetNuisanceParameters());
728 if (sbModel->GetSnapshot())
729 slrts.SetNullParameters(nullParams);
730 RooArgSet altParams(*bModel->GetSnapshot());
731 if (bModel->GetNuisanceParameters())
732 altParams.add(*bModel->GetNuisanceParameters());
733 if (bModel->GetSnapshot())
734 slrts.SetAltParameters(altParams);
735 if (mEnableDetOutput)
736 slrts.EnableDetailedOutput();
737
738 // ratio of profile likelihood - need to pass snapshot for the alt
739 RatioOfProfiledLikelihoodsTestStat ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
740 ropl.SetSubtractMLE(false);
741 if (testStatType == 11)
742 ropl.SetSubtractMLE(true);
743 ropl.SetPrintLevel(mPrintLevel);
744 ropl.SetMinimizer(mMinimizerType.c_str());
745 if (mEnableDetOutput)
746 ropl.EnableDetailedOutput();
747
748 ProfileLikelihoodTestStat profll(*sbModel->GetPdf());
749 if (testStatType == 3)
750 profll.SetOneSided(true);
751 if (testStatType == 4)
752 profll.SetSigned(true);
753 profll.SetMinimizer(mMinimizerType.c_str());
754 profll.SetPrintLevel(mPrintLevel);
755 if (mEnableDetOutput)
756 profll.EnableDetailedOutput();
757
758 profll.SetReuseNLL(mOptimize);
759 slrts.SetReuseNLL(mOptimize);
760 ropl.SetReuseNLL(mOptimize);
761
762 if (mOptimize) {
763 profll.SetStrategy(0);
764 ropl.SetStrategy(0);
766 }
767
768 if (mMaxPoi > 0)
769 poi->setMax(mMaxPoi); // increase limit
770
771 MaxLikelihoodEstimateTestStat maxll(*sbModel->GetPdf(), *poi);
772 NumEventsTestStat nevtts;
773
774 AsymptoticCalculator::SetPrintLevel(mPrintLevel);
775
776 // create the HypoTest calculator class
777 HypoTestCalculatorGeneric *hc = nullptr;
778 if (type == 0)
779 hc = new FrequentistCalculator(*data, *bModel, *sbModel);
780 else if (type == 1)
781 hc = new HybridCalculator(*data, *bModel, *sbModel);
782 // else if (type == 2 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false, mAsimovBins);
783 // else if (type == 3 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, true, mAsimovBins); // for using
784 // Asimov data generated with nominal values
785 else if (type == 2)
786 hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false);
787 else if (type == 3)
788 hc = new AsymptoticCalculator(*data, *bModel, *sbModel,
789 true); // for using Asimov data generated with nominal values
790 else {
791 Error("StandardHypoTestInvDemo", "Invalid - calculator type = %d supported values are only :\n\t\t\t 0 "
792 "(Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",
793 type);
794 return nullptr;
795 }
796
797 // set the test statistic
798 TestStatistic *testStat = nullptr;
799 if (testStatType == 0)
800 testStat = &slrts;
801 if (testStatType == 1 || testStatType == 11)
802 testStat = &ropl;
803 if (testStatType == 2 || testStatType == 3 || testStatType == 4)
804 testStat = &profll;
805 if (testStatType == 5)
806 testStat = &maxll;
807 if (testStatType == 6)
808 testStat = &nevtts;
809
810 if (testStat == nullptr) {
811 Error("StandardHypoTestInvDemo", "Invalid - test statistic type = %d supported values are only :\n\t\t\t 0 (SLR) "
812 ", 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",
813 testStatType);
814 return nullptr;
815 }
816
818 if (toymcs && (type == 0 || type == 1)) {
819 // look if pdf is number counting or extended
820 if (sbModel->GetPdf()->canBeExtended()) {
821 if (useNumberCounting)
822 Warning("StandardHypoTestInvDemo", "Pdf is extended: but number counting flag is set: ignore it ");
823 } else {
824 // for not extended pdf
825 if (!useNumberCounting) {
826 int nEvents = data->numEntries();
827 Info("StandardHypoTestInvDemo",
828 "Pdf is not extended: number of events to generate taken from observed data set is %d", nEvents);
829 toymcs->SetNEventsPerToy(nEvents);
830 } else {
831 Info("StandardHypoTestInvDemo", "using a number counting pdf");
832 toymcs->SetNEventsPerToy(1);
833 }
834 }
835
836 toymcs->SetTestStatistic(testStat);
837
838 if (data->isWeighted() && !mGenerateBinned) {
839 Info("StandardHypoTestInvDemo", "Data set is weighted, nentries = %d and sum of weights = %8.1f but toy "
840 "generation is unbinned - it would be faster to set mGenerateBinned to true\n",
841 data->numEntries(), data->sumEntries());
842 }
843 toymcs->SetGenerateBinned(mGenerateBinned);
844
845 toymcs->SetUseMultiGen(mOptimize);
846
847 if (mGenerateBinned && sbModel->GetObservables()->getSize() > 2) {
848 Warning("StandardHypoTestInvDemo", "generate binned is activated but the number of observable is %d. Too much "
849 "memory could be needed for allocating all the bins",
850 sbModel->GetObservables()->getSize());
851 }
852
853 // set the random seed if needed
854 if (mRandomSeed >= 0)
855 RooRandom::randomGenerator()->SetSeed(mRandomSeed);
856 }
857
858 // specify if need to re-use same toys
859 if (mReuseAltToys) {
860 hc->UseSameAltToys();
861 }
862
863 if (type == 1) {
864 HybridCalculator *hhc = dynamic_cast<HybridCalculator *>(hc);
865 assert(hhc);
866
867 hhc->SetToys(ntoys, ntoys / mNToysRatio); // can use less ntoys for b hypothesis
868
869 // remove global observables from ModelConfig (this is probably not needed anymore in 5.32)
872
873 // check for nuisance prior pdf in case of nuisance parameters
874 if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters()) {
875
876 // fix for using multigen (does not work in this case)
877 toymcs->SetUseMultiGen(false);
878 ToyMCSampler::SetAlwaysUseMultiGen(false);
879
880 RooAbsPdf *nuisPdf = nullptr;
881 if (nuisPriorName)
882 nuisPdf = w->pdf(nuisPriorName);
883 // use prior defined first in bModel (then in SbModel)
884 if (!nuisPdf) {
885 Info("StandardHypoTestInvDemo",
886 "No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
887 if (bModel->GetPdf() && bModel->GetObservables())
888 nuisPdf = RooStats::MakeNuisancePdf(*bModel, "nuisancePdf_bmodel");
889 else
890 nuisPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisancePdf_sbmodel");
891 }
892 if (!nuisPdf) {
893 if (bModel->GetPriorPdf()) {
894 nuisPdf = bModel->GetPriorPdf();
895 Info("StandardHypoTestInvDemo",
896 "No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
897 nuisPdf->GetName());
898 } else {
899 Error("StandardHypoTestInvDemo", "Cannot run Hybrid calculator because no prior on the nuisance "
900 "parameter is specified or can be derived");
901 return nullptr;
902 }
903 }
904 assert(nuisPdf);
905 Info("StandardHypoTestInvDemo", "Using as nuisance Pdf ... ");
906 nuisPdf->Print();
907
908 const RooArgSet *nuisParams =
909 (bModel->GetNuisanceParameters()) ? bModel->GetNuisanceParameters() : sbModel->GetNuisanceParameters();
910 std::unique_ptr<RooArgSet> np{nuisPdf->getObservables(*nuisParams)};
911 if (np->getSize() == 0) {
912 Warning("StandardHypoTestInvDemo",
913 "Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
914 }
915
916 hhc->ForcePriorNuisanceAlt(*nuisPdf);
917 hhc->ForcePriorNuisanceNull(*nuisPdf);
918 }
919 } else if (type == 2 || type == 3) {
920 if (testStatType == 3)
921 ((AsymptoticCalculator *)hc)->SetOneSided(true);
922 if (testStatType != 2 && testStatType != 3)
923 Warning("StandardHypoTestInvDemo",
924 "Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
925 } else if (type == 0) {
926 ((FrequentistCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
927 // store also the fit information for each poi point used by calculator based on toys
928 if (mEnableDetOutput)
929 ((FrequentistCalculator *)hc)->StoreFitInfo(true);
930 } else if (type == 1) {
931 ((HybridCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
932 // store also the fit information for each poi point used by calculator based on toys
933 // if (mEnableDetOutput) ((HybridCalculator*) hc)->StoreFitInfo(true);
934 }
935
936 // Get the result
938
939 HypoTestInverter calc(*hc);
940 calc.SetConfidenceLevel(optHTInv.confLevel);
941
942 calc.UseCLs(useCLs);
943 calc.SetVerbose(true);
944
945 // can speed up using proof-lite
946 if (mUseProof) {
947 ProofConfig pc(*w, mNWorkers, "", kFALSE);
948 toymcs->SetProofConfig(&pc); // enable proof
949 }
950
951 if (npoints > 0) {
952 if (poimin > poimax) {
953 // if no min/max given scan between MLE and +4 sigma
954 poimin = int(poihat);
955 poimax = int(poihat + 4 * poi->getError());
956 }
957 std::cout << "Doing a fixed scan in interval : " << poimin << " , " << poimax << std::endl;
958 calc.SetFixedScan(npoints, poimin, poimax);
959 } else {
960 // poi->setMax(10*int( (poihat+ 10 *poi->getError() )/10 ) );
961 std::cout << "Doing an automatic scan in interval : " << poi->getMin() << " , " << poi->getMax() << std::endl;
962 }
963
964 tw.Start();
965 HypoTestInverterResult *r = calc.GetInterval();
966 std::cout << "Time to perform limit scan \n";
967 tw.Print();
968
969 if (mRebuild) {
970
971 std::cout << "\n***************************************************************\n";
972 std::cout << "Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute "
973 "for each of them a new upper limit\n\n";
974
975 allParams = std::unique_ptr<RooArgSet>{sbModel->GetPdf()->getParameters(*data)};
976
977 // define on which value of nuisance parameters to do the rebuild
978 // default is best fit value for bmodel snapshot
979
980 if (mRebuildParamValues != 0) {
981 // set all parameters to their initial workspace values
982 allParams->assign(initialParameters);
983 }
984 if (mRebuildParamValues == 0 || mRebuildParamValues == 1) {
985 RooArgSet constrainParams;
986 if (sbModel->GetNuisanceParameters())
987 constrainParams.add(*sbModel->GetNuisanceParameters());
988 RooStats::RemoveConstantParameters(&constrainParams);
989
990 const RooArgSet *poiModel = sbModel->GetParametersOfInterest();
991 bModel->LoadSnapshot();
992
993 // do a profile using the B model snapshot
994 if (mRebuildParamValues == 0) {
995
996 RooStats::SetAllConstant(*poiModel, true);
997
998 sbModel->GetPdf()->fitTo(*data, InitialHesse(false), Hesse(false),
999 Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0), PrintLevel(mPrintLevel),
1000 Constrain(constrainParams), Offset(RooStats::IsNLLOffset()));
1001
1002 std::cout << "rebuild using fitted parameter value for B-model snapshot" << std::endl;
1003 constrainParams.Print("v");
1004
1005 RooStats::SetAllConstant(*poiModel, false);
1006 }
1007 }
1008 std::cout << "StandardHypoTestInvDemo: Initial parameters used for rebuilding: ";
1009 RooStats::PrintListContent(*allParams, std::cout);
1010
1011 calc.SetCloseProof(true);
1012 tw.Start();
1013 SamplingDistribution *limDist = calc.GetUpperLimitDistribution(true, mNToyToRebuild);
1014 std::cout << "Time to rebuild distributions " << std::endl;
1015 tw.Print();
1016
1017 if (limDist) {
1018 std::cout << "Expected limits after rebuild distribution " << std::endl;
1019 std::cout << "expected upper limit (median of limit distribution) " << limDist->InverseCDF(0.5) << std::endl;
1020 std::cout << "expected -1 sig limit (0.16% quantile of limit dist) "
1021 << limDist->InverseCDF(ROOT::Math::normal_cdf(-1)) << std::endl;
1022 std::cout << "expected +1 sig limit (0.84% quantile of limit dist) "
1023 << limDist->InverseCDF(ROOT::Math::normal_cdf(1)) << std::endl;
1024 std::cout << "expected -2 sig limit (.025% quantile of limit dist) "
1025 << limDist->InverseCDF(ROOT::Math::normal_cdf(-2)) << std::endl;
1026 std::cout << "expected +2 sig limit (.975% quantile of limit dist) "
1027 << limDist->InverseCDF(ROOT::Math::normal_cdf(2)) << std::endl;
1028
1029 // Plot the upper limit distribution
1030 SamplingDistPlot limPlot((mNToyToRebuild < 200) ? 50 : 100);
1031 limPlot.AddSamplingDistribution(limDist);
1032 limPlot.GetTH1F()->SetStats(true); // display statistics
1033 limPlot.SetLineColor(kBlue);
1034 new TCanvas("limPlot", "Upper Limit Distribution");
1035 limPlot.Draw();
1036
1037 /// save result in a file
1038 limDist->SetName("RULDist");
1039 TFile *fileOut = new TFile("RULDist.root", "RECREATE");
1040 limDist->Write();
1041 fileOut->Close();
1042
1043 // update r to a new updated result object containing the rebuilt expected p-values distributions
1044 // (it will not recompute the expected limit)
1045 if (r)
1046 delete r; // need to delete previous object since GetInterval will return a cloned copy
1047 r = calc.GetInterval();
1048
1049 } else
1050 std::cout << "ERROR : failed to re-build distributions " << std::endl;
1051 }
1052
1053 return r;
1054}
1055
1056void ReadResult(const char *fileName, const char *resultName = "", bool useCLs = true)
1057{
1058 // read a previous stored result from a file given the result name
1059
1060 StandardHypoTestInvDemo(fileName, resultName, "", "", "", 0, 0, useCLs);
1061}
1062
1063#ifdef USE_AS_MAIN
1064int main()
1065{
1066 StandardHypoTestInvDemo();
1067}
1068#endif
int main()
Definition Prototype.cxx:12
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
@ kBlue
Definition Rtypes.h:66
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
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:561
#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()
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:294
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * first() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
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
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:218
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:157
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.
StreamConfig & getStream(Int_t id)
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:48
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.
virtual void ForcePriorNuisanceNull(RooAbsPdf &priorNuisance)
Override the distribution used for marginalizing nuisance parameters that is inferred from ModelConfi...
virtual void ForcePriorNuisanceAlt(RooAbsPdf &priorNuisance)
void SetToys(int toysNull, int toysAlt)
set number of toys
Common base class for the Hypothesis Test Calculators.
TestStatSampler * GetTestStatSampler(void) const
Returns instance of TestStatSampler.
void UseSameAltToys()
Set this for re-using always the same toys for alternate hypothesis in case of calls at different nul...
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
virtual void SetSnapshot(const RooArgSet &set)
Set parameter values for a particular hypothesis if using a common PDF by saving a snapshot in the wo...
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
ModelConfig * Clone(const char *name="") const override
clone
Definition ModelConfig.h:57
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
void LoadSnapshot() const
load the snapshot from ws if it exists
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
RooAbsPdf * GetPriorPdf() const
get parameters prior pdf (return nullptr if not existing)
virtual void SetGlobalObservables(const RooArgSet &set)
Specify the global observables.
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...
Holds configuration options for proof and proof-lite.
Definition ProofConfig.h:45
TestStatistic that returns the ratio of profiled likelihoods.
This class provides simple and straightforward utilities to plot SamplingDistribution objects.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
void SetLineColor(Color_t color, const SamplingDistribution *samplDist=nullptr)
Sets line color for given sampling distribution and fill color for the associated shaded TH1F.
void SetLogYaxis(bool ly)
changes plot to log scale on y axis
double AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST")
adds the sampling distribution and returns the scale factor
TH1F * GetTH1F(const SamplingDistribution *samplDist=nullptr)
Returns the TH1F associated with the give SamplingDistribution.
This class simply holds a sampling distribution of some test statistic.
double InverseCDF(double pvalue)
get the inverse of the Cumulative distribution function
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.
void SetProofConfig(ProofConfig *pc=nullptr)
calling with argument or nullptr deactivates proof
virtual void SetTestStatistic(TestStatistic *testStatistic, unsigned int i)
Set the TestStatistic (want the argument to be a function of the data & parameter points.
void SetGenerateBinned(bool binned=true)
control to use bin data generation (=> see RooFit::AllBinned() option)
virtual void SetNEventsPerToy(const Int_t nevents)
Forces the generation of exactly n events even for extended PDFs.
void SetUseMultiGen(bool flag)
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:53
void ls(Option_t *option="") const override
List file contents.
Definition TFile.cxx:1457
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:4089
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8990
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
Mother of all ROOT objects.
Definition TObject.h:41
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:886
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:280
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
Basic string class.
Definition TString.h:139
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition TString.h:694
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).
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 JSONIO.h:26
@ NumIntegration
Namespace for the RooStats classes.
Definition Asimov.h:19
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
void UseNLLOffset(bool on)
function to set a global flag in RooStats to use NLL offset when performing nll computations Note tha...
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set
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:662
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:674
void removeTopic(RooFit::MsgTopic oldTopic)