Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
AsymptoticCalculator.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, Sven Kreiss 23/05/10
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11/** \class RooStats::AsymptoticCalculator
12 \ingroup Roostats
13
14Hypothesis Test Calculator based on the asymptotic formulae for the profile
15likelihood ratio.
16
17It performs hypothesis tests using the asymptotic formula for the profile likelihood, and
18uses the Asimov data set to compute expected significances or limits.
19
20See G. Cowan, K. Cranmer, E. Gross and O. Vitells: Asymptotic formulae for
21likelihood- based tests of new physics. Eur. Phys. J., C71:1–19, 2011.
22It provides methods to perform hypothesis tests using the likelihood function,
23and computes the \f$p\f$-values for the null and the alternate hypothesis using the asymptotic
24formulae for the profile likelihood ratio described in the given paper.
25
26The calculator provides methods to produce the Asimov dataset, *i.e.* a dataset
27generated where the observed values are equal to the expected ones.
28The Asimov data set is then used to compute the observed asymptotic \f$p\f$-value for
29the alternate hypothesis and the asymptotic expected \f$p\f$-values.
30
31The asymptotic formulae are valid only for one POI (parameter of interest). So
32the calculator works only for one-dimensional (one POI) models.
33If more than one POI exists, only the first one is used.
34
35The calculator can generate Asimov datasets from two kinds of PDFs:
36- "Counting" distributions: RooPoisson, RooGaussian, or products of RooPoissons.
37- Extended, *i.e.* number of events can be read off from extended likelihood term.
38*/
39
40
45
46#include "RooArgSet.h"
47#include "RooArgList.h"
48#include "RooProdPdf.h"
49#include "RooSimultaneous.h"
50#include "RooDataSet.h"
51#include "RooCategory.h"
52#include "RooRealVar.h"
53#include "RooMinimizer.h"
54#include "RooFitResult.h"
56#include "RooPoisson.h"
57#include "RooUniform.h"
58#include "RooGamma.h"
59#include "RooGaussian.h"
60#include "RooMultiVarGaussian.h"
61#include "RooBifurGauss.h"
62#include "RooLognormal.h"
63#include "RooDataHist.h"
64#include <cmath>
65#include <typeinfo>
66
69
70#include "TStopwatch.h"
71
72#include <ROOT/RSpan.hxx>
73
74using namespace RooStats;
75using std::string, std::unique_ptr;
76
77
78namespace {
79
80/// Control print level (0 minimal, 1 normal, 2 debug).
81int &fgPrintLevel()
82{
83
84 static int val = 1;
85 return val;
86}
87
88// Forward declaration.
90
91} // namespace
92
93////////////////////////////////////////////////////////////////////////////////
94/// set print level (static function)
95///
96/// - 0 minimal,
97/// - 1 normal,
98/// - 2 debug
99
101 fgPrintLevel() = level;
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// constructor for asymptotic calculator from Data set and ModelConfig
106
109 const ModelConfig &altModel,
110 const ModelConfig &nullModel, bool nominalAsimov) :
112 fOneSided(false), fOneSidedDiscovery(false), fNominalAsimov(nominalAsimov),
113 fUseQTilde(-1),
114 fNLLObs(0), fNLLAsimov(0),
115 fAsimovData(nullptr)
116{
117 if (!Initialize()) return;
118
119 int verbose = fgPrintLevel();
120 // try to guess default configuration
121 // (this part should be only in constructor because the null snapshot might change during HypoTestInversion
122 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
124 RooRealVar * muNull = dynamic_cast<RooRealVar*>(nullSnapshot->first() );
125 assert(muNull);
126 if (muNull->getVal() == muNull->getMin()) {
127 fOneSidedDiscovery = true;
128 if (verbose > 0)
129 oocoutI(nullptr,InputArguments) << "AsymptotiCalculator: Minimum of POI is " << muNull->getMin() << " corresponds to null snapshot - default configuration is one-sided discovery formulae " << std::endl;
130 }
131
132}
133
134////////////////////////////////////////////////////////////////////////////////
135/// Initialize the calculator
136/// The initialization will perform a global fit of the model to the data
137/// and build an Asimov data set.
138/// It will then also fit the model to the Asimov data set to find the likelihood value
139/// of the Asimov data set
140/// nominalAsimov is an option for using Asimov data set obtained using nominal nuisance parameter values
141/// By default the nuisance parameters are fitted to the data
142/// NOTE: If a fit has been done before, one for speeding up could set all the initial parameters
143/// to the fit value and in addition set the null snapshot to the best fit
144
146
147 int verbose = fgPrintLevel();
148 if (verbose >= 0)
149 oocoutP(nullptr,Eval) << "AsymptoticCalculator::Initialize...." << std::endl;
150
151
152 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
153 if (!nullPdf) {
154 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has not a pdf defined" << std::endl;
155 return false;
156 }
157 RooAbsData * obsData = const_cast<RooAbsData *>(GetData() );
158 if (!obsData ) {
159 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::Initialize - data set has not been defined" << std::endl;
160 return false;
161 }
163
164
165
166 const RooArgSet * poi = GetNullModel()->GetParametersOfInterest();
167 if (!poi || poi->empty()) {
168 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has not POI defined." << std::endl;
169 return false;
170 }
171 if (poi->size() > 1) {
172 oocoutW(nullptr,InputArguments) << "AsymptoticCalculator::Initialize - ModelConfig has more than one POI defined \n\t"
173 << "The asymptotic calculator works for only one POI - consider as POI only the first parameter"
174 << std::endl;
175 }
176
177
178 // This will set the poi value to the null snapshot value in the ModelConfig
179 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
180 if(nullSnapshot == nullptr || nullSnapshot->empty()) {
181 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::Initialize - Null model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << std::endl;
182 return false;
183 }
184
185 // GetNullModel()->Print();
186 // std::cout << "ASymptotic calc: null snapshot\n";
187 // nullSnapshot->Print("v");
188 // std::cout << "PDF variables ";
189 // nullPdf->getVariables()->Print("v");
190
191 // keep snapshot for the initial parameter values (need for nominal Asimov)
193 std::unique_ptr<RooArgSet> allParams{nullPdf->getParameters(data)};
195 if (fNominalAsimov) {
196 allParams->snapshot(nominalParams);
197 }
201
202 // evaluate the unconditional nll for the full model on the observed data
203 if (verbose >= 0)
204 oocoutP(nullptr,Eval) << "AsymptoticCalculator::Initialize - Find best unconditional NLL on observed data" << std::endl;
206 // fill also snapshot of best poi
207 poi->snapshot(fBestFitPoi);
208 RooRealVar * muBest = dynamic_cast<RooRealVar*>(fBestFitPoi.first());
209 assert(muBest);
210 if (verbose >= 0)
211 oocoutP(nullptr,Eval) << "Best fitted POI value = " << muBest->getVal() << " +/- " << muBest->getError() << std::endl;
212 // keep snapshot of all best fit parameters
213 allParams->snapshot(fBestFitParams);
214
215 // compute Asimov data set for the background (alt poi ) value
216 const RooArgSet * altSnapshot = GetAlternateModel()->GetSnapshot();
217 if(altSnapshot == nullptr || altSnapshot->empty()) {
218 oocoutE(nullptr,InputArguments) << "Alt (Background) model needs a snapshot. Set using modelconfig->SetSnapshot(poi)." << std::endl;
219 return false;
220 }
221
222 RooArgSet poiAlt(*altSnapshot); // this is the poi snapshot of B (i.e. for mu=0)
223
224 oocoutP(nullptr,Eval) << "AsymptoticCalculator: Building Asimov data Set" << std::endl;
225
226 // check that in case of binned models the n number of bins of the observables are consistent
227 // with the number of bins in the observed data
228 // This number will be used for making the Asimov data set so it will be more consistent with the
229 // observed data
230 int prevBins = 0;
231 RooRealVar * xobs = nullptr;
232 if (GetNullModel()->GetObservables() && GetNullModel()->GetObservables()->size() == 1 ) {
233 xobs = static_cast<RooRealVar*>((GetNullModel()->GetObservables())->first());
234 if (data.IsA() == RooDataHist::Class() ) {
235 if (data.numEntries() != xobs->getBins() ) {
236 prevBins = xobs->getBins();
237 oocoutW(nullptr,InputArguments) << "AsymptoticCalculator: number of bins in " << xobs->GetName() << " are different than data bins "
238 << " set the same data bins " << data.numEntries() << " in range "
239 << " [ " << xobs->getMin() << " , " << xobs->getMax() << " ]" << std::endl;
240 xobs->setBins(data.numEntries());
241 }
242 }
243 }
244
245 if (!fNominalAsimov) {
246 if (verbose >= 0)
247 oocoutI(nullptr,InputArguments) << "AsymptoticCalculator: Asimov data will be generated using fitted nuisance parameter values" << std::endl;
248 RooArgSet * tmp = (RooArgSet*) poiAlt.snapshot();
250 }
251
252 else {
253 // assume use current value of nuisance as nominal ones
254 if (verbose >= 0)
255 oocoutI(nullptr,InputArguments) << "AsymptoticCalculator: Asimovdata set will be generated using nominal (current) nuisance parameter values" << std::endl;
256 nominalParams.assign(poiAlt); // set poi to alt value but keep nuisance at the nominal one
258 }
259
260 if (!fAsimovData) {
261 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator: Error : Asimov data set could not be generated " << std::endl;
262 return false;
263 }
264
265 // set global observables to their Asimov values
268 if (GetNullModel()->GetGlobalObservables() ) {
269 globObs.add(*GetNullModel()->GetGlobalObservables());
270 assert(globObs.size() == fAsimovGlobObs.size() );
271 // store previous snapshot value
272 globObs.snapshot(globObsSnapshot);
273 globObs.assign(fAsimovGlobObs);
274 }
275
276
277 // evaluate the likelihood. Since we use on Asimov data , conditional and unconditional values should be the same
278 // do conditional fit since is faster
279
280 RooRealVar * muAlt = static_cast<RooRealVar*>(poiAlt.first());
281 assert(muAlt);
282 if (verbose >= 0) {
283 oocoutP(nullptr, Eval)
284 << "AsymptoticCalculator::Initialize Find best conditional NLL on ASIMOV data set for given alt POI ( "
285 << muAlt->GetName() << " ) = " << muAlt->getVal() << std::endl;
286 }
287
289 // for unconditional fit
290 //fNLLAsimov = EvaluateNLL( *nullPdf, *fAsimovData);
291 //poi->Print("v");
292
293 // restore previous value
294 globObs.assign(globObsSnapshot);
295
296 // restore number of bins
297 if (prevBins > 0 && xobs) xobs->setBins(prevBins);
298
299 fIsInitialized = true;
300 return true;
301}
302
303namespace {
304
306{
307 int verbose = fgPrintLevel();
308
309 RooAbsPdf &pdf = *modelConfig.GetPdf();
310
312 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
313
314
315 std::unique_ptr<RooArgSet> allParams{pdf.getParameters(data)};
317 // add constraint terms for all non-constant parameters
318
319 // need to call constrain for RooSimultaneous until stripDisconnected problem fixed
320 auto& config = GetGlobalRooStatsConfig();
321 std::unique_ptr<RooAbsReal> nll{modelConfig.createNLL(data, RooFit::Constrain(*allParams), RooFit::Offset(config.useLikelihoodOffset))};
322
323 std::unique_ptr<RooArgSet> attachedSet{nll->getVariables()};
324
325 // if poi are specified - do a conditional fit
327 // support now only one POI
328 if (poiSet && !poiSet->empty()) {
329 RooRealVar * muTest = static_cast<RooRealVar*> (poiSet->first());
330 RooRealVar * poiVar = dynamic_cast<RooRealVar*>(attachedSet->find( muTest->GetName() ) );
331 if (poiVar && !poiVar->isConstant() ) {
332 poiVar->setVal( muTest->getVal() );
333 poiVar->setConstant();
335 }
336 if (poiSet->size() > 1)
337 oocoutW(nullptr,InputArguments) << "Model with more than one POI are not supported - ignore extra parameters, consider only first one" << std::endl;
338
339
340
341 // This for more than one POI (not yet supported)
342 //
343 // RooLinkedListIter it = poiSet->iterator();
344 // RooRealVar* tmpPar = nullptr, *tmpParA=nullptr;
345 // while((tmpPar = (RooRealVar*)it.Next())){
346 // tmpParA = ((RooRealVar*)attachedSet->find(tmpPar->GetName()));
347 // tmpParA->setVal( tmpPar->getVal() );
348 // if (!tmpParA->isConstant() ) {
349 // tmpParA->setConstant();
350 // paramsSetConstant.add(*tmpParA);
351 // }
352 // }
353
354 // check if there are non-const parameters so it is worth to do the minimization
355
356 }
357
359 tw.Start();
360 double val = -1;
361
362 //check if needed to skip the fit
365 bool skipFit = (nllParams.empty());
366
367 if (skipFit) {
368 val = nll->getVal(); // just evaluate nll in conditional fits with model without nuisance params
369 } else {
370
372
373 RooMinimizer minim(*nll);
375 minim.setStrategy( strategy);
376 minim.setEvalErrorWall(config.useEvalErrorWall);
377 // use tolerance - but never smaller than 1 (default in RooMinimizer)
379 tol = std::max(tol,1.0); // 1.0 is the minimum value used in RooMinimizer
380 minim.setEps( tol );
381 //LM: RooMinimizer.setPrintLevel has +1 offset - so subtract here -1
382 minim.setPrintLevel(minimPrintLevel-1);
383 int status = -1;
384 minim.optimizeConst(2);
385 TString minimizer = ""; // empty string to take RooMinimizer default initially
387
388 if (verbose > 0) {
389 oocoutP(nullptr,Eval) << "AsymptoticCalculator::EvaluateNLL ........ using " << minimizer << " / " << algorithm
390 << " with strategy " << strategy << " and tolerance " << tol << std::endl;
391 }
392
393 for (int tries = 1, maxtries = 4; tries <= maxtries; ++tries) {
394 // status = minim.minimize(fMinimizer, ROOT::Math::MinimizerOptions::DefaultMinimizerAlgo().c_str());
395 status = minim.minimize(minimizer, algorithm);
396 // RooMinimizer::minimize returns -1 when the fit fails
397 if (status >= 0) {
398 break;
399 } else {
400 if (tries == 1) {
401 oocoutW(nullptr,Minimization) << " ----> Doing a re-scan first" << std::endl;
402 minim.minimize(minimizer,"Scan");
403 }
404 if (tries == 2) {
406 oocoutW(nullptr,Minimization) << " ----> trying with strategy = 1" << std::endl;
407 minim.setStrategy(1);
408 }
409 else
410 tries++; // skip this trial if strategy is already 1
411 }
412 if (tries == 3) {
413 oocoutW(nullptr,Minimization) << " ----> trying with improve" << std::endl;
414 minimizer = "Minuit";
415 algorithm = "migradimproved";
416 }
417 }
418 }
419
420 std::unique_ptr<RooFitResult> result;
421
422 // ignore errors in Hesse or in Improve and also when matrix was made pos def (status returned = 1)
423 if (status >= 0) {
424 result = std::unique_ptr<RooFitResult>{minim.save()};
425 }
426 if (result){
427 if (RooStats::NLLOffsetMode() != "initial") {
428 val = result->minNll();
429 } else {
432 val = nll->getVal();
434 }
435
436 }
437 else {
438 oocoutE(nullptr,Fitting) << "FIT FAILED !- return a NaN NLL " << std::endl;
439 val = TMath::QuietNaN();
440 }
441
442 minim.optimizeConst(false);
443 }
444
445 double muTest = 0;
446 if (verbose > 0) {
447 oocoutP(nullptr,Eval) << "AsymptoticCalculator::EvaluateNLL - value = " << val;
448 if (poiSet) {
449 muTest = ( static_cast<RooRealVar*>(poiSet->first()) )->getVal();
450 ooccoutP(nullptr,Eval) << " for poi fixed at = " << muTest;
451 }
452 if (!skipFit) {
453 tw.Stop();
454 ooccoutP(nullptr,Eval) << "\tfit time : " << tw.RealTime() << " s (real) " << tw.CpuTime() << " s (cpu)" << std::endl;
455 } else {
456 ooccoutP(nullptr,Eval) << std::endl;
457 }
458 }
459
460 // reset the parameter free which where set as constant
462
463
464 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
465
466 return val;
467}
468
469} // namespace
470
471////////////////////////////////////////////////////////////////////////////////
472/// It performs an hypothesis tests using the likelihood function
473/// and computes the p values for the null and the alternate using the asymptotic
474/// formulae for the profile likelihood ratio.
475/// See G. Cowan, K. Cranmer, E. Gross and O. Vitells.
476/// Asymptotic formulae for likelihood- based tests of new physics. Eur. Phys. J., C71:1–19, 2011.
477/// The formulae are valid only for one POI. If more than one POI exists consider as POI only the
478/// first one
479
481 int verbose = fgPrintLevel();
482
483 // re-initialized the calculator in case it is needed (pdf or data modified)
484 if (!fIsInitialized) {
485 if (!Initialize() ) {
486 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::GetHypoTest - Error initializing Asymptotic calculator - return nullptr result " << std::endl;
487 return nullptr;
488 }
489 }
490
491 if (!fAsimovData) {
492 oocoutE(nullptr,InputArguments) << "AsymptoticCalculator::GetHypoTest - Asimov data set has not been generated - return nullptr result " << std::endl;
493 return nullptr;
494 }
495
497 assert(GetData() );
498
499 RooAbsPdf * nullPdf = GetNullModel()->GetPdf();
501
502 // make conditional fit on null snapshot of poi
503
504 const RooArgSet * nullSnapshot = GetNullModel()->GetSnapshot();
505 assert(nullSnapshot && !nullSnapshot->empty());
506
507 // use as POI the nullSnapshot
508 // if more than one POI exists, consider only the first one
510
511 if (poiTest.size() > 1) {
512 oocoutW(nullptr,InputArguments) << "AsymptoticCalculator::GetHypoTest: snapshot has more than one POI - assume as POI first parameter " << std::endl;
513 }
514
515 std::unique_ptr<RooArgSet> allParams{nullPdf->getParameters(*GetData() )};
516 allParams->assign(fBestFitParams);
517
518 // set the one-side condition
519 // (this works when we have only one params of interest
520 RooRealVar * muHat = dynamic_cast<RooRealVar*> ( fBestFitPoi.first() );
521 assert(muHat && "no best fit parameter defined");
522 RooRealVar * muTest = dynamic_cast<RooRealVar*> ( nullSnapshot->find(muHat->GetName() ) );
523 assert(muTest && "poi snapshot is not existing");
524
525
526
527 if (verbose> 0) {
528 oocoutI(nullptr,Eval) << "\nAsymptoticCalculator::GetHypoTest: - perform an hypothesis test for POI ( " << muTest->GetName() << " ) = " << muTest->getVal() << std::endl;
529 oocoutP(nullptr,Eval) << "AsymptoticCalculator::GetHypoTest - Find best conditional NLL on OBSERVED data set ..... " << std::endl;
530 }
531
532 // evaluate the conditional NLL on the observed data for the snapshot value
533 double condNLL = EvaluateNLL(*GetNullModel(), const_cast<RooAbsData&>(*GetData()), &poiTest);
534
535 double qmu = 2.*(condNLL - fNLLObs);
536
537
538
539 if (verbose > 0)
540 oocoutP(nullptr,Eval) << "\t OBSERVED DATA : qmu = " << qmu << " condNLL = " << condNLL << " uncond " << fNLLObs << std::endl;
541
542
543 // this tolerance is used to avoid having negative qmu due to numerical errors
544 double tol = 2.E-3 * std::max(1.,ROOT::Math::MinimizerOptions::DefaultTolerance());
545 if (qmu < -tol || TMath::IsNaN(fNLLObs) ) {
546
547 if (qmu < 0) {
548 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a negative value of the qmu - retry to do the unconditional fit "
549 << std::endl;
550 } else {
551 oocoutW(nullptr, Minimization)
552 << "AsymptoticCalculator: unconditional fit failed before - retry to do it now " << std::endl;
553 }
554
555 double nll = EvaluateNLL(*GetNullModel(), const_cast<RooAbsData&>(*GetData()));
556
557 if (nll < fNLLObs || (TMath::IsNaN(fNLLObs) && !TMath::IsNaN(nll) ) ) {
558 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum "
559 << " old NLL = " << fNLLObs << " old muHat " << muHat->getVal() << std::endl;
560
561 // update values
562 fNLLObs = nll;
563 const RooArgSet * poi = GetNullModel()->GetParametersOfInterest();
564 assert(poi);
566 poi->snapshot(fBestFitPoi);
567 // restore also muHad since previous pointer has been deleted
568 muHat = dynamic_cast<RooRealVar*> ( fBestFitPoi.first() );
569 assert(muHat);
570
571 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: New minimum found for "
572 << " NLL = " << fNLLObs << " muHat " << muHat->getVal() << std::endl;
573
574
575 qmu = 2.*(condNLL - fNLLObs);
576
577 if (verbose > 0)
578 oocoutP(nullptr,Eval) << "After unconditional refit, new qmu value is " << qmu << std::endl;
579
580 }
581 }
582
583 if (qmu < -tol ) {
584 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: qmu is still < 0 for mu = "
585 << muTest->getVal() << " return a dummy result "
586 << std::endl;
587 return new HypoTestResult();
588 }
589 if (TMath::IsNaN(qmu) ) {
590 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
591 << muTest->getVal() << " return a dummy result "
592 << std::endl;
593 return new HypoTestResult();
594 }
595
596
597
598
599
600 // compute conditional ML on Asimov data set
601 // (need to const cast because it uses fitTo which is a non const method
602 // RooArgSet asimovGlobObs;
603 // RooAbsData * asimovData = (const_cast<AsymptoticCalculator*>(this))->MakeAsimovData( poi, asimovGlobObs);
604 // set global observables to their Asimov values
607 if (GetNullModel()->GetGlobalObservables() ) {
608 globObs.add(*GetNullModel()->GetGlobalObservables());
609 // store previous snapshot value
610 globObs.snapshot(globObsSnapshot);
611 globObs.assign(fAsimovGlobObs);
612 }
613
614
615 if (verbose > 0) oocoutP(nullptr,Eval) << "AsymptoticCalculator::GetHypoTest -- Find best conditional NLL on ASIMOV data set .... " << std::endl;
616
618
619
620 double qmu_A = 2.*(condNLL_A - fNLLAsimov );
621
622 if (verbose > 0)
623 oocoutP(nullptr,Eval) << "\t ASIMOV data qmu_A = " << qmu_A << " condNLL = " << condNLL_A << " uncond " << fNLLAsimov << std::endl;
624
625 if (qmu_A < -tol || TMath::IsNaN(fNLLAsimov) ) {
626
627 if (qmu_A < 0) {
628 oocoutW(nullptr, Minimization)
629 << "AsymptoticCalculator: Found a negative value of the qmu Asimov- retry to do the unconditional fit "
630 << std::endl;
631 } else {
632 oocoutW(nullptr, Minimization)
633 << "AsymptoticCalculator: Fit failed for unconditional the qmu Asimov- retry unconditional fit "
634 << std::endl;
635 }
636
637 double nll = EvaluateNLL(*GetNullModel(), *fAsimovData);
638
639 if (nll < fNLLAsimov || (TMath::IsNaN(fNLLAsimov) && !TMath::IsNaN(nll) )) {
640 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum for Asimov data set"
641 << " old NLL = " << fNLLAsimov << std::endl;
642
643 // update values
644 fNLLAsimov = nll;
645
646 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: New minimum found for "
647 << " NLL = " << fNLLAsimov << std::endl;
648 qmu_A = 2.*(condNLL_A - fNLLAsimov);
649
650 if (verbose > 0)
651 oocoutP(nullptr,Eval) << "After unconditional Asimov refit, new qmu_A value is " << qmu_A << std::endl;
652
653 }
654 }
655
656 if (qmu_A < - tol) {
657 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: qmu_A is still < 0 for mu = "
658 << muTest->getVal() << " return a dummy result "
659 << std::endl;
660 return new HypoTestResult();
661 }
662 if (TMath::IsNaN(qmu) ) {
663 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
664 << muTest->getVal() << " return a dummy result "
665 << std::endl;
666 return new HypoTestResult();
667 }
668
669
670 // restore previous value of global observables
671 globObs.assign(globObsSnapshot);
672
673 // now we compute p-values using the asymptotic formulae
674 // described in the paper
675 // Cowan et al, Eur.Phys.J. C (2011) 71:1554
676
677 // first try to guess automatically if needed to use qtilde (or ttilde in case of two sided)
678 // if explicitly fUseQTilde this was not set
679 // qtilde is in this case used if poi is bounded at the value of the alt hypothesis
680 // for Qtilde (need to distinguish case when qmu > qmuA = mu^2/ sigma^2)
681 // (see Cowan et al, Eur.Phys.J. C(2011) 71:1554 paper equations 64 and 65
682 // (remember qmu_A = mu^2/sigma^2 )
683 bool useQTilde = false;
684 // default case (check if poi is limited or not to a zero value)
685 if (!fOneSidedDiscovery) { // qtilde is not a discovery test
686 if (fUseQTilde == -1 && !fOneSidedDiscovery) {
687 // alternate snapshot is value for which background is zero (for limits)
688 RooRealVar * muAlt = dynamic_cast<RooRealVar*>(GetAlternateModel()->GetSnapshot()->first() );
689 // null snapshot is value for which background is zero (for discovery)
690 //RooRealVar * muNull = dynamic_cast<RooRealVar*>(GetNullModel()->GetSnapshot()->first() );
691 assert(muAlt != nullptr );
692 if (muTest->getMin() == muAlt->getVal() ) {
693 fUseQTilde = 1;
694 oocoutI(nullptr,InputArguments) << "Minimum of POI is " << muTest->getMin() << " corresponds to alt snapshot - using qtilde asymptotic formulae " << std::endl;
695 } else {
696 fUseQTilde = 0;
697 oocoutI(nullptr,InputArguments) << "Minimum of POI is " << muTest->getMin() << " is different to alt snapshot " << muAlt->getVal()
698 << " - using standard q asymptotic formulae " << std::endl;
699 }
700 }
702 }
703
704
705 //check for one side condition (remember this is valid only for one poi)
706 if (fOneSided ) {
707 if ( muHat->getVal() > muTest->getVal() ) {
708 oocoutI(nullptr,Eval) << "Using one-sided qmu - setting qmu to zero muHat = " << muHat->getVal()
709 << " muTest = " << muTest->getVal() << std::endl;
710 qmu = 0;
711 }
712 }
713 if (fOneSidedDiscovery ) {
714 if ( muHat->getVal() < muTest->getVal() ) {
715 oocoutI(nullptr,Eval) << "Using one-sided discovery qmu - setting qmu to zero muHat = " << muHat->getVal()
716 << " muTest = " << muTest->getVal() << std::endl;
717 qmu = 0;
718 }
719 }
720
721 // fix for negative qmu values due to numerical errors
722 if (qmu < 0 && qmu > -tol) qmu = 0;
724
725 // asymptotic formula for pnull and from paper Eur.Phys.J C 2011 71:1554
726 // we have 4 different cases:
727 // t(mu), t_tilde(mu) for the 2-sided
728 // q(mu) and q_tilde(mu) for the one -sided test statistics
729
730 double pnull = -1;
731 double palt = -1;
732
733 // asymptotic formula for pnull (for only one POI)
734 // From fact that qmu is a chi2 with ndf=1
735
736 double sqrtqmu = (qmu > 0) ? std::sqrt(qmu) : 0;
737 double sqrtqmu_A = (qmu_A > 0) ? std::sqrt(qmu_A) : 0;
738
739
741 // for one-sided PL (q_mu : equations 56,57)
742 if (verbose>2) {
743 if (fOneSided) {
744 oocoutI(nullptr,Eval) << "Using one-sided limit asymptotic formula (qmu)" << std::endl;
745 } else {
746 oocoutI(nullptr, Eval) << "Using one-sided discovery asymptotic formula (q0)" << std::endl;
747 }
748 }
751 }
752 else {
753 // for 2-sided PL (t_mu : equations 35,36 in asymptotic paper)
754 if (verbose > 2) oocoutI(nullptr,Eval) << "Using two-sided asymptotic formula (tmu)" << std::endl;
758
759 }
760
761 if (useQTilde ) {
762 if (fOneSided) {
763 // for bounded one-sided (q_mu_tilde: equations 64,65)
764 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) { // to avoid case 0/0
765 if (verbose > 2) oocoutI(nullptr,Eval) << "Using qmu_tilde (qmu is greater than qmu_A)" << std::endl;
766 pnull = ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
767 palt = ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
768 }
769 }
770 else {
771 // for 2 sided bounded test statistic (N.B there is no one sided discovery qtilde)
772 // t_mu_tilde: equations 43,44 in asymptotic paper
773 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) {
774 if (verbose > 2) oocoutI(nullptr,Eval) << "Using tmu_tilde (qmu is greater than qmu_A)" << std::endl;
776 ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
778 ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
779 }
780 }
781 }
782
783
784
785 // create an HypoTest result but where the sampling distributions are set to zero
786 string resultname = "HypoTestAsymptotic_result";
787 HypoTestResult* res = new HypoTestResult(resultname.c_str(), pnull, palt);
788
789 if (verbose > 0) {
790 oocoutP(nullptr, Eval) << "poi = " << muTest->getVal() << " qmu = " << qmu << " qmu_A = " << qmu_A
791 << " sigma = " << muTest->getVal() / sqrtqmu_A << " CLsplusb = " << pnull
792 << " CLb = " << palt << " CLs = " << res->CLs() << std::endl;
793 }
794
795 return res;
796
797}
798
800 PaltFunction( double offset, double pval, int icase) :
801 fOffset(offset), fPval(pval), fCase(icase) {}
802 double operator() (double x) const {
803 return ROOT::Math::normal_cdf_c(x + fOffset) + ROOT::Math::normal_cdf_c(fCase*(x - fOffset)) - fPval;
804 }
805 double fOffset;
806 double fPval;
807 int fCase;
808};
809
810////////////////////////////////////////////////////////////////////////////////
811/// function given the null and the alt p value - return the expected one given the N - sigma value
812
813double AsymptoticCalculator::GetExpectedPValues(double pnull, double palt, double nsigma, bool useCls, bool oneSided ) {
814 if (oneSided) {
818 if (!useCls) return clsplusb;
819 double clb = ROOT::Math::normal_cdf( nsigma, 1.);
820 return (clb == 0) ? -1 : clsplusb / clb;
821 }
822
823 // case of 2 sided test statistic
824 // need to compute numerically
826 if (sqrttmu == 0) {
827 // here cannot invert the function - skip the point
828 return -1;
829 }
830 // invert formula for palt to get sqrttmu_A
831 PaltFunction f( sqrttmu, palt, -1);
834 brf.SetFunction( wf, 0, 20);
835 bool ret = brf.Solve();
836 if (!ret) {
837 oocoutE(nullptr,Eval) << "Error finding expected p-values - return -1" << std::endl;
838 return -1;
839 }
840 double sqrttmu_A = brf.Root();
841
842 // now invert for expected value
845 brf.SetFunction(wf2,0,20);
846 ret = brf.Solve();
847 if (!ret) {
848 oocoutE(nullptr,Eval) << "Error finding expected p-values - return -1" << std::endl;
849 return -1;
850 }
851 return 2*ROOT::Math::normal_cdf_c( brf.Root(),1.);
852}
853
854namespace {
855
856////////////////////////////////////////////////////////////////////////////////
857/// Fill bins by looping recursively on observables.
858
859void FillBins(const RooAbsPdf & pdf, const RooArgList &obs, RooAbsData & data, int &index, double &binVolume, int &ibin) {
860
861 bool debug = (fgPrintLevel() >= 2);
862
863 RooRealVar * v = dynamic_cast<RooRealVar*>(&(obs[index]) );
864 if (!v) return;
865
866 RooArgSet obstmp(obs);
867 double expectedEvents = pdf.expectedEvents(obstmp);
868
869 if (debug) oocoutI(nullptr,Generation) << "looping on observable " << v->GetName() << std::endl;
870 for (int i = 0; i < v->getBins(); ++i) {
871 v->setBin(i);
872 if (index < int(obs.size()) -1) {
873 index++; // increase index
874 double prevBinVolume = binVolume;
875 binVolume *= v->getBinWidth(i); // increase bin volume
876 FillBins(pdf, obs, data, index, binVolume, ibin);
877 index--; // decrease index
878 binVolume = prevBinVolume; // decrease also bin volume
879 }
880 else {
881
882 // this is now a new bin - compute the pdf in this bin
883 double totBinVolume = binVolume * v->getBinWidth(i);
884 double fval = pdf.getVal(&obstmp)*totBinVolume;
885
886 if (fval*expectedEvents <= 0)
887 {
888 if (fval*expectedEvents < 0) {
889 oocoutW(nullptr,InputArguments)
890 << "AsymptoticCalculator::" << __func__
891 << "(): Bin " << i << " of " << v->GetName() << " has negative expected events! Please check your inputs." << std::endl;
892 }
893 else {
894 oocoutW(nullptr,InputArguments)
895 << "AsymptoticCalculator::" << __func__
896 << "(): Bin " << i << " of " << v->GetName() << " has zero expected events - skip it" << std::endl;
897 }
898 }
899 // have a cut off for overflows ??
900 else {
901 data.add(obs, fval*expectedEvents);
902 }
903
904 if (debug) {
905 oocoutI(nullptr,Generation) << "bin " << ibin << "\t";
906 for (std::size_t j=0; j < obs.size(); ++j) { ooccoutI(nullptr,Generation) << " " << (static_cast<RooRealVar&>( obs[j])).getVal(); }
907 ooccoutI(nullptr,Generation) << " w = " << fval*expectedEvents;
908 ooccoutI(nullptr,Generation) << std::endl;
909 }
910 ibin++;
911 }
912 }
913 //reset bin values
914 if (debug) {
915 oocoutI(nullptr,Generation) << "ending loop on .. " << v->GetName() << std::endl;
916 }
917
918 v->setBin(0);
919
920}
921
922bool setObsToExpected(std::span<RooAbsArg *> servers, const RooArgSet &obs, std::string const &errPrefix)
923{
924 RooRealVar *myobs = nullptr;
925 RooAbsReal *myexp = nullptr;
926 for (RooAbsArg *a : servers) {
927 if (obs.contains(*a)) {
928 if (myobs != nullptr) {
929 oocoutF(nullptr,Generation) << errPrefix << "Has two observables ?? " << std::endl;
930 return false;
931 }
932 myobs = dynamic_cast<RooRealVar *>(a);
933 if (myobs == nullptr) {
934 oocoutF(nullptr,Generation) << errPrefix << "Observable is not a RooRealVar??" << std::endl;
935 return false;
936 }
937 } else {
938 if (!a->isConstant() ) {
939 if (myexp != nullptr) {
940 oocoutE(nullptr,Generation) << errPrefix << "Has two non-const arguments " << std::endl;
941 return false;
942 }
943 myexp = dynamic_cast<RooAbsReal *>(a);
944 if (myexp == nullptr) {
945 oocoutF(nullptr,Generation) << errPrefix << "Expected is not a RooAbsReal??" << std::endl;
946 return false;
947 }
948 }
949 }
950 }
951 if (myobs == nullptr) {
952 oocoutF(nullptr,Generation) << errPrefix << "No observable?" << std::endl;
953 return false;
954 }
955 if (myexp == nullptr) {
956 oocoutF(nullptr,Generation) << errPrefix << "No observable?" << std::endl;
957 return false;
958 }
959
960 myobs->setVal(myexp->getVal());
961
962 if (fgPrintLevel() > 2) {
963 oocoutI(nullptr,Generation) << "SetObsToExpected : setting " << myobs->GetName() << " to expected value " << myexp->getVal() << " of " << myexp->GetName() << std::endl;
964 }
965
966 return true;
967}
968
969////////////////////////////////////////////////////////////////////////////////
970/// set observed value to the expected one
971/// works for Gaussian, Poisson or LogNormal
972/// assumes mean parameter value is the argument not constant and not depending on observables
973/// (if more than two arguments are not constant will use first one but print a warning !)
974/// need to iterate on the components of the Poisson to get n and nu (nu can be a RooAbsReal)
975/// (code from G. Petrucciani and extended by L.M.)
976
977bool SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs)
978{
979 std::string const &errPrefix = "AsymptoticCalculator::SetObsExpected( " + std::string{pdf.ClassName()} + " ) : ";
980 std::vector<RooAbsArg *> servers;
981 for (RooAbsArg *a : pdf.servers()) {
982 servers.emplace_back(a);
983 }
984 return setObsToExpected(servers, obs, errPrefix);
985}
986
988{
989 // In the case of the multi-variate Gaussian, we need to iterate over the
990 // dimensions and treat the servers for each dimension separately.
991
992 std::string const &errPrefix = "AsymptoticCalculator::SetObsExpected( " + std::string{mvgauss.ClassName()} + " ) : ";
993 std::vector<RooAbsArg *> servers{nullptr, nullptr};
994 bool ret = true;
995 for (std::size_t iDim = 0; iDim < mvgauss.xVec().size(); ++iDim) {
996 servers[0] = &mvgauss.xVec()[iDim];
997 servers[1] = &mvgauss.muVec()[iDim];
998 ret &= setObsToExpected(servers, obs, errPrefix + " : dim " + std::to_string(iDim) + " ");
999 }
1000 return ret;
1001}
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// Inpspect a product pdf to find all the Poisson or Gaussian parts to set the observed
1005/// values to expected ones.
1006
1007bool setObsToExpectedProdPdf(RooProdPdf &prod, const RooArgSet &obs)
1008{
1009 bool ret = true;
1010 for (auto *a : prod.pdfList()) {
1011 if (!a->dependsOn(obs)) continue;
1012 RooPoisson *pois = nullptr;
1013 RooGaussian *gauss = nullptr;
1014 RooMultiVarGaussian *mvgauss = nullptr;
1015 // should try to add also lognormal case ?
1016 if ((pois = dynamic_cast<RooPoisson *>(a)) != nullptr) {
1017 ret &= SetObsToExpected(*pois, obs);
1018 pois->setNoRounding(true); //needed since expected value is not an integer
1019 } else if ((gauss = dynamic_cast<RooGaussian *>(a)) != nullptr) {
1020 ret &= SetObsToExpected(*gauss, obs);
1021 } else if ((mvgauss = dynamic_cast<RooMultiVarGaussian *>(a)) != nullptr) {
1023 } else if (RooProdPdf *subprod = dynamic_cast<RooProdPdf *>(a)) {
1025 } else {
1026 oocoutE(nullptr, InputArguments)
1027 << "Illegal term in counting model: "
1028 << "the PDF " << a->GetName() << " depends on the observables, but is not a Poisson, Gaussian or Product"
1029 << std::endl;
1030 return false;
1031 }
1032 }
1033
1034 return ret;
1035}
1036
1037////////////////////////////////////////////////////////////////////////////////
1038/// Generate counting Asimov data for the case when the pdf cannot be extended.
1039/// This function assumes that the pdf is a RooPoisson or can be decomposed in a product of RooPoisson,
1040/// or is a RooGaussian. Otherwise, we cannot know how to make the Asimov data sets.
1041
1043 RooArgSet obs(observables);
1044 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
1045 RooPoisson *pois = nullptr;
1046 RooGaussian *gauss = nullptr;
1047 RooMultiVarGaussian *mvgauss = nullptr;
1048
1049 if (fgPrintLevel() > 1)
1050 oocoutI(nullptr,Generation) << "generate counting Asimov data for pdf of type " << pdf.ClassName() << std::endl;
1051
1052 bool r = false;
1053 if (prod != nullptr) {
1054 r = setObsToExpectedProdPdf(*prod, observables);
1055 } else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != nullptr) {
1056 r = SetObsToExpected(*pois, observables);
1057 // we need in this case to set Poisson to real values
1058 pois->setNoRounding(true);
1059 } else if ((gauss = dynamic_cast<RooGaussian *>(&pdf)) != nullptr) {
1060 r = SetObsToExpected(*gauss, observables);
1061 } else if ((mvgauss = dynamic_cast<RooMultiVarGaussian *>(&pdf)) != nullptr) {
1062 r = setObsToExpectedMultiVarGauss(*mvgauss, observables);
1063 } else {
1064 oocoutE(nullptr,InputArguments) << "A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << std::endl;
1065 }
1066 if (!r) return nullptr;
1067 int icat = 0;
1068 if (channelCat) {
1069 icat = channelCat->getCurrentIndex();
1070 }
1071
1072 RooDataSet *ret = new RooDataSet("CountingAsimovData" + std::to_string(icat),
1073 "CountingAsimovData" + std::to_string(icat), obs);
1074 ret->add(obs);
1075 return ret;
1076}
1077
1078////////////////////////////////////////////////////////////////////////////////
1079/// Compute the asimov data set for an observable of a pdf.
1080/// It generates binned data following the binning of the observables.
1081// TODO: (possibility to change number of bins)
1082// TODO: implement integration over bin content
1083
1085
1086 int printLevel = fgPrintLevel();
1087
1088 // Get observables defined by the pdf associated with this state
1089 std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );
1090
1091
1092 // if pdf cannot be extended assume is then a counting experiment
1093 if (!pdf.canBeExtended() ) return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);
1094
1095 RooArgSet obsAndWeight(*obs);
1096 obsAndWeight.add(weightVar);
1097
1098 std::unique_ptr<RooDataSet> asimovData;
1099 if (channelCat) {
1100 int icat = channelCat->getCurrentIndex();
1101 asimovData = std::make_unique<RooDataSet>("AsimovData" + std::to_string(icat),
1102 "combAsimovData" + std::to_string(icat),
1104 }
1105 else {
1106 asimovData = std::make_unique<RooDataSet>("AsimovData","AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));
1107 }
1108
1109 // This works only for 1D observables
1110 //RooRealVar* thisObs = ((RooRealVar*)obstmp->first());
1111
1112 RooArgList obsList(*obs);
1113
1114 // loop on observables and on the bins
1115 if (printLevel >= 2) {
1116 oocoutI(nullptr,Generation) << "Generating Asimov data for pdf " << pdf.GetName() << std::endl;
1117 oocoutI(nullptr,Generation) << "list of observables " << std::endl;
1118 obsList.Print();
1119 }
1120
1121 int obsIndex = 0;
1122 double binVolume = 1;
1123 int nbins = 0;
1124 FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
1125 if (printLevel >= 2)
1126 oocoutI(nullptr,Generation) << "filled from " << pdf.GetName() << " " << nbins << " nbins " << " volume is " << binVolume << std::endl;
1127
1128 // for (int iobs = 0; iobs < obsList.size(); ++iobs) {
1129 // RooRealVar * thisObs = dynamic_cast<RooRealVar*> &obsList[i];
1130 // if (thisObs == 0) continue;
1131 // // loop on the bin contents
1132 // for(int ibin=0; ibin<thisObs->numBins(); ++ibin){
1133 // thisObs->setBin(ibin);
1134
1135 // thisNorm=pdftmp->getVal(obstmp)*thisObs->getBinWidth(jj);
1136 // if (thisNorm*expectedEvents <= 0)
1137 // {
1138 // std::cout << "WARNING::Detected bin with zero expected events! Please check your inputs." << std::endl;
1139 // }
1140 // // have a cut off for overflows ??
1141 // obsDataUnbinned->add(*mc->GetObservables(), thisNorm*expectedEvents);
1142 // }
1143
1144 if (printLevel >= 1)
1145 {
1146 asimovData->Print();
1147 }
1148 if( TMath::IsNaN(asimovData->sumEntries()) ){
1149 oocoutE(nullptr,Generation) << "sum entries is nan"<< std::endl;
1150 assert(0);
1151 asimovData = nullptr;
1152 }
1153
1154 return asimovData.release();
1155
1156}
1157
1158} // namespace
1159
1160////////////////////////////////////////////////////////////////////////////////
1161/// generate the asimov data for the observables (not the global ones)
1162/// need to deal with the case of a sim pdf
1163/// \param pdf a RooAbsPdf. Can be also a RooSimultaneous, or a RooProPdf of a RooSimultaneous
1164
1166
1167 int printLevel = fgPrintLevel();
1168
1169 RooRealVar weightVar{"binWeightAsimov", "binWeightAsimov", 1, 0, 1.e30};
1170
1171 if (printLevel > 1) oocoutI(nullptr,Generation) <<" Generate Asimov data for observables"<< std::endl;
1172 // RooDataSet *simData = nullptr;
1173 const RooProdPdf *prodPdf = dynamic_cast<const RooProdPdf *>(&pdf);
1174 RooArgList strippedPdfSet("strippedPdfSet");
1175 if (prodPdf) {
1176 RooArgList list(prodPdf->pdfList());
1177 for (auto ele : list) {
1178 const RooAbsPdf *pdfi = dynamic_cast<const RooAbsPdf *>(ele);
1179 if (pdfi->dependsOn(observables))
1180 strippedPdfSet.add(*pdfi);
1181 }
1182 }
1183 RooProdPdf observableProdPdf("observableProdPdf", "observableProdPdf", strippedPdfSet);
1184 const RooAbsPdf *observablePdf;
1185 if (prodPdf) {
1186 if (strippedPdfSet.getSize() == 1)
1187 observablePdf = dynamic_cast<const RooAbsPdf *>(strippedPdfSet.at(0));
1188 else
1190 } else {
1191 observablePdf = &pdf;
1192 }
1193 const RooSimultaneous *simPdf = dynamic_cast<const RooSimultaneous *>(observablePdf);
1194 if (!simPdf) {
1195 // generate data for non sim pdf
1196 return GenerateAsimovDataSinglePdf(*observablePdf, observables, weightVar, nullptr);
1197 }
1198
1199 std::map<std::string, std::unique_ptr<RooDataSet>> asimovDataMap;
1200
1201 //look at category of simpdf
1202 RooCategory& channelCat = const_cast<RooCategory&>(dynamic_cast<const RooCategory&>(simPdf->indexCat()));
1203 int nrIndices = channelCat.numTypes();
1204 if( nrIndices == 0 ) {
1205 oocoutW(nullptr,Generation) << "Simultaneous pdf does not contain any categories." << std::endl;
1206 }
1207 for (int i=0;i<nrIndices;i++){
1208 channelCat.setIndex(i);
1209 //iFrame++;
1210 // Get pdf associated with state from simpdf
1211 RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getCurrentLabel()) ;
1212 assert(pdftmp != nullptr);
1213
1214 if (printLevel > 1)
1215 {
1216 oocoutI(nullptr,Generation) << "on type " << channelCat.getCurrentLabel() << " " << channelCat.getCurrentIndex() << std::endl;
1217 }
1218
1219 std::unique_ptr<RooDataSet> dataSinglePdf{static_cast<RooDataSet*>(GenerateAsimovDataSinglePdf( *pdftmp, observables, weightVar, &channelCat))};
1220 if (!dataSinglePdf) {
1221 oocoutE(nullptr,Generation) << "Error generating an Asimov data set for pdf " << pdftmp->GetName() << std::endl;
1222 return nullptr;
1223 }
1224
1225 if (asimovDataMap.count(string(channelCat.getCurrentLabel())) != 0) {
1226 oocoutE(nullptr,Generation) << "AsymptoticCalculator::GenerateAsimovData(): The PDF for " << channelCat.getCurrentLabel()
1227 << " was already defined. It will be overridden. The faulty category definitions follow:" << std::endl;
1228 channelCat.Print("V");
1229 }
1230
1231 if (printLevel > 1)
1232 {
1233 oocoutI(nullptr,Generation) << "channel: " << channelCat.getCurrentLabel() << ", data: ";
1234 dataSinglePdf->Print();
1235 ooccoutI(nullptr,Generation) << std::endl;
1236 }
1237
1238 asimovDataMap[string(channelCat.getCurrentLabel())] = std::move(dataSinglePdf);
1239 }
1240
1241 RooArgSet obsAndWeight(observables);
1242 obsAndWeight.add(weightVar);
1243
1244
1245 return new RooDataSet("asimovDataFullModel","asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
1247}
1248
1249////////////////////////////////////////////////////////////////////////////////
1250/// Make the Asimov data from the ModelConfig and list of poi
1251/// \param realData Real data
1252/// \param model Model config defining the pdf and the parameters
1253/// \param paramValues The snapshot of POI and parameters used for finding the best nuisance parameter values (conditioned at these values)
1254/// \param[out] asimovGlobObs Global observables set to values satisfying the constraints
1255/// \param genPoiValues Optional. A different set of POI values used for generating. By default the same POI are used for generating and for finding the nuisance parameters
1256/// given an observed data set, a model and a snapshot of the poi.
1257/// \return The asimov data set. The user takes ownership.
1258///
1259
1261
1262 int verbose = fgPrintLevel();
1263
1264
1265 RooArgSet poi(*model.GetParametersOfInterest());
1266 poi.assign(paramValues);
1267
1268 // set poi constant for conditional MLE
1269 // need to fit nuisance parameters at their conditional MLE value
1271 for (auto *tmpPar : static_range_cast<RooRealVar *>(poi)) {
1272 tmpPar->setConstant();
1273 if (verbose>0)
1274 oocoutI(nullptr,Generation) << "MakeAsimov: Setting poi " << tmpPar->GetName() << " to a constant value = " << tmpPar->getVal() << std::endl;
1276 }
1277
1278 // find conditional value of the nuisance parameters
1279 bool hasFloatParams = false;
1281 if (model.GetNuisanceParameters()) {
1284 if (!constrainParams.empty()) hasFloatParams = true;
1285
1286 } else {
1287 // Do we have free parameters anyway that need fitting?
1288 std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
1289 for (auto const *rrv : dynamic_range_cast<RooRealVar *>(*params)) {
1290 if ( rrv != nullptr && rrv->isConstant() == false ) { hasFloatParams = true; break; }
1291 }
1292 }
1293 if (hasFloatParams) {
1294 // models need to be fitted to find best nuisance parameter values
1295
1296 TStopwatch tw2; tw2.Start();
1298 if (verbose>0) {
1299 oocoutP(nullptr,Generation) << "MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
1300 minimPrintLevel = verbose;
1301 if (verbose>1) {
1302 oocoutI(nullptr,Generation) << "POI values:\n"; poi.Print("v");
1303 if (verbose > 2) {
1304 oocoutI(nullptr,Generation) << "Nuis param values:\n";
1305 constrainParams.Print("v");
1306 }
1307 }
1308 }
1310 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
1311
1313 std::vector<RooCmdArg> args{
1314 RooFit::Minimizer("",minimizerAlgo.c_str()), // empty mimimizer type to select default
1317 RooFit::Hesse(false),
1319 RooFit::Offset(GetGlobalRooStatsConfig().useLikelihoodOffset),
1321 };
1322
1323 RooLinkedList argList;
1324 for (auto& arg : args) {
1325 argList.Add(&arg);
1326 }
1327 model.fitTo(realData, argList);
1328 if (verbose>0) {
1329 tw2.Stop();
1330 oocoutP(nullptr,Generation) << "fit time : " << tw2.RealTime() << " s (real) " << tw2.CpuTime() << " s (cpu)" << std::endl;
1331 }
1332 if (verbose > 1) {
1333 // after the fit the nuisance parameters will have their best fit value
1334 if (model.GetNuisanceParameters() ) {
1335 oocoutI(nullptr,Generation) << "Nuisance parameters after fit for asimov dataset: " << std::endl;
1336 model.GetNuisanceParameters()->Print("V");
1337 }
1338 }
1339
1340 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
1341
1342 }
1343
1344 // restore the parameters which were set constant
1346
1347 std::unique_ptr<RooArgSet> allParams{model.GetPdf()->getParameters(realData)};
1348
1350
1351 // if a RooArgSet of poi is passed , different poi will be used for generating the Asimov data set
1352 if (genPoiValues) {
1353 allParams->assign(*genPoiValues);
1354 }
1355
1356 // now do the actual generation of the AsimovData Set
1357 // no need to pass parameters values since we have set them before
1358 return MakeAsimovData(model, *allParams, asimovGlobObs);
1359}
1360
1361////////////////////////////////////////////////////////////////////////////////
1362/// \param model ModelConfig that contains the model pdf and the model parameters
1363/// \param allParamValues The parameters of the model will be set to the values given in this set
1364/// \param[out] asimovGlobObs Global observables set to values satisfying the constraints
1365/// \return Asimov data set. The user takes ownership.
1366///
1367/// The parameter values (including the nuisance parameter) can result from a fit to data or be at the nominal values.
1368///
1369
1371
1372 int verbose = fgPrintLevel();
1373
1374 TStopwatch tw;
1375 tw.Start();
1376
1377 // set the parameter values (do I need the poi to be constant ? )
1378 // the nuisance parameter values could be set at their fitted value (the MLE)
1379 if (!allParamValues.empty()) {
1380 std::unique_ptr<RooArgSet> allVars{model.GetPdf()->getVariables()};
1381 allVars->assign(allParamValues);
1382 }
1383
1384
1385 // generate the Asimov data set for the observables
1386 RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
1387
1388 if (verbose>0) {
1389 oocoutI(nullptr,Generation) << "Generated Asimov data for observables "; (model.GetObservables() )->Print();
1390 if (verbose > 1) {
1391 if (asimov->numEntries() == 1 ) {
1392 oocoutI(nullptr,Generation) << "--- Asimov data values \n";
1393 asimov->get()->Print("v");
1394 }
1395 else {
1396 oocoutI(nullptr,Generation) << "--- Asimov data numEntries = " << asimov->numEntries() << " sumOfEntries = " << asimov->sumEntries() << std::endl;
1397 }
1398 tw.Stop();
1399 oocoutI(nullptr,Generation) << "\ttime for generating : " << tw.RealTime() << " s (real) " << tw.CpuTime() << " s (cpu)" << std::endl;
1400 }
1401 }
1402
1403
1404 // Now need to have in ASIMOV the data sets also the global observables
1405 // Their values must be the one satisfying the constraint.
1406 // to do it make a nuisance pdf with all product of constraints and then
1407 // assign to each constraint a glob observable value = to the current fitted nuisance parameter value
1408 // IN general one should solve in general the system of equations f( gobs| nuispar ) = 0 where f are the
1409 // derivatives of the constraint with respect the nuisance parameter and they are evaluated at the best fit nuisance
1410 // parameter points
1411 // As simple solution assume that constrain has a direct dependence on the nuisance parameter, i.e.
1412 // Constraint (gobs, func( nuispar) ) and the condition is satisfied for
1413 // gobs = func( nuispar) where nunispar is at the MLE value
1414
1415
1416 if (model.GetGlobalObservables() && !model.GetGlobalObservables()->empty()) {
1417
1418 if (verbose>1) {
1419 oocoutI(nullptr,Generation) << "Generating Asimov data for global observables " << std::endl;
1420 }
1421
1423
1424 // snapshot data global observables
1426 SetAllConstant(gobs, true);
1427 gobs.snapshot(snapGlobalObsData);
1428
1429
1431 if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
1432 if (nuis.empty()) {
1433 oocoutW(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
1434 << " set global observables to model values " << std::endl;
1435 asimovGlobObs.assign(gobs);
1436 return asimov;
1437 }
1438
1439 // part 1: create the nuisance pdf
1440 std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,"TempNuisPdf") );
1441 if (nuispdf == nullptr) {
1442 oocoutF(nullptr, Generation) << "AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and "
1443 "global obs but no nuisance pdf "
1444 << std::endl;
1445 }
1446 // unfold the nuisance pdf if it is a prod pdf
1447 RooArgList pdfList;
1448 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(nuispdf.get());
1449 if (prod ) {
1450 pdfList.add(prod->pdfList());
1451 } else {
1452 // nothing to unfold - just use the pdf
1453 pdfList.add(*nuispdf);
1454 }
1455
1456 for (auto *cterm : static_range_cast<RooAbsPdf *>(pdfList)) {
1457 assert(dynamic_cast<RooAbsPdf *>(static_cast<RooAbsArg *>(cterm)) &&
1458 "AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
1459
1460 if (!cterm->dependsOn(nuis)) continue; // dummy constraints
1461 // skip also the case of uniform components
1462 if (typeid(*cterm) == typeid(RooUniform)) continue;
1463
1464 std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
1465 std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
1466 if (cgobs->size() > 1) {
1467 oocoutE(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1468 << " has multiple global observables -cannot generate - skip it" << std::endl;
1469 continue;
1470 }
1471 else if (cgobs->empty()) {
1472 oocoutW(nullptr, Generation)
1473 << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1474 << " has no global observables - skip it" << std::endl;
1475 continue;
1476 }
1477 // the variable representing the global observable
1478 RooRealVar &rrv = dynamic_cast<RooRealVar &>(*cgobs->first());
1479
1480 // remove the constant parameters in cpars
1482 if (cpars->size() != 1) {
1483 oocoutE(nullptr, Generation)
1484 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1485 << cterm->GetName() << " has multiple floating params - cannot generate - skip it " << std::endl;
1486 continue;
1487 }
1488
1489 bool foundServer = false;
1490 // note : this will work only for this type of constraints
1491 // expressed as RooPoisson, RooGaussian, RooLognormal, RooGamma
1492 TClass * cClass = cterm->IsA();
1493 if (verbose > 2) oocoutI(nullptr,Generation) << "Constraint " << cterm->GetName() << " of type " << cClass->GetName() << std::endl;
1497 TString className = (cClass) ? cClass->GetName() : "undefined";
1498 oocoutW(nullptr, Generation)
1499 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1500 << cterm->GetName() << " of type " << className
1501 << " is a non-supported type - result might be not correct " << std::endl;
1502 }
1503
1504 // in case of a Poisson constraint make sure the rounding is not set
1505 if (cClass == RooPoisson::Class() ) {
1506 RooPoisson * pois = static_cast<RooPoisson*>(cterm);
1507 assert(dynamic_cast<RooPoisson *>(cterm));
1508 pois->setNoRounding(true);
1509 }
1510
1511 // look at server of the constraint term and check if the global observable is part of the server
1512 RooAbsArg * arg = cterm->findServer(rrv);
1513 if (!arg) {
1514 // special case is for the Gamma where one might define the global observable n and you have a Gamma(b, n+1, ...._
1515 // in this case n+1 is the server and we don;t have a direct dependency, but we want to set n to the b value
1516 // so in case of the Gamma ignore this test
1517 if ( cClass != RooGamma::Class() ) {
1518 oocoutE(nullptr, Generation)
1519 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1520 << cterm->GetName() << " has no direct dependence on global observable- cannot generate it " << std::endl;
1521 continue;
1522 }
1523 }
1524
1525 // loop on the server of the constraint term
1526 // need to treat the Gamma as a special case
1527 // the mode of the Gamma is (k-1)*theta where theta is the inverse of the rate parameter.
1528 // we assume that the global observable is defined as ngobs = k-1 and the theta parameter has the name theta otherwise we use other procedure which might be wrong
1529 RooAbsReal * thetaGamma = nullptr;
1530 if ( cClass == RooGamma::Class() ) {
1531 for (RooAbsArg *a2 : cterm->servers()) {
1532 if (TString(a2->GetName()).Contains("theta") ) {
1533 thetaGamma = dynamic_cast<RooAbsReal*>(a2);
1534 break;
1535 }
1536 }
1537 if (thetaGamma == nullptr) {
1538 oocoutI(nullptr, Generation)
1539 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1540 << cterm->GetName() << " is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is 1 " << std::endl;
1541 }
1542 else if (verbose>2) {
1543 oocoutI(nullptr,Generation) << "Gamma constraint has a scale " << thetaGamma->GetName() << " = " << thetaGamma->getVal() << std::endl;
1544 }
1545 }
1546 for (RooAbsArg *a2 : cterm->servers()) {
1547 RooAbsReal * rrv2 = dynamic_cast<RooAbsReal *>(a2);
1548 if (verbose > 2) oocoutI(nullptr,Generation) << "Loop on constraint server term " << a2->GetName() << std::endl;
1549 if (rrv2 && rrv2->dependsOn(nuis) ) {
1550
1551
1552 // found server depending on nuisance
1553 if (foundServer) {
1554 oocoutE(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1555 << cterm->GetName() << " constraint term has more server depending on nuisance- cannot generate it " <<
1556 std::endl;
1557 foundServer = false;
1558 break;
1559 }
1560 if (thetaGamma && thetaGamma->getVal() > 0) {
1561 rrv.setVal( rrv2->getVal() / thetaGamma->getVal() );
1562 } else {
1563 rrv.setVal(rrv2->getVal());
1564 }
1565 foundServer = true;
1566
1567 if (verbose > 2) {
1568 oocoutI(nullptr,Generation) << "setting global observable " << rrv.GetName() << " to value " << rrv.getVal()
1569 << " which comes from " << rrv2->GetName() << std::endl;
1570 }
1571 }
1572 }
1573
1574 if (!foundServer) {
1575 oocoutE(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData - can't find nuisance for constraint term - global observables will not be set to Asimov value " << cterm->GetName() << std::endl;
1576 oocoutE(nullptr,Generation) << "Parameters: " << std::endl;
1577 cpars->Print("V");
1578 oocoutE(nullptr,Generation) << "Observables: " << std::endl;
1579 cgobs->Print("V");
1580 }
1581 }
1582
1583 // make a snapshot of global observables
1584 // needed this ?? (LM)
1585
1586 asimovGlobObs.removeAll();
1587 SetAllConstant(gobs, true);
1588 gobs.snapshot(asimovGlobObs);
1589
1590 // revert global observables to the data value
1591 gobs.assign(snapGlobalObsData);
1592
1593 if (verbose>0) {
1594 oocoutI(nullptr,Generation) << "Generated Asimov data for global observables ";
1595 if (verbose == 1) gobs.Print();
1596 }
1597
1598 if (verbose > 1) {
1599 oocoutI(nullptr,Generation) << "\nGlobal observables for data: " << std::endl;
1600 gobs.Print("V");
1601 oocoutI(nullptr,Generation) << "\nGlobal observables for asimov: " << std::endl;
1602 asimovGlobObs.Print("V");
1603 }
1604
1605
1606 }
1607
1608 return asimov;
1609
1610}
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define oocoutW(o, a)
#define oocoutE(o, a)
#define oocoutF(o, a)
#define oocoutI(o, a)
#define ooccoutI(o, a)
#define ooccoutP(o, a)
#define oocoutP(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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 char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
TRObject operator()(const T1 &t1) const
Class for finding the root of a one dimensional function using the Brent algorithm.
static const std::string & DefaultMinimizerAlgo()
Template class to wrap any C++ callable object which takes one argument i.e.
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
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.
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
bool contains(const char *name) const
Check if collection contains an argument with a specific name.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
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:56
virtual double sumEntries() const =0
Return effective number of entries in dataset, i.e., sum all weights.
virtual const RooArgSet * get() const
Definition RooAbsData.h:100
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:214
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
static bool hideOffset()
static void setHideOffset(bool flag)
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
static TClass * Class()
Object to represent discrete states.
Definition RooCategory.h:28
static TClass * Class()
Container class to hold unbinned data.
Definition RooDataSet.h:32
static TClass * Class()
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
static TClass * Class()
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
virtual void Add(TObject *arg)
static TClass * Class()
Wrapper class around ROOT::Math::Minimizer that provides a seamless interface between the minimizer f...
static RooMsgService & instance()
Return reference to singleton instance.
Multivariate Gaussian p.d.f.
Poisson pdf.
Definition RooPoisson.h:19
static TClass * Class()
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:36
const RooArgList & pdfList() const
Definition RooProdPdf.h:70
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
static double GetExpectedPValues(double pnull, double palt, double nsigma, bool usecls, bool oneSided=true)
function given the null and the alt p value - return the expected one given the N - sigma value
static void SetPrintLevel(int level)
set print level (static function)
RooArgSet fAsimovGlobObs
snapshot of Asimov global observables
static RooAbsData * GenerateAsimovData(const RooAbsPdf &pdf, const RooArgSet &observables)
generate the asimov data for the observables (not the global ones) need to deal with the case of a si...
int fUseQTilde
flag to indicate if using qtilde or not (-1 (default based on RooRealVar)), 0 false,...
bool fIsInitialized
! flag to check if calculator is initialized
HypoTestResult * GetHypoTest() const override
re-implement HypoTest computation using the asymptotic
bool fOneSided
for one sided PL test statistic (upper limits)
RooArgSet fBestFitParams
snapshot of all best fitted Parameter values
AsymptoticCalculator(RooAbsData &data, const ModelConfig &altModel, const ModelConfig &nullModel, bool nominalAsimov=false)
constructor for asymptotic calculator from Data set and ModelConfig
bool fOneSidedDiscovery
for one sided PL test statistic (for discovery)
RooAbsData * fAsimovData
asimov data set
RooArgSet fBestFitPoi
snapshot of best fitted POI values
static RooAbsData * MakeAsimovData(RooAbsData &data, const ModelConfig &model, const RooArgSet &poiValues, RooArgSet &globObs, const RooArgSet *genPoiValues=nullptr)
Make Asimov data.
bool fNominalAsimov
make Asimov at nominal parameter values
bool Initialize() const
initialize the calculator by performing a global fit and make the Asimov data set
Common base class for the Hypothesis Test Calculators.
const ModelConfig * GetNullModel(void) const
const ModelConfig * GetAlternateModel(void) const
HypoTestResult is a base class for results from hypothesis tests.
virtual double CLs() const
is simply (not a method, but a quantity)
< A class that holds configuration information for a model using a workspace as a store
Definition ModelConfig.h:34
std::unique_ptr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &...cmdArgs) const
Wrapper around RooAbsPdf::fitTo(), where the pdf and some configuration options are retrieved from th...
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
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)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
Flat p.d.f.
Definition RooUniform.h:24
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:224
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:138
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
RooCmdArg Index(RooCategory &icat)
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Import(const char *state, TH1 &histo)
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 EvalErrorWall(bool flag)
RooCmdArg PrintLevel(Int_t code)
double normal_cdf_c(double x, double sigma=1, double x0=0)
Complement of the cumulative distribution function of the normal (Gaussian) distribution (upper tail)...
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
double normal_quantile(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the lower tail of the normal (Gaussian) distri...
double normal_quantile_c(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the upper tail of the normal (Gaussian) distri...
Double_t x[n]
Definition legend1.C:17
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:404
MsgLevel
Verbosity level for RooMsgService::StreamConfig in RooMsgService.
Namespace for the RooStats classes.
Definition CodegenImpl.h:66
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
std::string const & NLLOffsetMode()
Test what offsetting mode RooStats should use by default.
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
RooStatsConfig & GetGlobalRooStatsConfig()
Retrieve the config object which can be used to set flags for things like offsetting the likelihood o...
Bool_t IsNaN(Double_t x)
Definition TMath.h:903
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:913
PaltFunction(double offset, double pval, int icase)