69 using namespace RooFit;
70 using namespace RooStats;
74 struct HypoTestInvOptions {
76 bool plotHypoTestResult =
true;
77 bool writeResult =
true;
78 TString resultFileName;
80 bool useVectorStore =
true;
81 bool generateBinned =
false;
82 bool noSystematics =
false;
84 double nToysRatio = 2;
86 bool useProof =
false;
88 bool enableDetailedOutput =
92 int nToyToRebuild = 100;
93 int rebuildParamValues = 0;
105 bool reuseAltToys =
false;
106 double confLevel = 0.95;
108 std::string minimizerType =
110 std::string massValue =
"";
113 bool useNLLOffset =
false;
116 HypoTestInvOptions optHTInv;
122 class HypoTestInvTool {
126 ~HypoTestInvTool(){};
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 = 0);
133 void AnalyzeResult(HypoTestInverterResult *r,
int calculatorType,
int testStatType,
bool useCLs,
int npoints,
134 const char *fileNameBase = 0);
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);
142 bool mPlotHypoTestResult;
145 bool mUseVectorStore;
146 bool mGenerateBinned;
150 bool mEnableDetOutput;
153 int mRebuildParamValues;
160 std::string mMassValue;
163 TString mResultFileName;
168 RooStats::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()
176 void RooStats::HypoTestInvTool::SetParameter(
const char *name,
bool value)
182 std::string s_name(name);
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)
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)
196 if (s_name.find(
"EnableDetailedOutput") != std::string::npos)
197 mEnableDetOutput = value;
198 if (s_name.find(
"Rebuild") != std::string::npos)
200 if (s_name.find(
"ReuseAltToys") != std::string::npos)
201 mReuseAltToys = value;
206 void RooStats::HypoTestInvTool::SetParameter(
const char *name,
int value)
212 std::string s_name(name);
214 if (s_name.find(
"NWorkers") != std::string::npos)
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)
222 if (s_name.find(
"InitialFit") != std::string::npos)
224 if (s_name.find(
"RandomSeed") != std::string::npos)
226 if (s_name.find(
"AsimovBins") != std::string::npos)
232 void RooStats::HypoTestInvTool::SetParameter(
const char *name,
double value)
238 std::string s_name(name);
240 if (s_name.find(
"NToysRatio") != std::string::npos)
242 if (s_name.find(
"MaxPOI") != std::string::npos)
248 void RooStats::HypoTestInvTool::SetParameter(
const char *name,
const char *value)
254 std::string s_name(name);
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;
266 void StandardHypoTestInvDemo(
const char *infile = 0,
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 = 0)
318 TString filename(infile);
319 if (filename.IsNull()) {
320 filename =
"results/example_combined_GaussExample_model.root";
321 bool fileExist = !gSystem->AccessPathName(filename);
325 cout <<
"HistFactory file cannot be generated on Windows - exit" << endl;
329 cout <<
"will run standard hist2workspace example" << endl;
330 gROOT->ProcessLine(
".! prepareHistFactory .");
331 gROOT->ProcessLine(
".! hist2workspace config/example.xml");
332 cout <<
"\n\n---------------------" << endl;
333 cout <<
"Done creating example input" << endl;
334 cout <<
"---------------------\n\n" << endl;
341 TFile *file = TFile::Open(filename);
345 cout <<
"StandardRooStatsDemoMacro: Input file " << filename <<
" is not found" << endl;
349 HypoTestInvTool calc;
352 calc.SetParameter(
"PlotHypoTestResult", optHTInv.plotHypoTestResult);
353 calc.SetParameter(
"WriteResult", optHTInv.writeResult);
354 calc.SetParameter(
"Optimize", optHTInv.optimize);
355 calc.SetParameter(
"UseVectorStore", optHTInv.useVectorStore);
356 calc.SetParameter(
"GenerateBinned", optHTInv.generateBinned);
357 calc.SetParameter(
"NToysRatio", optHTInv.nToysRatio);
358 calc.SetParameter(
"MaxPOI", optHTInv.maxPOI);
359 calc.SetParameter(
"UseProof", optHTInv.useProof);
360 calc.SetParameter(
"EnableDetailedOutput", optHTInv.enableDetailedOutput);
361 calc.SetParameter(
"NWorkers", optHTInv.nworkers);
362 calc.SetParameter(
"Rebuild", optHTInv.rebuild);
363 calc.SetParameter(
"ReuseAltToys", optHTInv.reuseAltToys);
364 calc.SetParameter(
"NToyToRebuild", optHTInv.nToyToRebuild);
365 calc.SetParameter(
"RebuildParamValues", optHTInv.rebuildParamValues);
366 calc.SetParameter(
"MassValue", optHTInv.massValue.c_str());
367 calc.SetParameter(
"MinimizerType", optHTInv.minimizerType.c_str());
368 calc.SetParameter(
"PrintLevel", optHTInv.printLevel);
369 calc.SetParameter(
"InitialFit", optHTInv.initialFit);
370 calc.SetParameter(
"ResultFileName", optHTInv.resultFileName);
371 calc.SetParameter(
"RandomSeed", optHTInv.randomSeed);
372 calc.SetParameter(
"AsimovBins", optHTInv.nAsimovBins);
375 if (optHTInv.useNLLOffset)
376 RooStats::UseNLLOffset(
true);
378 RooWorkspace *w =
dynamic_cast<RooWorkspace *
>(file->Get(wsName));
379 HypoTestInverterResult *r = 0;
380 std::cout << w <<
"\t" << filename << std::endl;
382 r = calc.RunInverter(w, modelSBName, modelBName, dataName, calculatorType, testStatType, useCLs, npoints, poimin,
383 poimax, ntoys, useNumberCounting, nuisPriorName);
385 std::cerr <<
"Error running the HypoTestInverter - Exit " << std::endl;
390 std::cout <<
"Reading an HypoTestInverterResult with name " << wsName <<
" from file " << filename << std::endl;
391 r =
dynamic_cast<HypoTestInverterResult *
>(file->Get(wsName));
393 std::cerr <<
"File " << filename <<
" does not contain a workspace or an HypoTestInverterResult - Exit "
400 calc.AnalyzeResult(r, calculatorType, testStatType, useCLs, npoints, infile);
405 void RooStats::HypoTestInvTool::AnalyzeResult(HypoTestInverterResult *r,
int calculatorType,
int testStatType,
406 bool useCLs,
int npoints,
const char *fileNameBase)
411 double lowerLimit = 0;
413 #if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
414 if (r->IsTwoSided()) {
415 lowerLimit = r->LowerLimit();
416 llError = r->LowerLimitEstimatedError();
419 lowerLimit = r->LowerLimit();
420 llError = r->LowerLimitEstimatedError();
423 double upperLimit = r->UpperLimit();
424 double ulError = r->UpperLimitEstimatedError();
428 if (lowerLimit < upperLimit * (1. - 1.E-4) && lowerLimit != 0)
429 std::cout <<
"The computed lower limit is: " << lowerLimit <<
" +/- " << llError << std::endl;
430 std::cout <<
"The computed upper limit is: " << upperLimit <<
" +/- " << ulError << std::endl;
433 std::cout <<
"Expected upper limits, using the B (alternate) model : " << std::endl;
434 std::cout <<
" expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
435 std::cout <<
" expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
436 std::cout <<
" expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
437 std::cout <<
" expected limit (-2 sig) " << r->GetExpectedUpperLimit(-2) << std::endl;
438 std::cout <<
" expected limit (+2 sig) " << r->GetExpectedUpperLimit(2) << std::endl;
441 if (mEnableDetOutput) {
443 Info(
"StandardHypoTestInvDemo",
"detailed output will be written in output result file");
447 if (r != NULL && mWriteResult) {
450 const char *calcType = (calculatorType == 0) ?
"Freq" : (calculatorType == 1) ?
"Hybr" :
"Asym";
451 const char *limitType = (useCLs) ?
"CLs" :
"Cls+b";
452 const char *scanType = (npoints < 0) ?
"auto" :
"grid";
453 if (mResultFileName.IsNull()) {
454 mResultFileName = TString::Format(
"%s_%s_%s_ts%d_", calcType, limitType, scanType, testStatType);
456 if (mMassValue.size() > 0) {
457 mResultFileName += mMassValue.c_str();
458 mResultFileName +=
"_";
461 TString name = fileNameBase;
462 name.Replace(0, name.Last(
'/') + 1,
"");
463 mResultFileName += name;
467 TString uldistFile =
"RULDist.root";
469 bool existULDist = !gSystem->AccessPathName(uldistFile);
471 TFile *fileULDist = TFile::Open(uldistFile);
473 ulDist = fileULDist->Get(
"RULDist");
476 TFile *fileOut =
new TFile(mResultFileName,
"RECREATE");
480 Info(
"StandardHypoTestInvDemo",
"HypoTestInverterResult has been written in the file %s", mResultFileName.Data());
486 std::string typeName =
"";
487 if (calculatorType == 0)
488 typeName =
"Frequentist";
489 if (calculatorType == 1)
491 else if (calculatorType == 2 || calculatorType == 3) {
492 typeName =
"Asymptotic";
493 mPlotHypoTestResult =
false;
496 const char *resultName = r->GetName();
497 TString plotTitle = TString::Format(
"%s CL Scan for workspace %s", typeName.c_str(), resultName);
498 HypoTestInverterPlot *plot =
new HypoTestInverterPlot(
"HTI_Result_Plot", plotTitle, r);
501 TString c1Name = TString::Format(
"%s_Scan", typeName.c_str());
502 TCanvas *c1 =
new TCanvas(c1Name);
505 plot->Draw(
"CLb 2CL");
512 const int nEntries = r->ArraySize();
515 if (mPlotHypoTestResult) {
516 TCanvas *c2 =
new TCanvas(
"c2");
518 int ny = TMath::CeilNint(TMath::Sqrt(nEntries));
519 int nx = TMath::CeilNint(
double(nEntries) / ny);
522 for (
int i = 0; i < nEntries; i++) {
525 SamplingDistPlot *pl = plot->MakeTestStatPlot(i);
526 pl->SetLogYaxis(
true);
534 HypoTestInverterResult *RooStats::HypoTestInvTool::RunInverter(RooWorkspace *w,
const char *modelSBName,
535 const char *modelBName,
const char *dataName,
int type,
536 int testStatType,
bool useCLs,
int npoints,
537 double poimin,
double poimax,
int ntoys,
538 bool useNumberCounting,
const char *nuisPriorName)
541 std::cout <<
"Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
545 RooAbsData *data = w->data(dataName);
547 Error(
"StandardHypoTestDemo",
"Not existing data %s", dataName);
550 std::cout <<
"Using data set " << dataName << std::endl;
552 if (mUseVectorStore) {
553 RooAbsData::setDefaultStorageType(RooAbsData::Vector);
554 data->convertToVectorStore();
559 ModelConfig *bModel = (ModelConfig *)w->obj(modelBName);
560 ModelConfig *sbModel = (ModelConfig *)w->obj(modelSBName);
563 Error(
"StandardHypoTestDemo",
"Not existing ModelConfig %s", modelSBName);
567 if (!sbModel->GetPdf()) {
568 Error(
"StandardHypoTestDemo",
"Model %s has no pdf ", modelSBName);
571 if (!sbModel->GetParametersOfInterest()) {
572 Error(
"StandardHypoTestDemo",
"Model %s has no poi ", modelSBName);
575 if (!sbModel->GetObservables()) {
576 Error(
"StandardHypoTestInvDemo",
"Model %s has no observables ", modelSBName);
579 if (!sbModel->GetSnapshot()) {
580 Info(
"StandardHypoTestInvDemo",
"Model %s has no snapshot - make one using model poi", modelSBName);
581 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
586 if (optHTInv.noSystematics) {
587 const RooArgSet *nuisPar = sbModel->GetNuisanceParameters();
588 if (nuisPar && nuisPar->getSize() > 0) {
589 std::cout <<
"StandardHypoTestInvDemo"
590 <<
" - Switch off all systematics by setting them constant to their initial values" << std::endl;
591 RooStats::SetAllConstant(*nuisPar);
594 const RooArgSet *bnuisPar = bModel->GetNuisanceParameters();
596 RooStats::SetAllConstant(*bnuisPar);
600 if (!bModel || bModel == sbModel) {
601 Info(
"StandardHypoTestInvDemo",
"The background model %s does not exist", modelBName);
602 Info(
"StandardHypoTestInvDemo",
"Copy it from ModelConfig %s and set POI to zero", modelSBName);
603 bModel = (ModelConfig *)sbModel->Clone();
604 bModel->SetName(TString(modelSBName) + TString(
"_with_poi_0"));
605 RooRealVar *var =
dynamic_cast<RooRealVar *
>(bModel->GetParametersOfInterest()->first());
608 double oldval = var->getVal();
610 bModel->SetSnapshot(RooArgSet(*var));
613 if (!bModel->GetSnapshot()) {
614 Info(
"StandardHypoTestInvDemo",
"Model %s has no snapshot - make one using model poi and 0 values ",
616 RooRealVar *var =
dynamic_cast<RooRealVar *
>(bModel->GetParametersOfInterest()->first());
618 double oldval = var->getVal();
620 bModel->SetSnapshot(RooArgSet(*var));
623 Error(
"StandardHypoTestInvDemo",
"Model %s has no valid poi", modelBName);
632 bool hasNuisParam = (sbModel->GetNuisanceParameters() && sbModel->GetNuisanceParameters()->getSize() > 0);
633 bool hasGlobalObs = (sbModel->GetGlobalObservables() && sbModel->GetGlobalObservables()->getSize() > 0);
634 if (hasNuisParam && !hasGlobalObs) {
636 RooAbsPdf *constrPdf = RooStats::MakeNuisancePdf(*sbModel,
"nuisanceConstraintPdf_sbmodel");
638 Warning(
"StandardHypoTestInvDemo",
"Model %s has nuisance parameters but no global observables associated",
640 Warning(
"StandardHypoTestInvDemo",
641 "\tThe effect of the nuisance parameters will not be treated correctly ");
647 RooArgSet initialParameters;
648 RooArgSet *allParams = sbModel->GetPdf()->getParameters(*data);
649 allParams->snapshot(initialParameters);
654 const RooArgSet *poiSet = sbModel->GetParametersOfInterest();
655 RooRealVar *poi = (RooRealVar *)poiSet->first();
657 std::cout <<
"StandardHypoTestInvDemo : POI initial value: " << poi->GetName() <<
" = " << poi->getVal()
663 bool doFit = mInitialFit;
664 if (testStatType == 0 && mInitialFit == -1)
666 if (type == 3 && mInitialFit == -1)
670 if (mMinimizerType.size() == 0)
671 mMinimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
673 ROOT::Math::MinimizerOptions::SetDefaultMinimizer(mMinimizerType.c_str());
675 Info(
"StandardHypoTestInvDemo",
"Using %s as minimizer for computing the test statistic",
676 ROOT::Math::MinimizerOptions::DefaultMinimizerType().c_str());
684 Info(
"StandardHypoTestInvDemo",
" Doing a first fit to the observed data ");
685 RooArgSet constrainParams;
686 if (sbModel->GetNuisanceParameters())
687 constrainParams.add(*sbModel->GetNuisanceParameters());
688 RooStats::RemoveConstantParameters(&constrainParams);
690 RooFitResult *fitres = sbModel->GetPdf()->fitTo(
691 *data, InitialHesse(
false), Hesse(
false), Minimizer(mMinimizerType.c_str(),
"Migrad"), Strategy(0),
692 PrintLevel(mPrintLevel), Constrain(constrainParams), Save(
true), Offset(RooStats::IsNLLOffset()));
693 if (fitres->status() != 0) {
694 Warning(
"StandardHypoTestInvDemo",
695 "Fit to the model failed - try with strategy 1 and perform first an Hesse computation");
696 fitres = sbModel->GetPdf()->fitTo(
697 *data, InitialHesse(
true), Hesse(
false), Minimizer(mMinimizerType.c_str(),
"Migrad"), Strategy(1),
698 PrintLevel(mPrintLevel + 1), Constrain(constrainParams), Save(
true), Offset(RooStats::IsNLLOffset()));
700 if (fitres->status() != 0)
701 Warning(
"StandardHypoTestInvDemo",
" Fit still failed - continue anyway.....");
703 poihat = poi->getVal();
704 std::cout <<
"StandardHypoTestInvDemo - Best Fit value : " << poi->GetName() <<
" = " << poihat <<
" +/- "
705 << poi->getError() << std::endl;
706 std::cout <<
"Time for fitting : ";
710 sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
711 std::cout <<
"StandardHypoTestInvo: snapshot of S+B Model " << sbModel->GetName()
712 <<
" is set to the best fit value" << std::endl;
716 if (testStatType == 0) {
718 Info(
"StandardHypoTestInvDemo",
"Using LEP test statistic - an initial fit is not done and the TS will use "
719 "the nuisances at the model value");
721 Info(
"StandardHypoTestInvDemo",
"Using LEP test statistic - an initial fit has been done and the TS will use "
722 "the nuisances at the best fit value");
727 SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(), *bModel->GetPdf());
730 RooArgSet nullParams(*sbModel->GetSnapshot());
731 if (sbModel->GetNuisanceParameters())
732 nullParams.add(*sbModel->GetNuisanceParameters());
733 if (sbModel->GetSnapshot())
734 slrts.SetNullParameters(nullParams);
735 RooArgSet altParams(*bModel->GetSnapshot());
736 if (bModel->GetNuisanceParameters())
737 altParams.add(*bModel->GetNuisanceParameters());
738 if (bModel->GetSnapshot())
739 slrts.SetAltParameters(altParams);
740 if (mEnableDetOutput)
741 slrts.EnableDetailedOutput();
744 RatioOfProfiledLikelihoodsTestStat ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
745 ropl.SetSubtractMLE(
false);
746 if (testStatType == 11)
747 ropl.SetSubtractMLE(
true);
748 ropl.SetPrintLevel(mPrintLevel);
749 ropl.SetMinimizer(mMinimizerType.c_str());
750 if (mEnableDetOutput)
751 ropl.EnableDetailedOutput();
753 ProfileLikelihoodTestStat profll(*sbModel->GetPdf());
754 if (testStatType == 3)
755 profll.SetOneSided(
true);
756 if (testStatType == 4)
757 profll.SetSigned(
true);
758 profll.SetMinimizer(mMinimizerType.c_str());
759 profll.SetPrintLevel(mPrintLevel);
760 if (mEnableDetOutput)
761 profll.EnableDetailedOutput();
763 profll.SetReuseNLL(mOptimize);
764 slrts.SetReuseNLL(mOptimize);
765 ropl.SetReuseNLL(mOptimize);
768 profll.SetStrategy(0);
770 ROOT::Math::MinimizerOptions::SetDefaultStrategy(0);
774 poi->setMax(mMaxPoi);
776 MaxLikelihoodEstimateTestStat maxll(*sbModel->GetPdf(), *poi);
777 NumEventsTestStat nevtts;
779 AsymptoticCalculator::SetPrintLevel(mPrintLevel);
782 HypoTestCalculatorGeneric *hc = 0;
784 hc =
new FrequentistCalculator(*data, *bModel, *sbModel);
786 hc =
new HybridCalculator(*data, *bModel, *sbModel);
791 hc =
new AsymptoticCalculator(*data, *bModel, *sbModel,
false);
793 hc =
new AsymptoticCalculator(*data, *bModel, *sbModel,
796 Error(
"StandardHypoTestInvDemo",
"Invalid - calculator type = %d supported values are only :\n\t\t\t 0 "
797 "(Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",
803 TestStatistic *testStat = 0;
804 if (testStatType == 0)
806 if (testStatType == 1 || testStatType == 11)
808 if (testStatType == 2 || testStatType == 3 || testStatType == 4)
810 if (testStatType == 5)
812 if (testStatType == 6)
816 Error(
"StandardHypoTestInvDemo",
"Invalid - test statistic type = %d supported values are only :\n\t\t\t 0 (SLR) "
817 ", 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",
822 ToyMCSampler *toymcs = (ToyMCSampler *)hc->GetTestStatSampler();
823 if (toymcs && (type == 0 || type == 1)) {
825 if (sbModel->GetPdf()->canBeExtended()) {
826 if (useNumberCounting)
827 Warning(
"StandardHypoTestInvDemo",
"Pdf is extended: but number counting flag is set: ignore it ");
830 if (!useNumberCounting) {
831 int nEvents = data->numEntries();
832 Info(
"StandardHypoTestInvDemo",
833 "Pdf is not extended: number of events to generate taken from observed data set is %d", nEvents);
834 toymcs->SetNEventsPerToy(nEvents);
836 Info(
"StandardHypoTestInvDemo",
"using a number counting pdf");
837 toymcs->SetNEventsPerToy(1);
841 toymcs->SetTestStatistic(testStat);
843 if (data->isWeighted() && !mGenerateBinned) {
844 Info(
"StandardHypoTestInvDemo",
"Data set is weighted, nentries = %d and sum of weights = %8.1f but toy "
845 "generation is unbinned - it would be faster to set mGenerateBinned to true\n",
846 data->numEntries(), data->sumEntries());
848 toymcs->SetGenerateBinned(mGenerateBinned);
850 toymcs->SetUseMultiGen(mOptimize);
852 if (mGenerateBinned && sbModel->GetObservables()->getSize() > 2) {
853 Warning(
"StandardHypoTestInvDemo",
"generate binned is activated but the number of observable is %d. Too much "
854 "memory could be needed for allocating all the bins",
855 sbModel->GetObservables()->getSize());
859 if (mRandomSeed >= 0)
860 RooRandom::randomGenerator()->SetSeed(mRandomSeed);
865 hc->UseSameAltToys();
869 HybridCalculator *hhc =
dynamic_cast<HybridCalculator *
>(hc);
872 hhc->SetToys(ntoys, ntoys / mNToysRatio);
875 bModel->SetGlobalObservables(RooArgSet());
876 sbModel->SetGlobalObservables(RooArgSet());
879 if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters()) {
882 toymcs->SetUseMultiGen(
false);
883 ToyMCSampler::SetAlwaysUseMultiGen(
false);
885 RooAbsPdf *nuisPdf = 0;
887 nuisPdf = w->pdf(nuisPriorName);
890 Info(
"StandardHypoTestInvDemo",
891 "No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
892 if (bModel->GetPdf() && bModel->GetObservables())
893 nuisPdf = RooStats::MakeNuisancePdf(*bModel,
"nuisancePdf_bmodel");
895 nuisPdf = RooStats::MakeNuisancePdf(*sbModel,
"nuisancePdf_sbmodel");
898 if (bModel->GetPriorPdf()) {
899 nuisPdf = bModel->GetPriorPdf();
900 Info(
"StandardHypoTestInvDemo",
901 "No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
904 Error(
"StandardHypoTestInvDemo",
"Cannot run Hybrid calculator because no prior on the nuisance "
905 "parameter is specified or can be derived");
910 Info(
"StandardHypoTestInvDemo",
"Using as nuisance Pdf ... ");
913 const RooArgSet *nuisParams =
914 (bModel->GetNuisanceParameters()) ? bModel->GetNuisanceParameters() : sbModel->GetNuisanceParameters();
915 RooArgSet *np = nuisPdf->getObservables(*nuisParams);
916 if (np->getSize() == 0) {
917 Warning(
"StandardHypoTestInvDemo",
918 "Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
922 hhc->ForcePriorNuisanceAlt(*nuisPdf);
923 hhc->ForcePriorNuisanceNull(*nuisPdf);
925 }
else if (type == 2 || type == 3) {
926 if (testStatType == 3)
927 ((AsymptoticCalculator *)hc)->SetOneSided(
true);
928 if (testStatType != 2 && testStatType != 3)
929 Warning(
"StandardHypoTestInvDemo",
930 "Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
931 }
else if (type == 0) {
932 ((FrequentistCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
934 if (mEnableDetOutput)
935 ((FrequentistCalculator *)hc)->StoreFitInfo(
true);
936 }
else if (type == 1) {
937 ((HybridCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
943 RooMsgService::instance().getStream(1).removeTopic(RooFit::NumIntegration);
945 HypoTestInverter calc(*hc);
946 calc.SetConfidenceLevel(optHTInv.confLevel);
949 calc.SetVerbose(
true);
953 ProofConfig pc(*w, mNWorkers,
"", kFALSE);
954 toymcs->SetProofConfig(&pc);
958 if (poimin > poimax) {
960 poimin = int(poihat);
961 poimax = int(poihat + 4 * poi->getError());
963 std::cout <<
"Doing a fixed scan in interval : " << poimin <<
" , " << poimax << std::endl;
964 calc.SetFixedScan(npoints, poimin, poimax);
967 std::cout <<
"Doing an automatic scan in interval : " << poi->getMin() <<
" , " << poi->getMax() << std::endl;
971 HypoTestInverterResult *r = calc.GetInterval();
972 std::cout <<
"Time to perform limit scan \n";
977 std::cout <<
"\n***************************************************************\n";
978 std::cout <<
"Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute "
979 "for each of them a new upper limit\n\n";
981 allParams = sbModel->GetPdf()->getParameters(*data);
986 if (mRebuildParamValues != 0) {
988 *allParams = initialParameters;
990 if (mRebuildParamValues == 0 || mRebuildParamValues == 1) {
991 RooArgSet constrainParams;
992 if (sbModel->GetNuisanceParameters())
993 constrainParams.add(*sbModel->GetNuisanceParameters());
994 RooStats::RemoveConstantParameters(&constrainParams);
996 const RooArgSet *poiModel = sbModel->GetParametersOfInterest();
997 bModel->LoadSnapshot();
1000 if (mRebuildParamValues == 0) {
1002 RooStats::SetAllConstant(*poiModel,
true);
1004 sbModel->GetPdf()->fitTo(*data, InitialHesse(
false), Hesse(
false),
1005 Minimizer(mMinimizerType.c_str(),
"Migrad"), Strategy(0), PrintLevel(mPrintLevel),
1006 Constrain(constrainParams), Offset(RooStats::IsNLLOffset()));
1008 std::cout <<
"rebuild using fitted parameter value for B-model snapshot" << std::endl;
1009 constrainParams.Print(
"v");
1011 RooStats::SetAllConstant(*poiModel,
false);
1014 std::cout <<
"StandardHypoTestInvDemo: Initial parameters used for rebuilding: ";
1015 RooStats::PrintListContent(*allParams, std::cout);
1018 calc.SetCloseProof(1);
1020 SamplingDistribution *limDist = calc.GetUpperLimitDistribution(
true, mNToyToRebuild);
1021 std::cout <<
"Time to rebuild distributions " << std::endl;
1025 std::cout <<
"Expected limits after rebuild distribution " << std::endl;
1026 std::cout <<
"expected upper limit (median of limit distribution) " << limDist->InverseCDF(0.5) << std::endl;
1027 std::cout <<
"expected -1 sig limit (0.16% quantile of limit dist) "
1028 << limDist->InverseCDF(ROOT::Math::normal_cdf(-1)) << std::endl;
1029 std::cout <<
"expected +1 sig limit (0.84% quantile of limit dist) "
1030 << limDist->InverseCDF(ROOT::Math::normal_cdf(1)) << std::endl;
1031 std::cout <<
"expected -2 sig limit (.025% quantile of limit dist) "
1032 << limDist->InverseCDF(ROOT::Math::normal_cdf(-2)) << std::endl;
1033 std::cout <<
"expected +2 sig limit (.975% quantile of limit dist) "
1034 << limDist->InverseCDF(ROOT::Math::normal_cdf(2)) << std::endl;
1037 SamplingDistPlot limPlot((mNToyToRebuild < 200) ? 50 : 100);
1038 limPlot.AddSamplingDistribution(limDist);
1039 limPlot.GetTH1F()->SetStats(
true);
1040 limPlot.SetLineColor(kBlue);
1041 new TCanvas(
"limPlot",
"Upper Limit Distribution");
1045 limDist->SetName(
"RULDist");
1046 TFile *fileOut =
new TFile(
"RULDist.root",
"RECREATE");
1054 r = calc.GetInterval();
1057 std::cout <<
"ERROR : failed to re-build distributions " << std::endl;
1063 void ReadResult(
const char *fileName,
const char *resultName =
"",
bool useCLs =
true)
1067 StandardHypoTestInvDemo(fileName, resultName,
"",
"",
"", 0, 0, useCLs);
1073 StandardHypoTestInvDemo();