72 using namespace RooStats;
76 ClassImp(RooStats::AsymptoticCalculator);
78 int AsymptoticCalculator::fgPrintLevel = 1;
87 void AsymptoticCalculator::SetPrintLevel(
int level) {
94 AsymptoticCalculator::AsymptoticCalculator(
96 const ModelConfig &altModel,
97 const ModelConfig &nullModel,
bool nominalAsimov) :
98 HypoTestCalculatorGeneric(data, altModel, nullModel, 0),
99 fOneSided(false), fOneSidedDiscovery(false), fNominalAsimov(nominalAsimov),
101 fNLLObs(0), fNLLAsimov(0),
104 if (!Initialize())
return;
106 int verbose = fgPrintLevel;
109 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
110 assert(nullSnapshot);
111 RooRealVar * muNull =
dynamic_cast<RooRealVar*
>( nullSnapshot->first() );
113 if (muNull->getVal() == muNull->getMin()) {
114 fOneSidedDiscovery =
true;
116 oocoutI((TObject*)0,InputArguments) <<
"AsymptotiCalculator: Minimum of POI is " << muNull->getMin() <<
" corresponds to null snapshot - default configuration is one-sided discovery formulae " << std::endl;
132 bool AsymptoticCalculator::Initialize()
const {
134 int verbose = fgPrintLevel;
136 oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator::Initialize...." << std::endl;
139 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
141 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::Initialize - ModelConfig has not a pdf defined" << std::endl;
144 RooAbsData * obsData =
const_cast<RooAbsData *
>(GetData() );
146 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::Initialize - data set has not been defined" << std::endl;
149 RooAbsData & data = *obsData;
153 const RooArgSet * poi = GetNullModel()->GetParametersOfInterest();
154 if (!poi || poi->getSize() == 0) {
155 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::Initialize - ModelConfig has not POI defined." << endl;
158 if (poi->getSize() > 1) {
159 oocoutW((TObject*)0,InputArguments) <<
"AsymptoticCalculator::Initialize - ModelConfig has more than one POI defined \n\t"
160 <<
"The asymptotic calculator works for only one POI - consider as POI only the first parameter"
166 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
167 if(nullSnapshot == NULL || nullSnapshot->getSize() == 0) {
168 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::Initialize - Null model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
179 RooArgSet nominalParams;
180 RooArgSet * allParams = nullPdf->getParameters(data);
181 RemoveConstantParameters(allParams);
182 if (fNominalAsimov) {
183 allParams->snapshot(nominalParams);
185 fBestFitPoi.removeAll();
186 fBestFitParams.removeAll();
187 fAsimovGlobObs.removeAll();
191 oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator::Initialize - Find best unconditional NLL on observed data" << endl;
192 fNLLObs = EvaluateNLL( *nullPdf, data, GetNullModel()->GetConditionalObservables(),GetNullModel()->GetGlobalObservables());
194 poi->snapshot(fBestFitPoi);
195 RooRealVar * muBest =
dynamic_cast<RooRealVar*
>(fBestFitPoi.first());
198 oocoutP((TObject*)0,Eval) <<
"Best fitted POI value = " << muBest->getVal() <<
" +/- " << muBest->getError() << std::endl;
200 allParams->snapshot(fBestFitParams);
204 const RooArgSet * altSnapshot = GetAlternateModel()->GetSnapshot();
205 if(altSnapshot == NULL || altSnapshot->getSize() == 0) {
206 oocoutE((TObject*)0,InputArguments) <<
"Alt (Background) model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << endl;
210 RooArgSet poiAlt(*altSnapshot);
212 oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator: Building Asimov data Set" << endl;
219 RooRealVar * xobs = 0;
220 if (GetNullModel()->GetObservables() && GetNullModel()->GetObservables()->getSize() == 1 ) {
221 xobs = (RooRealVar*) (GetNullModel()->GetObservables())->first();
222 if (data.IsA() == RooDataHist::Class() ) {
223 if (data.numEntries() != xobs->getBins() ) {
224 prevBins = xobs->getBins();
225 oocoutW((TObject*)0,InputArguments) <<
"AsymptoticCalculator: number of bins in " << xobs->GetName() <<
" are different than data bins "
226 <<
" set the same data bins " << data.numEntries() <<
" in range "
227 <<
" [ " << xobs->getMin() <<
" , " << xobs->getMax() <<
" ]" << std::endl;
228 xobs->setBins(data.numEntries());
233 if (!fNominalAsimov) {
235 oocoutI((TObject*)0,InputArguments) <<
"AsymptoticCalculator: Asimov data will be generated using fitted nuisance parameter values" << endl;
236 RooArgSet * tmp = (RooArgSet*) poiAlt.snapshot();
237 fAsimovData = MakeAsimovData( data, *GetNullModel(), poiAlt, fAsimovGlobObs,tmp);
243 oocoutI((TObject*)0,InputArguments) <<
"AsymptoticCalculator: Asimovdata set will be generated using nominal (current) nuisance parameter values" << endl;
244 nominalParams = poiAlt;
245 fAsimovData = MakeAsimovData( *GetNullModel(), nominalParams, fAsimovGlobObs);
249 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator: Error : Asimov data set could not be generated " << endl;
255 RooArgSet globObsSnapshot;
256 if (GetNullModel()->GetGlobalObservables() ) {
257 globObs.add(*GetNullModel()->GetGlobalObservables());
258 assert(globObs.getSize() == fAsimovGlobObs.getSize() );
260 globObs.snapshot(globObsSnapshot);
261 globObs = fAsimovGlobObs;
268 RooRealVar * muAlt = (RooRealVar*) poiAlt.first();
271 oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator::Initialize Find best conditional NLL on ASIMOV data set for given alt POI ( " <<
272 muAlt->GetName() <<
" ) = " << muAlt->getVal() << std::endl;
274 fNLLAsimov = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiAlt );
280 globObs = globObsSnapshot;
283 if (prevBins > 0 && xobs) xobs->setBins(prevBins);
285 fIsInitialized =
true;
291 Double_t AsymptoticCalculator::EvaluateNLL(RooAbsPdf & pdf, RooAbsData& data,
const RooArgSet * condObs,
const RooArgSet * globObs,
const RooArgSet *poiSet) {
292 int verbose = fgPrintLevel;
295 RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
296 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
299 RooArgSet* allParams = pdf.getParameters(data);
300 RooStats::RemoveConstantParameters(allParams);
303 RooArgSet conditionalObs;
304 if (condObs) conditionalObs.add(*condObs);
306 if (globObs) globalObs.add(*globObs);
309 auto& config = GetGlobalRooStatsConfig();
310 RooAbsReal* nll = pdf.createNLL(data, RooFit::CloneData(kFALSE),RooFit::Constrain(*allParams),RooFit::ConditionalObservables(conditionalObs), RooFit::GlobalObservables(globalObs),
311 RooFit::Offset(config.useLikelihoodOffset));
313 RooArgSet* attachedSet = nll->getVariables();
316 RooArgSet paramsSetConstant;
318 if (poiSet && poiSet->getSize() > 0) {
319 RooRealVar * muTest = (RooRealVar*) (poiSet->first());
320 RooRealVar * poiVar =
dynamic_cast<RooRealVar*
>( attachedSet->find( muTest->GetName() ) );
321 if (poiVar && !poiVar->isConstant() ) {
322 poiVar->setVal( muTest->getVal() );
323 poiVar->setConstant();
324 paramsSetConstant.add(*poiVar);
326 if (poiSet->getSize() > 1)
327 std::cout <<
"Model with more than one POI are not supported - ignore extra parameters, consider only first one" << std::endl;
353 RooArgSet nllParams(*attachedSet);
354 RooStats::RemoveConstantParameters(&nllParams);
356 bool skipFit = (nllParams.getSize() == 0);
363 int minimPrintLevel = verbose;
365 RooMinimizer minim(*nll);
366 int strategy = ROOT::Math::MinimizerOptions::DefaultStrategy();
367 minim.setStrategy( strategy);
368 minim.setEvalErrorWall(config.useEvalErrorWall);
370 double tol = ROOT::Math::MinimizerOptions::DefaultTolerance();
371 tol = std::max(tol,1.0);
374 minim.setPrintLevel(minimPrintLevel-1);
376 minim.optimizeConst(2);
377 TString minimizer = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
378 TString algorithm = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo();
381 std::cout <<
"AsymptoticCalculator::EvaluateNLL ........ using " << minimizer <<
" / " << algorithm
382 <<
" with strategy " << strategy <<
" and tolerance " << tol << std::endl;
385 for (
int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
387 status = minim.minimize(minimizer, algorithm);
393 printf(
" ----> Doing a re-scan first\n");
394 minim.minimize(minimizer,
"Scan");
397 if (ROOT::Math::MinimizerOptions::DefaultStrategy() == 0 ) {
398 printf(
" ----> trying with strategy = 1\n");
399 minim.setStrategy(1);
405 printf(
" ----> trying with improve\n");
406 minimizer =
"Minuit";
407 algorithm =
"migradimproved";
412 RooFitResult * result = 0;
416 result = minim.save();
419 if (!RooStats::IsNLLOffset() )
420 val = result->minNll();
422 bool previous = RooAbsReal::hideOffset();
423 RooAbsReal::setHideOffset(kTRUE) ;
425 if (!previous) RooAbsReal::setHideOffset(kFALSE) ;
430 oocoutE((TObject*)0,Fitting) <<
"FIT FAILED !- return a NaN NLL " << std::endl;
431 val = TMath::QuietNaN();
434 minim.optimizeConst(
false);
435 if (result)
delete result;
442 std::cout <<
"AsymptoticCalculator::EvaluateNLL - value = " << val;
444 muTest = ( (RooRealVar*) poiSet->first() )->getVal();
445 std::cout <<
" for poi fixed at = " << muTest;
448 std::cout <<
"\tfit time : ";
452 std::cout << std::endl;
456 if (poiSet && paramsSetConstant.getSize() > 0) SetAllConstant(paramsSetConstant,
false);
459 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
476 HypoTestResult* AsymptoticCalculator::GetHypoTest()
const {
477 int verbose = fgPrintLevel;
480 if (!fIsInitialized) {
481 if (!Initialize() ) {
482 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::GetHypoTest - Error initializing Asymptotic calculator - return NULL result " << endl;
488 oocoutE((TObject*)0,InputArguments) <<
"AsymptoticCalculator::GetHypoTest - Asimov data set has not been generated - return NULL result " << endl;
492 assert(GetNullModel() );
495 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
500 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
501 assert(nullSnapshot && nullSnapshot->getSize() > 0);
505 RooArgSet poiTest(*nullSnapshot);
507 if (poiTest.getSize() > 1) {
508 oocoutW((TObject*)0,InputArguments) <<
"AsymptoticCalculator::GetHypoTest: snapshot has more than one POI - assume as POI first parameter " << std::endl;
511 RooArgSet * allParams = nullPdf->getParameters(*GetData() );
512 *allParams = fBestFitParams;
517 RooRealVar * muHat =
dynamic_cast<RooRealVar*
> ( fBestFitPoi.first() );
518 assert(muHat &&
"no best fit parameter defined");
519 RooRealVar * muTest =
dynamic_cast<RooRealVar*
> ( nullSnapshot->find(muHat->GetName() ) );
520 assert(muTest &&
"poi snapshot is not existing");
525 std::cout << std::endl;
526 oocoutI((TObject*)0,Eval) <<
"AsymptoticCalculator::GetHypoTest: - perform an hypothesis test for POI ( " << muTest->GetName() <<
" ) = " << muTest->getVal() << std::endl;
527 oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator::GetHypoTest - Find best conditional NLL on OBSERVED data set ..... " << std::endl;
531 double condNLL = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()), GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiTest);
533 double qmu = 2.*(condNLL - fNLLObs);
538 oocoutP((TObject*)0,Eval) <<
"\t OBSERVED DATA : qmu = " << qmu <<
" condNLL = " << condNLL <<
" uncond " << fNLLObs << std::endl;
542 double tol = 2.E-3 * std::max(1.,ROOT::Math::MinimizerOptions::DefaultTolerance());
543 if (qmu < -tol || TMath::IsNaN(fNLLObs) ) {
546 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: Found a negative value of the qmu - retry to do the unconditional fit "
549 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: unconditional fit failed before - retry to do it now "
553 double nll = EvaluateNLL( *nullPdf, const_cast<RooAbsData&>(*GetData()),GetNullModel()->GetConditionalObservables(),GetNullModel()->GetGlobalObservables());
555 if (nll < fNLLObs || (TMath::IsNaN(fNLLObs) && !TMath::IsNaN(nll) ) ) {
556 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: Found a better unconditional minimum "
557 <<
" old NLL = " << fNLLObs <<
" old muHat " << muHat->getVal() << std::endl;
561 const RooArgSet * poi = GetNullModel()->GetParametersOfInterest();
563 fBestFitPoi.removeAll();
564 poi->snapshot(fBestFitPoi);
566 muHat =
dynamic_cast<RooRealVar*
> ( fBestFitPoi.first() );
569 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: New minimum found for "
570 <<
" NLL = " << fNLLObs <<
" muHat " << muHat->getVal() << std::endl;
573 qmu = 2.*(condNLL - fNLLObs);
576 oocoutP((TObject*)0,Eval) <<
"After unconditional refit, new qmu value is " << qmu << std::endl;
582 oocoutE((TObject*)0,Minimization) <<
"AsymptoticCalculator: qmu is still < 0 for mu = "
583 << muTest->getVal() <<
" return a dummy result "
585 return new HypoTestResult();
587 if (TMath::IsNaN(qmu) ) {
588 oocoutE((TObject*)0,Minimization) <<
"AsymptoticCalculator: failure in fitting for qmu or qmuA "
589 << muTest->getVal() <<
" return a dummy result "
591 return new HypoTestResult();
604 RooArgSet globObsSnapshot;
605 if (GetNullModel()->GetGlobalObservables() ) {
606 globObs.add(*GetNullModel()->GetGlobalObservables());
608 globObs.snapshot(globObsSnapshot);
609 globObs = fAsimovGlobObs;
613 if (verbose > 0) oocoutP((TObject*)0,Eval) <<
"AsymptoticCalculator::GetHypoTest -- Find best conditional NLL on ASIMOV data set .... " << std::endl;
615 double condNLL_A = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables(), &poiTest);
618 double qmu_A = 2.*(condNLL_A - fNLLAsimov );
621 oocoutP((TObject*)0,Eval) <<
"\t ASIMOV data qmu_A = " << qmu_A <<
" condNLL = " << condNLL_A <<
" uncond " << fNLLAsimov << std::endl;
623 if (qmu_A < -tol || TMath::IsNaN(fNLLAsimov) ) {
626 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: Found a negative value of the qmu Asimov- retry to do the unconditional fit "
629 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: Fit failed for unconditional the qmu Asimov- retry unconditional fit "
633 double nll = EvaluateNLL( *nullPdf, *fAsimovData, GetNullModel()->GetConditionalObservables(), GetNullModel()->GetGlobalObservables() );
635 if (nll < fNLLAsimov || (TMath::IsNaN(fNLLAsimov) && !TMath::IsNaN(nll) )) {
636 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: Found a better unconditional minimum for Asimov data set"
637 <<
" old NLL = " << fNLLAsimov << std::endl;
642 oocoutW((TObject*)0,Minimization) <<
"AsymptoticCalculator: New minimum found for "
643 <<
" NLL = " << fNLLAsimov << std::endl;
644 qmu_A = 2.*(condNLL_A - fNLLAsimov);
647 oocoutP((TObject*)0,Eval) <<
"After unconditional Asimov refit, new qmu_A value is " << qmu_A << std::endl;
653 oocoutE((TObject*)0,Minimization) <<
"AsymptoticCalculator: qmu_A is still < 0 for mu = "
654 << muTest->getVal() <<
" return a dummy result "
656 return new HypoTestResult();
658 if (TMath::IsNaN(qmu) ) {
659 oocoutE((TObject*)0,Minimization) <<
"AsymptoticCalculator: failure in fitting for qmu or qmuA "
660 << muTest->getVal() <<
" return a dummy result "
662 return new HypoTestResult();
667 globObs = globObsSnapshot;
679 bool useQTilde =
false;
681 if (!fOneSidedDiscovery) {
682 if (fUseQTilde == -1 && !fOneSidedDiscovery) {
684 RooRealVar * muAlt =
dynamic_cast<RooRealVar*
>( GetAlternateModel()->GetSnapshot()->first() );
688 if (muTest->getMin() == muAlt->getVal() ) {
690 oocoutI((TObject*)0,InputArguments) <<
"Minimum of POI is " << muTest->getMin() <<
" corresponds to alt snapshot - using qtilde asymptotic formulae " << std::endl;
693 oocoutI((TObject*)0,InputArguments) <<
"Minimum of POI is " << muTest->getMin() <<
" is different to alt snapshot " << muAlt->getVal()
694 <<
" - using standard q asymptotic formulae " << std::endl;
697 useQTilde = fUseQTilde;
703 if ( muHat->getVal() > muTest->getVal() ) {
704 oocoutI((TObject*)0,Eval) <<
"Using one-sided qmu - setting qmu to zero muHat = " << muHat->getVal()
705 <<
" muTest = " << muTest->getVal() << std::endl;
709 if (fOneSidedDiscovery ) {
710 if ( muHat->getVal() < muTest->getVal() ) {
711 oocoutI((TObject*)0,Eval) <<
"Using one-sided discovery qmu - setting qmu to zero muHat = " << muHat->getVal()
712 <<
" muTest = " << muTest->getVal() << std::endl;
718 if (qmu < 0 && qmu > -tol) qmu = 0;
719 if (qmu_A < 0 && qmu_A > -tol) qmu_A = 0;
726 double pnull = -1, palt = -1;
731 double sqrtqmu = (qmu > 0) ? std::sqrt(qmu) : 0;
732 double sqrtqmu_A = (qmu_A > 0) ? std::sqrt(qmu_A) : 0;
735 if (fOneSided || fOneSidedDiscovery) {
739 oocoutI((TObject*)0,Eval) <<
"Using one-sided limit asymptotic formula (qmu)" << endl;
741 oocoutI((TObject*)0,Eval) <<
"Using one-sided discovery asymptotic formula (q0)" << endl;
743 pnull = ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
744 palt = ROOT::Math::normal_cdf( sqrtqmu_A - sqrtqmu, 1.);
748 if (verbose > 2) oocoutI((TObject*)0,Eval) <<
"Using two-sided asymptotic formula (tmu)" << endl;
749 pnull = 2.*ROOT::Math::normal_cdf_c( sqrtqmu, 1.);
750 palt = ROOT::Math::normal_cdf_c( sqrtqmu + sqrtqmu_A, 1.) +
751 ROOT::Math::normal_cdf_c( sqrtqmu - sqrtqmu_A, 1.);
758 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) {
759 if (verbose > 2) oocoutI((TObject*)0,Eval) <<
"Using qmu_tilde (qmu is greater than qmu_A)" << endl;
760 pnull = ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
761 palt = ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
767 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) {
768 if (verbose > 2) oocoutI((TObject*)0,Eval) <<
"Using tmu_tilde (qmu is greater than qmu_A)" << endl;
769 pnull = ROOT::Math::normal_cdf_c(sqrtqmu,1.) +
770 ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
771 palt = ROOT::Math::normal_cdf_c( sqrtqmu_A + sqrtqmu, 1.) +
772 ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
780 string resultname =
"HypoTestAsymptotic_result";
781 HypoTestResult* res =
new HypoTestResult(resultname.c_str(), pnull, palt);
785 oocoutP((TObject*)0,Eval)
786 <<
"poi = " << muTest->getVal() <<
" qmu = " << qmu <<
" qmu_A = " << qmu_A
787 <<
" sigma = " << muTest->getVal()/sqrtqmu_A
788 <<
" CLsplusb = " << pnull <<
" CLb = " << palt <<
" CLs = " << res->CLs() << std::endl;
794 struct PaltFunction {
795 PaltFunction(
double offset,
double pval,
int icase) :
796 fOffset(offset), fPval(pval), fCase(icase) {}
797 double operator() (
double x)
const {
798 return ROOT::Math::normal_cdf_c(x + fOffset) + ROOT::Math::normal_cdf_c(fCase*(x - fOffset)) - fPval;
808 double AsymptoticCalculator::GetExpectedPValues(
double pnull,
double palt,
double nsigma,
bool useCls,
bool oneSided ) {
810 double sqrtqmu = ROOT::Math::normal_quantile_c( pnull,1.);
811 double sqrtqmu_A = ROOT::Math::normal_quantile( palt,1.) + sqrtqmu;
812 double clsplusb = ROOT::Math::normal_cdf_c( sqrtqmu_A - nsigma, 1.);
813 if (!useCls)
return clsplusb;
814 double clb = ROOT::Math::normal_cdf( nsigma, 1.);
815 return (clb == 0) ? -1 : clsplusb / clb;
820 double sqrttmu = ROOT::Math::normal_quantile_c( 0.5*pnull,1.);
826 PaltFunction f( sqrttmu, palt, -1);
827 ROOT::Math::BrentRootFinder brf;
828 ROOT::Math::WrappedFunction<PaltFunction> wf(f);
829 brf.SetFunction( wf, 0, 20);
830 bool ret = brf.Solve();
832 oocoutE((TObject*)0,Eval) <<
"Error finding expected p-values - return -1" << std::endl;
835 double sqrttmu_A = brf.Root();
838 PaltFunction f2( sqrttmu_A, ROOT::Math::normal_cdf( nsigma, 1.), 1);
839 ROOT::Math::WrappedFunction<PaltFunction> wf2(f2);
840 brf.SetFunction(wf2,0,20);
843 oocoutE((TObject*)0,Eval) <<
"Error finding expected p-values - return -1" << std::endl;
846 return 2*ROOT::Math::normal_cdf_c( brf.Root(),1.);
857 void AsymptoticCalculator::FillBins(
const RooAbsPdf & pdf,
const RooArgList &obs, RooAbsData & data,
int &index,
double &binVolume,
int &ibin) {
859 bool debug = (fgPrintLevel >= 2);
861 RooRealVar * v =
dynamic_cast<RooRealVar*
>( &(obs[index]) );
864 RooArgSet obstmp(obs);
865 double expectedEvents = pdf.expectedEvents(obstmp);
870 if (debug) cout <<
"looping on observable " << v->GetName() << endl;
871 for (
int i = 0; i < v->getBins(); ++i) {
873 if (index < obs.getSize() -1) {
875 double prevBinVolume = binVolume;
876 binVolume *= v->getBinWidth(i);
877 FillBins(pdf, obs, data, index, binVolume, ibin);
879 binVolume = prevBinVolume;
884 double totBinVolume = binVolume * v->getBinWidth(i);
885 double fval = pdf.getVal(&obstmp)*totBinVolume;
888 if (fval*expectedEvents <= 0)
890 if (fval*expectedEvents < 0)
891 cout <<
"WARNING::Detected a bin with negative expected events! Please check your inputs." << endl;
893 cout <<
"WARNING::Detected a bin with zero expected events- skip it" << endl;
897 data.add(obs, fval*expectedEvents);
900 cout <<
"bin " << ibin <<
"\t";
901 for (
int j=0; j < obs.getSize(); ++j) { cout <<
" " << ((RooRealVar&) obs[j]).getVal(); }
902 cout <<
" w = " << fval*expectedEvents;
913 cout <<
"ending loop on .. " << v->GetName() << endl;
923 bool AsymptoticCalculator::SetObsToExpected(RooProdPdf &prod,
const RooArgSet &obs)
925 RooLinkedListIter iter(prod.pdfList().iterator());
927 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
928 if (!a->dependsOn(obs))
continue;
929 RooPoisson *pois = 0;
930 RooGaussian * gaus = 0;
931 if ((pois = dynamic_cast<RooPoisson *>(a)) != 0) {
932 ret &= SetObsToExpected(*pois, obs);
933 pois->setNoRounding(
true);
934 }
else if ((gaus = dynamic_cast<RooGaussian *>(a)) != 0) {
935 ret &= SetObsToExpected(*gaus, obs);
938 RooProdPdf *subprod =
dynamic_cast<RooProdPdf *
>(a);
940 ret &= SetObsToExpected(*subprod, obs);
942 oocoutE((TObject*)0,InputArguments) <<
"Illegal term in counting model: "
943 <<
"the PDF " << a->GetName()
944 <<
" depends on the observables, but is not a Poisson, Gaussian or Product"
962 bool AsymptoticCalculator::SetObsToExpected(RooAbsPdf &pdf,
const RooArgSet &obs)
964 RooRealVar *myobs = 0;
965 RooAbsReal *myexp = 0;
966 const char * pdfName = pdf.IsA()->GetName();
967 RooFIter iter(pdf.serverMIterator());
968 for (RooAbsArg *a = iter.next(); a != 0; a = iter.next()) {
969 if (obs.contains(*a)) {
971 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : Has two observables ?? " << endl;
974 myobs =
dynamic_cast<RooRealVar *
>(a);
976 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : Observable is not a RooRealVar??" << endl;
980 if (!a->isConstant() ) {
982 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : Has two non-const arguments " << endl;
985 myexp =
dynamic_cast<RooAbsReal *
>(a);
987 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : Expected is not a RooAbsReal??" << endl;
994 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : No observable?" << endl;
998 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::SetObsExpected( " << pdfName <<
" ) : No observable?" << endl;
1002 myobs->setVal(myexp->getVal());
1004 if (fgPrintLevel > 2) {
1005 std::cout <<
"SetObsToExpected : setting " << myobs->GetName() <<
" to expected value " << myexp->getVal() <<
" of " << myexp->GetName() << std::endl;
1016 RooAbsData * AsymptoticCalculator::GenerateCountingAsimovData(RooAbsPdf & pdf,
const RooArgSet & observables,
const RooRealVar & , RooCategory * channelCat) {
1017 RooArgSet obs(observables);
1018 RooProdPdf *prod =
dynamic_cast<RooProdPdf *
>(&pdf);
1019 RooPoisson *pois = 0;
1020 RooGaussian *gaus = 0;
1022 if (fgPrintLevel > 1)
1023 std::cout <<
"generate counting Asimov data for pdf of type " << pdf.IsA()->GetName() << std::endl;
1027 r = SetObsToExpected(*prod, observables);
1028 }
else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != 0) {
1029 r = SetObsToExpected(*pois, observables);
1031 pois->setNoRounding(
true);
1032 }
else if ((gaus = dynamic_cast<RooGaussian *>(&pdf)) != 0) {
1033 r = SetObsToExpected(*gaus, observables);
1035 oocoutE((TObject*)0,InputArguments) <<
"A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << endl;
1040 icat = channelCat->getIndex();
1043 RooDataSet *ret =
new RooDataSet(TString::Format(
"CountingAsimovData%d",icat),TString::Format(
"CountingAsimovData%d",icat), obs);
1054 RooAbsData * AsymptoticCalculator::GenerateAsimovDataSinglePdf(
const RooAbsPdf & pdf,
const RooArgSet & allobs,
const RooRealVar & weightVar, RooCategory * channelCat) {
1056 int printLevel = fgPrintLevel;
1059 std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );
1063 if (!pdf.canBeExtended() )
return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);
1065 RooArgSet obsAndWeight(*obs);
1066 obsAndWeight.add(weightVar);
1068 RooDataSet* asimovData = 0;
1070 int icat = channelCat->getIndex();
1071 asimovData =
new RooDataSet(TString::Format(
"AsimovData%d",icat),TString::Format(
"combAsimovData%d",icat),
1072 RooArgSet(obsAndWeight,*channelCat),RooFit::WeightVar(weightVar));
1075 asimovData =
new RooDataSet(
"AsimovData",
"AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));
1080 RooArgList obsList(*obs);
1083 if (printLevel >= 2) {
1084 cout <<
"Generating Asimov data for pdf " << pdf.GetName() << endl;
1085 cout <<
"list of observables " << endl;
1090 double binVolume = 1;
1092 FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
1093 if (printLevel >= 2)
1094 cout <<
"filled from " << pdf.GetName() <<
" " << nbins <<
" nbins " <<
" volume is " << binVolume << endl;
1112 if (printLevel >= 1)
1114 asimovData->Print();
1117 if( TMath::IsNaN(asimovData->sumEntries()) ){
1118 cout <<
"sum entries is nan"<<endl;
1132 RooAbsData * AsymptoticCalculator::GenerateAsimovData(
const RooAbsPdf & pdf,
const RooArgSet & observables ) {
1134 int printLevel = fgPrintLevel;
1136 unique_ptr<RooRealVar> weightVar (
new RooRealVar(
"binWeightAsimov",
"binWeightAsimov", 1, 0, 1.E30 ));
1138 if (printLevel > 1) cout <<
" Generate Asimov data for observables"<<endl;
1140 const RooSimultaneous* simPdf =
dynamic_cast<const RooSimultaneous*
>(&pdf);
1143 return GenerateAsimovDataSinglePdf( pdf, observables, *weightVar, 0);
1146 std::map<std::string, RooDataSet*> asimovDataMap;
1149 RooCategory& channelCat = (RooCategory&)simPdf->indexCat();
1150 int nrIndices = channelCat.numTypes();
1151 if( nrIndices == 0 ) {
1152 oocoutW((TObject*)0,Generation) <<
"Simultaneous pdf does not contain any categories." << endl;
1154 for (
int i=0;i<nrIndices;i++){
1155 channelCat.setIndex(i);
1158 RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getLabel()) ;
1159 assert(pdftmp != 0);
1163 cout <<
"on type " << channelCat.getLabel() <<
" " << channelCat.getIndex() << endl;
1166 RooAbsData * dataSinglePdf = GenerateAsimovDataSinglePdf( *pdftmp, observables, *weightVar, &channelCat);
1169 if (!dataSinglePdf) {
1170 oocoutE((TObject*)0,Generation) <<
"Error generating an Asimov data set for pdf " << pdftmp->GetName() << endl;
1174 if (asimovDataMap.count(
string(channelCat.getLabel())) != 0) {
1175 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::GenerateAsimovData(): The PDF for " << channelCat.getLabel()
1176 <<
" was already defined. It will be overridden. The faulty category definitions follow:" << endl;
1177 channelCat.Print(
"V");
1180 asimovDataMap[string(channelCat.getLabel())] = (RooDataSet*) dataSinglePdf;
1184 cout <<
"channel: " << channelCat.getLabel() <<
", data: ";
1185 dataSinglePdf->Print();
1190 RooArgSet obsAndWeight(observables);
1191 obsAndWeight.add(*weightVar);
1194 RooDataSet* asimovData =
new RooDataSet(
"asimovDataFullModel",
"asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
1195 RooFit::Index(channelCat),RooFit::Import(asimovDataMap),RooFit::WeightVar(*weightVar));
1197 for (
auto element : asimovDataMap) {
1198 delete element.second;
1216 RooAbsData * AsymptoticCalculator::MakeAsimovData(RooAbsData & realData,
const ModelConfig & model,
const RooArgSet & paramValues, RooArgSet & asimovGlobObs,
const RooArgSet * genPoiValues ) {
1218 int verbose = fgPrintLevel;
1221 RooArgSet poi(*model.GetParametersOfInterest());
1226 RooLinkedListIter it = poi.iterator();
1227 RooRealVar* tmpPar = NULL;
1228 RooArgSet paramsSetConstant;
1229 while((tmpPar = (RooRealVar*)it.Next())){
1230 tmpPar->setConstant();
1232 std::cout <<
"MakeAsimov: Setting poi " << tmpPar->GetName() <<
" to a constant value = " << tmpPar->getVal() << std::endl;
1233 paramsSetConstant.add(*tmpPar);
1237 bool hasFloatParams =
false;
1238 RooArgSet constrainParams;
1239 if (model.GetNuisanceParameters()) {
1240 constrainParams.add(*model.GetNuisanceParameters());
1241 RooStats::RemoveConstantParameters(&constrainParams);
1242 if (constrainParams.getSize() > 0) hasFloatParams =
true;
1246 std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
1247 RooLinkedListIter iter(params->iterator());
1248 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1249 RooRealVar *rrv =
dynamic_cast<RooRealVar *
>(a);
1250 if ( rrv != 0 && rrv->isConstant() == false ) { hasFloatParams =
true;
break; }
1253 if (hasFloatParams) {
1256 TStopwatch tw2; tw2.Start();
1257 int minimPrintLevel = ROOT::Math::MinimizerOptions::DefaultPrintLevel();
1259 std::cout <<
"MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
1260 minimPrintLevel = verbose;
1262 std::cout <<
"POI values:\n"; poi.Print(
"v");
1264 std::cout <<
"Nuis param values:\n";
1265 constrainParams.Print(
"v");
1269 RooFit::MsgLevel msglevel = RooMsgService::instance().globalKillBelow();
1270 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
1272 RooArgSet conditionalObs;
1273 if (model.GetConditionalObservables()) conditionalObs.add(*model.GetConditionalObservables());
1274 RooArgSet globalObs;
1275 if (model.GetGlobalObservables()) globalObs.add(*model.GetGlobalObservables());
1277 std::string minimizerType = ROOT::Math::MinimizerOptions::DefaultMinimizerType();
1278 std::string minimizerAlgo = ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo();
1279 std::vector<RooCmdArg> args;
1280 args.push_back(RooFit::Minimizer(minimizerType.c_str(),minimizerAlgo.c_str()));
1281 args.push_back(RooFit::Strategy(ROOT::Math::MinimizerOptions::DefaultStrategy()));
1282 args.push_back(RooFit::PrintLevel(minimPrintLevel-1));
1283 args.push_back(RooFit::Hesse(
false));
1284 args.push_back(RooFit::Constrain(constrainParams));
1285 args.push_back(RooFit::GlobalObservables(globalObs));
1286 args.push_back(RooFit::ConditionalObservables(conditionalObs));
1287 args.push_back(RooFit::Offset(GetGlobalRooStatsConfig().useLikelihoodOffset));
1288 args.push_back(RooFit::EvalErrorWall(GetGlobalRooStatsConfig().useEvalErrorWall));
1290 RooLinkedList argList;
1291 for (
auto& arg : args) {
1294 model.GetPdf()->fitTo(realData, argList);
1295 if (verbose>0) { std::cout <<
"fit time "; tw2.Print();}
1298 if (model.GetNuisanceParameters() ) {
1299 std::cout <<
"Nuisance parameters after fit for asimov dataset: " << std::endl;
1300 model.GetNuisanceParameters()->Print(
"V");
1304 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
1309 SetAllConstant(paramsSetConstant,
false);
1312 RooArgSet * allParams = model.GetPdf()->getParameters(realData);
1313 RooStats::RemoveConstantParameters( allParams );
1316 if (genPoiValues) *allParams = *genPoiValues;
1320 RooAbsData * asimovData = MakeAsimovData(model, *allParams, asimovGlobObs);
1337 RooAbsData * AsymptoticCalculator::MakeAsimovData(
const ModelConfig & model,
const RooArgSet & allParamValues, RooArgSet & asimovGlobObs) {
1339 int verbose = fgPrintLevel;
1346 if (allParamValues.getSize() > 0) {
1347 RooArgSet * allVars = model.GetPdf()->getVariables();
1348 *allVars = allParamValues;
1354 RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
1357 std::cout <<
"Generated Asimov data for observables "; (model.GetObservables() )->Print();
1359 if (asimov->numEntries() == 1 ) {
1360 std::cout <<
"--- Asimov data values \n";
1361 asimov->get()->Print(
"v");
1364 std::cout <<
"--- Asimov data numEntries = " << asimov->numEntries() <<
" sumOfEntries = " << asimov->sumEntries() << std::endl;
1366 std::cout <<
"\ttime for generating : "; tw.Print();
1383 if (model.GetGlobalObservables() && model.GetGlobalObservables()->getSize() > 0) {
1386 std::cout <<
"Generating Asimov data for global observables " << std::endl;
1389 RooArgSet gobs(*model.GetGlobalObservables());
1392 RooArgSet snapGlobalObsData;
1393 SetAllConstant(gobs,
true);
1394 gobs.snapshot(snapGlobalObsData);
1398 if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
1399 if (nuis.getSize() == 0) {
1400 oocoutW((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
1401 <<
" set global observables to model values " << endl;
1402 asimovGlobObs = gobs;
1407 std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,
"TempNuisPdf") );
1408 if (nuispdf.get() == 0) {
1409 oocoutF((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and global obs but no nuisance pdf "
1414 RooProdPdf *prod =
dynamic_cast<RooProdPdf *
>(nuispdf.get());
1416 pdfList.add(prod->pdfList());
1420 pdfList.add(*nuispdf.get());
1422 RooLinkedListIter iter(pdfList.iterator());
1423 for (RooAbsArg *a = (RooAbsArg *) iter.Next(); a != 0; a = (RooAbsArg *) iter.Next()) {
1424 RooAbsPdf *cterm =
dynamic_cast<RooAbsPdf *
>(a);
1425 assert(cterm &&
"AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
1426 if (!cterm->dependsOn(nuis))
continue;
1428 if (
typeid(*cterm) ==
typeid(RooUniform))
continue;
1430 std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
1431 std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
1432 if (cgobs->getSize() > 1) {
1433 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1434 <<
" has multiple global observables -cannot generate - skip it" << std::endl;
1437 else if (cgobs->getSize() == 0) {
1438 oocoutW((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1439 <<
" has no global observables - skip it" << std::endl;
1443 RooRealVar &rrv =
dynamic_cast<RooRealVar &
>(*cgobs->first());
1446 RooStats::RemoveConstantParameters(cpars.get());
1447 if (cpars->getSize() != 1) {
1448 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData:constraint term "
1449 << cterm->GetName() <<
" has multiple floating params - cannot generate - skip it " << std::endl;
1453 bool foundServer =
false;
1456 TClass * cClass = cterm->IsA();
1457 if (verbose > 2) std::cout <<
"Constraint " << cterm->GetName() <<
" of type " << cClass->GetName() << std::endl;
1458 if ( cClass != RooGaussian::Class() && cClass != RooPoisson::Class() &&
1459 cClass != RooGamma::Class() && cClass != RooLognormal::Class() &&
1460 cClass != RooBifurGauss::Class() ) {
1461 TString className = (cClass) ? cClass->GetName() :
"undefined";
1462 oocoutW((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData:constraint term "
1463 << cterm->GetName() <<
" of type " << className
1464 <<
" is a non-supported type - result might be not correct " << std::endl;
1468 if (cClass == RooPoisson::Class() ) {
1469 RooPoisson * pois =
dynamic_cast<RooPoisson*
>(cterm);
1471 pois->setNoRounding(
true);
1475 RooAbsArg * arg = cterm->findServer(rrv);
1480 if ( cClass != RooGamma::Class() ) {
1481 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData:constraint term "
1482 << cterm->GetName() <<
" has no direct dependence on global observable- cannot generate it " << std::endl;
1491 RooAbsReal * thetaGamma = 0;
1492 if ( cClass == RooGamma::Class() ) {
1493 RooFIter itc(cterm->serverMIterator() );
1494 for (RooAbsArg *a2 = itc.next(); a2 != 0; a2 = itc.next()) {
1495 if (TString(a2->GetName()).Contains(
"theta") ) {
1496 thetaGamma =
dynamic_cast<RooAbsReal*
>(a2);
1500 if (thetaGamma == 0) {
1501 oocoutI((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData:constraint term "
1502 << cterm->GetName() <<
" is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is 1 " << std::endl;
1506 std::cout <<
"Gamma constraint has a scale " << thetaGamma->GetName() <<
" = " << thetaGamma->getVal() << std::endl;
1509 RooFIter iter2(cterm->serverMIterator() );
1510 for (RooAbsArg *a2 = iter2.next(); a2 != 0; a2 = iter2.next()) {
1511 RooAbsReal * rrv2 =
dynamic_cast<RooAbsReal *
>(a2);
1512 if (verbose > 2) std::cout <<
"Loop on constraint server term " << a2->GetName() << std::endl;
1513 if (rrv2 && rrv2->dependsOn(nuis) ) {
1518 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData:constraint term "
1519 << cterm->GetName() <<
" constraint term has more server depending on nuisance- cannot generate it " <<
1521 foundServer =
false;
1524 if (thetaGamma && thetaGamma->getVal() > 0)
1525 rrv.setVal( rrv2->getVal() / thetaGamma->getVal() );
1527 rrv.setVal( rrv2->getVal() );
1531 std::cout <<
"setting global observable "<< rrv.GetName() <<
" to value " << rrv.getVal()
1532 <<
" which comes from " << rrv2->GetName() << std::endl;
1537 oocoutE((TObject*)0,Generation) <<
"AsymptoticCalculator::MakeAsimovData - can't find nuisance for constraint term - global observables will not be set to Asimov value " << cterm->GetName() << std::endl;
1538 std::cerr <<
"Parameters: " << std::endl;
1540 std::cerr <<
"Observables: " << std::endl;
1549 asimovGlobObs.removeAll();
1550 SetAllConstant(gobs,
true);
1551 gobs.snapshot(asimovGlobObs);
1554 gobs = snapGlobalObsData;
1557 std::cout <<
"Generated Asimov data for global observables ";
1558 if (verbose == 1) gobs.Print();
1562 std::cout <<
"\nGlobal observables for data: " << std::endl;
1564 std::cout <<
"\nGlobal observables for asimov: " << std::endl;
1565 asimovGlobObs.Print(
"V");