Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 std::cout << "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 std::cout << "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 std::cout << " ----> Doing a re-scan first\n";
402 minim.minimize(minimizer,"Scan");
403 }
404 if (tries == 2) {
406 std::cout << " ----> trying with strategy = 1\n";
407 minim.setStrategy(1);
408 }
409 else
410 tries++; // skip this trial if strategy is already 1
411 }
412 if (tries == 3) {
413 std::cout << " ----> trying with improve\n";
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 std::cout << "AsymptoticCalculator::EvaluateNLL - value = " << val;
448 if (poiSet) {
449 muTest = ( static_cast<RooRealVar*>(poiSet->first()) )->getVal();
450 std::cout << " for poi fixed at = " << muTest;
451 }
452 if (!skipFit) {
453 std::cout << "\tfit time : ";
454 tw.Print();
455 } else {
456 std::cout << 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 std::cout << std::endl;
529 oocoutI(nullptr,Eval) << "AsymptoticCalculator::GetHypoTest: - perform an hypothesis test for POI ( " << muTest->GetName() << " ) = " << muTest->getVal() << std::endl;
530 oocoutP(nullptr,Eval) << "AsymptoticCalculator::GetHypoTest - Find best conditional NLL on OBSERVED data set ..... " << std::endl;
531 }
532
533 // evaluate the conditional NLL on the observed data for the snapshot value
534 double condNLL = EvaluateNLL(*GetNullModel(), const_cast<RooAbsData&>(*GetData()), &poiTest);
535
536 double qmu = 2.*(condNLL - fNLLObs);
537
538
539
540 if (verbose > 0)
541 oocoutP(nullptr,Eval) << "\t OBSERVED DATA : qmu = " << qmu << " condNLL = " << condNLL << " uncond " << fNLLObs << std::endl;
542
543
544 // this tolerance is used to avoid having negative qmu due to numerical errors
545 double tol = 2.E-3 * std::max(1.,ROOT::Math::MinimizerOptions::DefaultTolerance());
546 if (qmu < -tol || TMath::IsNaN(fNLLObs) ) {
547
548 if (qmu < 0) {
549 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a negative value of the qmu - retry to do the unconditional fit "
550 << std::endl;
551 } else {
552 oocoutW(nullptr, Minimization)
553 << "AsymptoticCalculator: unconditional fit failed before - retry to do it now " << std::endl;
554 }
555
556 double nll = EvaluateNLL(*GetNullModel(), const_cast<RooAbsData&>(*GetData()));
557
558 if (nll < fNLLObs || (TMath::IsNaN(fNLLObs) && !TMath::IsNaN(nll) ) ) {
559 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum "
560 << " old NLL = " << fNLLObs << " old muHat " << muHat->getVal() << std::endl;
561
562 // update values
563 fNLLObs = nll;
564 const RooArgSet * poi = GetNullModel()->GetParametersOfInterest();
565 assert(poi);
567 poi->snapshot(fBestFitPoi);
568 // restore also muHad since previous pointer has been deleted
569 muHat = dynamic_cast<RooRealVar*> ( fBestFitPoi.first() );
570 assert(muHat);
571
572 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: New minimum found for "
573 << " NLL = " << fNLLObs << " muHat " << muHat->getVal() << std::endl;
574
575
576 qmu = 2.*(condNLL - fNLLObs);
577
578 if (verbose > 0)
579 oocoutP(nullptr,Eval) << "After unconditional refit, new qmu value is " << qmu << std::endl;
580
581 }
582 }
583
584 if (qmu < -tol ) {
585 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: qmu is still < 0 for mu = "
586 << muTest->getVal() << " return a dummy result "
587 << std::endl;
588 return new HypoTestResult();
589 }
590 if (TMath::IsNaN(qmu) ) {
591 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
592 << muTest->getVal() << " return a dummy result "
593 << std::endl;
594 return new HypoTestResult();
595 }
596
597
598
599
600
601 // compute conditional ML on Asimov data set
602 // (need to const cast because it uses fitTo which is a non const method
603 // RooArgSet asimovGlobObs;
604 // RooAbsData * asimovData = (const_cast<AsymptoticCalculator*>(this))->MakeAsimovData( poi, asimovGlobObs);
605 // set global observables to their Asimov values
608 if (GetNullModel()->GetGlobalObservables() ) {
609 globObs.add(*GetNullModel()->GetGlobalObservables());
610 // store previous snapshot value
611 globObs.snapshot(globObsSnapshot);
612 globObs.assign(fAsimovGlobObs);
613 }
614
615
616 if (verbose > 0) oocoutP(nullptr,Eval) << "AsymptoticCalculator::GetHypoTest -- Find best conditional NLL on ASIMOV data set .... " << std::endl;
617
619
620
621 double qmu_A = 2.*(condNLL_A - fNLLAsimov );
622
623 if (verbose > 0)
624 oocoutP(nullptr,Eval) << "\t ASIMOV data qmu_A = " << qmu_A << " condNLL = " << condNLL_A << " uncond " << fNLLAsimov << std::endl;
625
626 if (qmu_A < -tol || TMath::IsNaN(fNLLAsimov) ) {
627
628 if (qmu_A < 0) {
629 oocoutW(nullptr, Minimization)
630 << "AsymptoticCalculator: Found a negative value of the qmu Asimov- retry to do the unconditional fit "
631 << std::endl;
632 } else {
633 oocoutW(nullptr, Minimization)
634 << "AsymptoticCalculator: Fit failed for unconditional the qmu Asimov- retry unconditional fit "
635 << std::endl;
636 }
637
638 double nll = EvaluateNLL(*GetNullModel(), *fAsimovData);
639
640 if (nll < fNLLAsimov || (TMath::IsNaN(fNLLAsimov) && !TMath::IsNaN(nll) )) {
641 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: Found a better unconditional minimum for Asimov data set"
642 << " old NLL = " << fNLLAsimov << std::endl;
643
644 // update values
645 fNLLAsimov = nll;
646
647 oocoutW(nullptr,Minimization) << "AsymptoticCalculator: New minimum found for "
648 << " NLL = " << fNLLAsimov << std::endl;
649 qmu_A = 2.*(condNLL_A - fNLLAsimov);
650
651 if (verbose > 0)
652 oocoutP(nullptr,Eval) << "After unconditional Asimov refit, new qmu_A value is " << qmu_A << std::endl;
653
654 }
655 }
656
657 if (qmu_A < - tol) {
658 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: qmu_A is still < 0 for mu = "
659 << muTest->getVal() << " return a dummy result "
660 << std::endl;
661 return new HypoTestResult();
662 }
663 if (TMath::IsNaN(qmu) ) {
664 oocoutE(nullptr,Minimization) << "AsymptoticCalculator: failure in fitting for qmu or qmuA "
665 << muTest->getVal() << " return a dummy result "
666 << std::endl;
667 return new HypoTestResult();
668 }
669
670
671 // restore previous value of global observables
672 globObs.assign(globObsSnapshot);
673
674 // now we compute p-values using the asymptotic formulae
675 // described in the paper
676 // Cowan et al, Eur.Phys.J. C (2011) 71:1554
677
678 // first try to guess automatically if needed to use qtilde (or ttilde in case of two sided)
679 // if explicitly fUseQTilde this was not set
680 // qtilde is in this case used if poi is bounded at the value of the alt hypothesis
681 // for Qtilde (need to distinguish case when qmu > qmuA = mu^2/ sigma^2)
682 // (see Cowan et al, Eur.Phys.J. C(2011) 71:1554 paper equations 64 and 65
683 // (remember qmu_A = mu^2/sigma^2 )
684 bool useQTilde = false;
685 // default case (check if poi is limited or not to a zero value)
686 if (!fOneSidedDiscovery) { // qtilde is not a discovery test
687 if (fUseQTilde == -1 && !fOneSidedDiscovery) {
688 // alternate snapshot is value for which background is zero (for limits)
689 RooRealVar * muAlt = dynamic_cast<RooRealVar*>(GetAlternateModel()->GetSnapshot()->first() );
690 // null snapshot is value for which background is zero (for discovery)
691 //RooRealVar * muNull = dynamic_cast<RooRealVar*>(GetNullModel()->GetSnapshot()->first() );
692 assert(muAlt != nullptr );
693 if (muTest->getMin() == muAlt->getVal() ) {
694 fUseQTilde = 1;
695 oocoutI(nullptr,InputArguments) << "Minimum of POI is " << muTest->getMin() << " corresponds to alt snapshot - using qtilde asymptotic formulae " << std::endl;
696 } else {
697 fUseQTilde = 0;
698 oocoutI(nullptr,InputArguments) << "Minimum of POI is " << muTest->getMin() << " is different to alt snapshot " << muAlt->getVal()
699 << " - using standard q asymptotic formulae " << std::endl;
700 }
701 }
703 }
704
705
706 //check for one side condition (remember this is valid only for one poi)
707 if (fOneSided ) {
708 if ( muHat->getVal() > muTest->getVal() ) {
709 oocoutI(nullptr,Eval) << "Using one-sided qmu - setting qmu to zero muHat = " << muHat->getVal()
710 << " muTest = " << muTest->getVal() << std::endl;
711 qmu = 0;
712 }
713 }
714 if (fOneSidedDiscovery ) {
715 if ( muHat->getVal() < muTest->getVal() ) {
716 oocoutI(nullptr,Eval) << "Using one-sided discovery qmu - setting qmu to zero muHat = " << muHat->getVal()
717 << " muTest = " << muTest->getVal() << std::endl;
718 qmu = 0;
719 }
720 }
721
722 // fix for negative qmu values due to numerical errors
723 if (qmu < 0 && qmu > -tol) qmu = 0;
725
726 // asymptotic formula for pnull and from paper Eur.Phys.J C 2011 71:1554
727 // we have 4 different cases:
728 // t(mu), t_tilde(mu) for the 2-sided
729 // q(mu) and q_tilde(mu) for the one -sided test statistics
730
731 double pnull = -1;
732 double palt = -1;
733
734 // asymptotic formula for pnull (for only one POI)
735 // From fact that qmu is a chi2 with ndf=1
736
737 double sqrtqmu = (qmu > 0) ? std::sqrt(qmu) : 0;
738 double sqrtqmu_A = (qmu_A > 0) ? std::sqrt(qmu_A) : 0;
739
740
742 // for one-sided PL (q_mu : equations 56,57)
743 if (verbose>2) {
744 if (fOneSided) {
745 oocoutI(nullptr,Eval) << "Using one-sided limit asymptotic formula (qmu)" << std::endl;
746 } else {
747 oocoutI(nullptr, Eval) << "Using one-sided discovery asymptotic formula (q0)" << std::endl;
748 }
749 }
752 }
753 else {
754 // for 2-sided PL (t_mu : equations 35,36 in asymptotic paper)
755 if (verbose > 2) oocoutI(nullptr,Eval) << "Using two-sided asymptotic formula (tmu)" << std::endl;
759
760 }
761
762 if (useQTilde ) {
763 if (fOneSided) {
764 // for bounded one-sided (q_mu_tilde: equations 64,65)
765 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) { // to avoid case 0/0
766 if (verbose > 2) oocoutI(nullptr,Eval) << "Using qmu_tilde (qmu is greater than qmu_A)" << std::endl;
767 pnull = ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
768 palt = ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
769 }
770 }
771 else {
772 // for 2 sided bounded test statistic (N.B there is no one sided discovery qtilde)
773 // t_mu_tilde: equations 43,44 in asymptotic paper
774 if ( qmu > qmu_A && (qmu_A > 0 || qmu > tol) ) {
775 if (verbose > 2) oocoutI(nullptr,Eval) << "Using tmu_tilde (qmu is greater than qmu_A)" << std::endl;
777 ROOT::Math::normal_cdf_c( (qmu + qmu_A)/(2 * sqrtqmu_A), 1.);
779 ROOT::Math::normal_cdf_c( (qmu - qmu_A)/(2 * sqrtqmu_A), 1.);
780 }
781 }
782 }
783
784
785
786 // create an HypoTest result but where the sampling distributions are set to zero
787 string resultname = "HypoTestAsymptotic_result";
788 HypoTestResult* res = new HypoTestResult(resultname.c_str(), pnull, palt);
789
790 if (verbose > 0) {
791 oocoutP(nullptr, Eval) << "poi = " << muTest->getVal() << " qmu = " << qmu << " qmu_A = " << qmu_A
792 << " sigma = " << muTest->getVal() / sqrtqmu_A << " CLsplusb = " << pnull
793 << " CLb = " << palt << " CLs = " << res->CLs() << std::endl;
794 }
795
796 return res;
797
798}
799
801 PaltFunction( double offset, double pval, int icase) :
802 fOffset(offset), fPval(pval), fCase(icase) {}
803 double operator() (double x) const {
804 return ROOT::Math::normal_cdf_c(x + fOffset) + ROOT::Math::normal_cdf_c(fCase*(x - fOffset)) - fPval;
805 }
806 double fOffset;
807 double fPval;
808 int fCase;
809};
810
811////////////////////////////////////////////////////////////////////////////////
812/// function given the null and the alt p value - return the expected one given the N - sigma value
813
814double AsymptoticCalculator::GetExpectedPValues(double pnull, double palt, double nsigma, bool useCls, bool oneSided ) {
815 if (oneSided) {
819 if (!useCls) return clsplusb;
820 double clb = ROOT::Math::normal_cdf( nsigma, 1.);
821 return (clb == 0) ? -1 : clsplusb / clb;
822 }
823
824 // case of 2 sided test statistic
825 // need to compute numerically
827 if (sqrttmu == 0) {
828 // here cannot invert the function - skip the point
829 return -1;
830 }
831 // invert formula for palt to get sqrttmu_A
832 PaltFunction f( sqrttmu, palt, -1);
835 brf.SetFunction( wf, 0, 20);
836 bool ret = brf.Solve();
837 if (!ret) {
838 oocoutE(nullptr,Eval) << "Error finding expected p-values - return -1" << std::endl;
839 return -1;
840 }
841 double sqrttmu_A = brf.Root();
842
843 // now invert for expected value
846 brf.SetFunction(wf2,0,20);
847 ret = brf.Solve();
848 if (!ret) {
849 oocoutE(nullptr,Eval) << "Error finding expected p-values - return -1" << std::endl;
850 return -1;
851 }
852 return 2*ROOT::Math::normal_cdf_c( brf.Root(),1.);
853}
854
855namespace {
856
857////////////////////////////////////////////////////////////////////////////////
858/// Fill bins by looping recursively on observables.
859
860void FillBins(const RooAbsPdf & pdf, const RooArgList &obs, RooAbsData & data, int &index, double &binVolume, int &ibin) {
861
862 bool debug = (fgPrintLevel() >= 2);
863
864 RooRealVar * v = dynamic_cast<RooRealVar*>(&(obs[index]) );
865 if (!v) return;
866
867 RooArgSet obstmp(obs);
868 double expectedEvents = pdf.expectedEvents(obstmp);
869
870 if (debug) std::cout << "looping on observable " << v->GetName() << std::endl;
871 for (int i = 0; i < v->getBins(); ++i) {
872 v->setBin(i);
873 if (index < int(obs.size()) -1) {
874 index++; // increase index
875 double prevBinVolume = binVolume;
876 binVolume *= v->getBinWidth(i); // increase bin volume
877 FillBins(pdf, obs, data, index, binVolume, ibin);
878 index--; // decrease index
879 binVolume = prevBinVolume; // decrease also bin volume
880 }
881 else {
882
883 // this is now a new bin - compute the pdf in this bin
884 double totBinVolume = binVolume * v->getBinWidth(i);
885 double fval = pdf.getVal(&obstmp)*totBinVolume;
886
887 if (fval*expectedEvents <= 0)
888 {
889 if (fval*expectedEvents < 0) {
890 oocoutW(nullptr,InputArguments)
891 << "AsymptoticCalculator::" << __func__
892 << "(): Bin " << i << " of " << v->GetName() << " has negative expected events! Please check your inputs." << std::endl;
893 }
894 else {
895 oocoutW(nullptr,InputArguments)
896 << "AsymptoticCalculator::" << __func__
897 << "(): Bin " << i << " of " << v->GetName() << " has zero expected events - skip it" << std::endl;
898 }
899 }
900 // have a cut off for overflows ??
901 else {
902 data.add(obs, fval*expectedEvents);
903 }
904
905 if (debug) {
906 std::cout << "bin " << ibin << "\t";
907 for (std::size_t j=0; j < obs.size(); ++j) { std::cout << " " << (static_cast<RooRealVar&>( obs[j])).getVal(); }
908 std::cout << " w = " << fval*expectedEvents;
909 std::cout << std::endl;
910 }
911 ibin++;
912 }
913 }
914 //reset bin values
915 if (debug) {
916 std::cout << "ending loop on .. " << v->GetName() << std::endl;
917 }
918
919 v->setBin(0);
920
921}
922
923bool setObsToExpected(std::span<RooAbsArg *> servers, const RooArgSet &obs, std::string const &errPrefix)
924{
925 RooRealVar *myobs = nullptr;
926 RooAbsReal *myexp = nullptr;
927 for (RooAbsArg *a : servers) {
928 if (obs.contains(*a)) {
929 if (myobs != nullptr) {
930 oocoutF(nullptr,Generation) << errPrefix << "Has two observables ?? " << std::endl;
931 return false;
932 }
933 myobs = dynamic_cast<RooRealVar *>(a);
934 if (myobs == nullptr) {
935 oocoutF(nullptr,Generation) << errPrefix << "Observable is not a RooRealVar??" << std::endl;
936 return false;
937 }
938 } else {
939 if (!a->isConstant() ) {
940 if (myexp != nullptr) {
941 oocoutE(nullptr,Generation) << errPrefix << "Has two non-const arguments " << std::endl;
942 return false;
943 }
944 myexp = dynamic_cast<RooAbsReal *>(a);
945 if (myexp == nullptr) {
946 oocoutF(nullptr,Generation) << errPrefix << "Expected is not a RooAbsReal??" << std::endl;
947 return false;
948 }
949 }
950 }
951 }
952 if (myobs == nullptr) {
953 oocoutF(nullptr,Generation) << errPrefix << "No observable?" << std::endl;
954 return false;
955 }
956 if (myexp == nullptr) {
957 oocoutF(nullptr,Generation) << errPrefix << "No observable?" << std::endl;
958 return false;
959 }
960
961 myobs->setVal(myexp->getVal());
962
963 if (fgPrintLevel() > 2) {
964 std::cout << "SetObsToExpected : setting " << myobs->GetName() << " to expected value " << myexp->getVal() << " of " << myexp->GetName() << std::endl;
965 }
966
967 return true;
968}
969
970////////////////////////////////////////////////////////////////////////////////
971/// set observed value to the expected one
972/// works for Gaussian, Poisson or LogNormal
973/// assumes mean parameter value is the argument not constant and not depending on observables
974/// (if more than two arguments are not constant will use first one but print a warning !)
975/// need to iterate on the components of the Poisson to get n and nu (nu can be a RooAbsReal)
976/// (code from G. Petrucciani and extended by L.M.)
977
978bool SetObsToExpected(RooAbsPdf &pdf, const RooArgSet &obs)
979{
980 std::string const &errPrefix = "AsymptoticCalculator::SetObsExpected( " + std::string{pdf.ClassName()} + " ) : ";
981 std::vector<RooAbsArg *> servers;
982 for (RooAbsArg *a : pdf.servers()) {
983 servers.emplace_back(a);
984 }
985 return setObsToExpected(servers, obs, errPrefix);
986}
987
989{
990 // In the case of the multi-variate Gaussian, we need to iterate over the
991 // dimensions and treat the servers for each dimension separately.
992
993 std::string const &errPrefix = "AsymptoticCalculator::SetObsExpected( " + std::string{mvgauss.ClassName()} + " ) : ";
994 std::vector<RooAbsArg *> servers{nullptr, nullptr};
995 bool ret = true;
996 for (std::size_t iDim = 0; iDim < mvgauss.xVec().size(); ++iDim) {
997 servers[0] = &mvgauss.xVec()[iDim];
998 servers[1] = &mvgauss.muVec()[iDim];
999 ret &= setObsToExpected(servers, obs, errPrefix + " : dim " + std::to_string(iDim) + " ");
1000 }
1001 return ret;
1002}
1003
1004////////////////////////////////////////////////////////////////////////////////
1005/// Inpspect a product pdf to find all the Poisson or Gaussian parts to set the observed
1006/// values to expected ones.
1007
1008bool setObsToExpectedProdPdf(RooProdPdf &prod, const RooArgSet &obs)
1009{
1010 bool ret = true;
1011 for (auto *a : prod.pdfList()) {
1012 if (!a->dependsOn(obs)) continue;
1013 RooPoisson *pois = nullptr;
1014 RooGaussian *gauss = nullptr;
1015 RooMultiVarGaussian *mvgauss = nullptr;
1016 // should try to add also lognormal case ?
1017 if ((pois = dynamic_cast<RooPoisson *>(a)) != nullptr) {
1018 ret &= SetObsToExpected(*pois, obs);
1019 pois->setNoRounding(true); //needed since expected value is not an integer
1020 } else if ((gauss = dynamic_cast<RooGaussian *>(a)) != nullptr) {
1021 ret &= SetObsToExpected(*gauss, obs);
1022 } else if ((mvgauss = dynamic_cast<RooMultiVarGaussian *>(a)) != nullptr) {
1024 } else if (RooProdPdf *subprod = dynamic_cast<RooProdPdf *>(a)) {
1026 } else {
1027 oocoutE(nullptr, InputArguments)
1028 << "Illegal term in counting model: "
1029 << "the PDF " << a->GetName() << " depends on the observables, but is not a Poisson, Gaussian or Product"
1030 << std::endl;
1031 return false;
1032 }
1033 }
1034
1035 return ret;
1036}
1037
1038////////////////////////////////////////////////////////////////////////////////
1039/// Generate counting Asimov data for the case when the pdf cannot be extended.
1040/// This function assumes that the pdf is a RooPoisson or can be decomposed in a product of RooPoisson,
1041/// or is a RooGaussian. Otherwise, we cannot know how to make the Asimov data sets.
1042
1044 RooArgSet obs(observables);
1045 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
1046 RooPoisson *pois = nullptr;
1047 RooGaussian *gauss = nullptr;
1048 RooMultiVarGaussian *mvgauss = nullptr;
1049
1050 if (fgPrintLevel() > 1)
1051 std::cout << "generate counting Asimov data for pdf of type " << pdf.ClassName() << std::endl;
1052
1053 bool r = false;
1054 if (prod != nullptr) {
1055 r = setObsToExpectedProdPdf(*prod, observables);
1056 } else if ((pois = dynamic_cast<RooPoisson *>(&pdf)) != nullptr) {
1057 r = SetObsToExpected(*pois, observables);
1058 // we need in this case to set Poisson to real values
1059 pois->setNoRounding(true);
1060 } else if ((gauss = dynamic_cast<RooGaussian *>(&pdf)) != nullptr) {
1061 r = SetObsToExpected(*gauss, observables);
1062 } else if ((mvgauss = dynamic_cast<RooMultiVarGaussian *>(&pdf)) != nullptr) {
1063 r = setObsToExpectedMultiVarGauss(*mvgauss, observables);
1064 } else {
1065 oocoutE(nullptr,InputArguments) << "A counting model pdf must be either a RooProdPdf or a RooPoisson or a RooGaussian" << std::endl;
1066 }
1067 if (!r) return nullptr;
1068 int icat = 0;
1069 if (channelCat) {
1070 icat = channelCat->getCurrentIndex();
1071 }
1072
1073 RooDataSet *ret = new RooDataSet("CountingAsimovData" + std::to_string(icat),
1074 "CountingAsimovData" + std::to_string(icat), obs);
1075 ret->add(obs);
1076 return ret;
1077}
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// Compute the asimov data set for an observable of a pdf.
1081/// It generates binned data following the binning of the observables.
1082// TODO: (possibility to change number of bins)
1083// TODO: implement integration over bin content
1084
1086
1087 int printLevel = fgPrintLevel();
1088
1089 // Get observables defined by the pdf associated with this state
1090 std::unique_ptr<RooArgSet> obs(pdf.getObservables(allobs) );
1091
1092
1093 // if pdf cannot be extended assume is then a counting experiment
1094 if (!pdf.canBeExtended() ) return GenerateCountingAsimovData(const_cast<RooAbsPdf&>(pdf), *obs, weightVar, channelCat);
1095
1096 RooArgSet obsAndWeight(*obs);
1097 obsAndWeight.add(weightVar);
1098
1099 std::unique_ptr<RooDataSet> asimovData;
1100 if (channelCat) {
1101 int icat = channelCat->getCurrentIndex();
1102 asimovData = std::make_unique<RooDataSet>("AsimovData" + std::to_string(icat),
1103 "combAsimovData" + std::to_string(icat),
1105 }
1106 else {
1107 asimovData = std::make_unique<RooDataSet>("AsimovData","AsimovData",RooArgSet(obsAndWeight),RooFit::WeightVar(weightVar));
1108 }
1109
1110 // This works only for 1D observables
1111 //RooRealVar* thisObs = ((RooRealVar*)obstmp->first());
1112
1113 RooArgList obsList(*obs);
1114
1115 // loop on observables and on the bins
1116 if (printLevel >= 2) {
1117 std::cout << "Generating Asimov data for pdf " << pdf.GetName() << std::endl;
1118 std::cout << "list of observables " << std::endl;
1119 obsList.Print();
1120 }
1121
1122 int obsIndex = 0;
1123 double binVolume = 1;
1124 int nbins = 0;
1125 FillBins(pdf, obsList, *asimovData, obsIndex, binVolume, nbins);
1126 if (printLevel >= 2)
1127 std::cout << "filled from " << pdf.GetName() << " " << nbins << " nbins " << " volume is " << binVolume << std::endl;
1128
1129 // for (int iobs = 0; iobs < obsList.size(); ++iobs) {
1130 // RooRealVar * thisObs = dynamic_cast<RooRealVar*> &obsList[i];
1131 // if (thisObs == 0) continue;
1132 // // loop on the bin contents
1133 // for(int ibin=0; ibin<thisObs->numBins(); ++ibin){
1134 // thisObs->setBin(ibin);
1135
1136 // thisNorm=pdftmp->getVal(obstmp)*thisObs->getBinWidth(jj);
1137 // if (thisNorm*expectedEvents <= 0)
1138 // {
1139 // std::cout << "WARNING::Detected bin with zero expected events! Please check your inputs." << std::endl;
1140 // }
1141 // // have a cut off for overflows ??
1142 // obsDataUnbinned->add(*mc->GetObservables(), thisNorm*expectedEvents);
1143 // }
1144
1145 if (printLevel >= 1)
1146 {
1147 asimovData->Print();
1148 }
1149 if( TMath::IsNaN(asimovData->sumEntries()) ){
1150 std::cout << "sum entries is nan"<< std::endl;
1151 assert(0);
1152 asimovData = nullptr;
1153 }
1154
1155 return asimovData.release();
1156
1157}
1158
1159} // namespace
1160
1161////////////////////////////////////////////////////////////////////////////////
1162/// generate the asimov data for the observables (not the global ones)
1163/// need to deal with the case of a sim pdf
1164
1166
1167 int printLevel = fgPrintLevel();
1168
1169 RooRealVar weightVar{"binWeightAsimov", "binWeightAsimov", 1, 0, 1.e30};
1170
1171 if (printLevel > 1) std::cout <<" Generate Asimov data for observables"<< std::endl;
1172 //RooDataSet* simData=nullptr;
1173 const RooSimultaneous* simPdf = dynamic_cast<const RooSimultaneous*>(&pdf);
1174 if (!simPdf) {
1175 // generate data for non sim pdf
1176 return GenerateAsimovDataSinglePdf( pdf, observables, weightVar, nullptr);
1177 }
1178
1179 std::map<std::string, std::unique_ptr<RooDataSet>> asimovDataMap;
1180
1181 //look at category of simpdf
1182 RooCategory& channelCat = const_cast<RooCategory&>(dynamic_cast<const RooCategory&>(simPdf->indexCat()));
1183 int nrIndices = channelCat.numTypes();
1184 if( nrIndices == 0 ) {
1185 oocoutW(nullptr,Generation) << "Simultaneous pdf does not contain any categories." << std::endl;
1186 }
1187 for (int i=0;i<nrIndices;i++){
1188 channelCat.setIndex(i);
1189 //iFrame++;
1190 // Get pdf associated with state from simpdf
1191 RooAbsPdf* pdftmp = simPdf->getPdf(channelCat.getCurrentLabel()) ;
1192 assert(pdftmp != nullptr);
1193
1194 if (printLevel > 1)
1195 {
1196 std::cout << "on type " << channelCat.getCurrentLabel() << " " << channelCat.getCurrentIndex() << std::endl;
1197 }
1198
1199 std::unique_ptr<RooDataSet> dataSinglePdf{static_cast<RooDataSet*>(GenerateAsimovDataSinglePdf( *pdftmp, observables, weightVar, &channelCat))};
1200 if (!dataSinglePdf) {
1201 oocoutE(nullptr,Generation) << "Error generating an Asimov data set for pdf " << pdftmp->GetName() << std::endl;
1202 return nullptr;
1203 }
1204
1205 if (asimovDataMap.count(string(channelCat.getCurrentLabel())) != 0) {
1206 oocoutE(nullptr,Generation) << "AsymptoticCalculator::GenerateAsimovData(): The PDF for " << channelCat.getCurrentLabel()
1207 << " was already defined. It will be overridden. The faulty category definitions follow:" << std::endl;
1208 channelCat.Print("V");
1209 }
1210
1211 if (printLevel > 1)
1212 {
1213 std::cout << "channel: " << channelCat.getCurrentLabel() << ", data: ";
1214 dataSinglePdf->Print();
1215 std::cout << std::endl;
1216 }
1217
1218 asimovDataMap[string(channelCat.getCurrentLabel())] = std::move(dataSinglePdf);
1219 }
1220
1221 RooArgSet obsAndWeight(observables);
1222 obsAndWeight.add(weightVar);
1223
1224
1225 return new RooDataSet("asimovDataFullModel","asimovDataFullModel",RooArgSet(obsAndWeight,channelCat),
1227}
1228
1229////////////////////////////////////////////////////////////////////////////////
1230/// Make the Asimov data from the ModelConfig and list of poi
1231/// \param realData Real data
1232/// \param model Model config defining the pdf and the parameters
1233/// \param paramValues The snapshot of POI and parameters used for finding the best nuisance parameter values (conditioned at these values)
1234/// \param[out] asimovGlobObs Global observables set to values satisfying the constraints
1235/// \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
1236/// given an observed data set, a model and a snapshot of the poi.
1237/// \return The asimov data set. The user takes ownership.
1238///
1239
1241
1242 int verbose = fgPrintLevel();
1243
1244
1245 RooArgSet poi(*model.GetParametersOfInterest());
1246 poi.assign(paramValues);
1247
1248 // set poi constant for conditional MLE
1249 // need to fit nuisance parameters at their conditional MLE value
1251 for (auto *tmpPar : static_range_cast<RooRealVar *>(poi)) {
1252 tmpPar->setConstant();
1253 if (verbose>0)
1254 std::cout << "MakeAsimov: Setting poi " << tmpPar->GetName() << " to a constant value = " << tmpPar->getVal() << std::endl;
1256 }
1257
1258 // find conditional value of the nuisance parameters
1259 bool hasFloatParams = false;
1261 if (model.GetNuisanceParameters()) {
1264 if (!constrainParams.empty()) hasFloatParams = true;
1265
1266 } else {
1267 // Do we have free parameters anyway that need fitting?
1268 std::unique_ptr<RooArgSet> params(model.GetPdf()->getParameters(realData));
1269 for (auto const *rrv : dynamic_range_cast<RooRealVar *>(*params)) {
1270 if ( rrv != nullptr && rrv->isConstant() == false ) { hasFloatParams = true; break; }
1271 }
1272 }
1273 if (hasFloatParams) {
1274 // models need to be fitted to find best nuisance parameter values
1275
1276 TStopwatch tw2; tw2.Start();
1278 if (verbose>0) {
1279 std::cout << "MakeAsimov: doing a conditional fit for finding best nuisance values " << std::endl;
1280 minimPrintLevel = verbose;
1281 if (verbose>1) {
1282 std::cout << "POI values:\n"; poi.Print("v");
1283 if (verbose > 2) {
1284 std::cout << "Nuis param values:\n";
1285 constrainParams.Print("v");
1286 }
1287 }
1288 }
1290 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
1291
1293 std::vector<RooCmdArg> args{
1294 RooFit::Minimizer("",minimizerAlgo.c_str()), // empty mimimizer type to select default
1297 RooFit::Hesse(false),
1299 RooFit::Offset(GetGlobalRooStatsConfig().useLikelihoodOffset),
1301 };
1302
1303 RooLinkedList argList;
1304 for (auto& arg : args) {
1305 argList.Add(&arg);
1306 }
1307 model.fitTo(realData, argList);
1308 if (verbose>0) { std::cout << "fit time "; tw2.Print();}
1309 if (verbose > 1) {
1310 // after the fit the nuisance parameters will have their best fit value
1311 if (model.GetNuisanceParameters() ) {
1312 std::cout << "Nuisance parameters after fit for asimov dataset: " << std::endl;
1313 model.GetNuisanceParameters()->Print("V");
1314 }
1315 }
1316
1317 if (verbose < 2) RooMsgService::instance().setGlobalKillBelow(msglevel);
1318
1319 }
1320
1321 // restore the parameters which were set constant
1323
1324 std::unique_ptr<RooArgSet> allParams{model.GetPdf()->getParameters(realData)};
1325
1327
1328 // if a RooArgSet of poi is passed , different poi will be used for generating the Asimov data set
1329 if (genPoiValues) {
1330 allParams->assign(*genPoiValues);
1331 }
1332
1333 // now do the actual generation of the AsimovData Set
1334 // no need to pass parameters values since we have set them before
1335 return MakeAsimovData(model, *allParams, asimovGlobObs);
1336}
1337
1338////////////////////////////////////////////////////////////////////////////////
1339/// \param model ModelConfig that contains the model pdf and the model parameters
1340/// \param allParamValues The parameters of the model will be set to the values given in this set
1341/// \param[out] asimovGlobObs Global observables set to values satisfying the constraints
1342/// \return Asimov data set. The user takes ownership.
1343///
1344/// The parameter values (including the nuisance parameter) can result from a fit to data or be at the nominal values.
1345///
1346
1348
1349 int verbose = fgPrintLevel();
1350
1351 TStopwatch tw;
1352 tw.Start();
1353
1354 // set the parameter values (do I need the poi to be constant ? )
1355 // the nuisance parameter values could be set at their fitted value (the MLE)
1356 if (!allParamValues.empty()) {
1357 std::unique_ptr<RooArgSet> allVars{model.GetPdf()->getVariables()};
1358 allVars->assign(allParamValues);
1359 }
1360
1361
1362 // generate the Asimov data set for the observables
1363 RooAbsData * asimov = GenerateAsimovData(*model.GetPdf() , *model.GetObservables() );
1364
1365 if (verbose>0) {
1366 std::cout << "Generated Asimov data for observables "; (model.GetObservables() )->Print();
1367 if (verbose > 1) {
1368 if (asimov->numEntries() == 1 ) {
1369 std::cout << "--- Asimov data values \n";
1370 asimov->get()->Print("v");
1371 }
1372 else {
1373 std::cout << "--- Asimov data numEntries = " << asimov->numEntries() << " sumOfEntries = " << asimov->sumEntries() << std::endl;
1374 }
1375 std::cout << "\ttime for generating : "; tw.Print();
1376 }
1377 }
1378
1379
1380 // Now need to have in ASIMOV the data sets also the global observables
1381 // Their values must be the one satisfying the constraint.
1382 // to do it make a nuisance pdf with all product of constraints and then
1383 // assign to each constraint a glob observable value = to the current fitted nuisance parameter value
1384 // IN general one should solve in general the system of equations f( gobs| nuispar ) = 0 where f are the
1385 // derivatives of the constraint with respect the nuisance parameter and they are evaluated at the best fit nuisance
1386 // parameter points
1387 // As simple solution assume that constrain has a direct dependence on the nuisance parameter, i.e.
1388 // Constraint (gobs, func( nuispar) ) and the condition is satisfied for
1389 // gobs = func( nuispar) where nunispar is at the MLE value
1390
1391
1392 if (model.GetGlobalObservables() && !model.GetGlobalObservables()->empty()) {
1393
1394 if (verbose>1) {
1395 std::cout << "Generating Asimov data for global observables " << std::endl;
1396 }
1397
1399
1400 // snapshot data global observables
1402 SetAllConstant(gobs, true);
1403 gobs.snapshot(snapGlobalObsData);
1404
1405
1407 if (model.GetNuisanceParameters()) nuis.add(*model.GetNuisanceParameters());
1408 if (nuis.empty()) {
1409 oocoutW(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData: model does not have nuisance parameters but has global observables"
1410 << " set global observables to model values " << std::endl;
1411 asimovGlobObs.assign(gobs);
1412 return asimov;
1413 }
1414
1415 // part 1: create the nuisance pdf
1416 std::unique_ptr<RooAbsPdf> nuispdf(RooStats::MakeNuisancePdf(model,"TempNuisPdf") );
1417 if (nuispdf == nullptr) {
1418 oocoutF(nullptr, Generation) << "AsymptoticCalculator::MakeAsimovData: model has nuisance parameters and "
1419 "global obs but no nuisance pdf "
1420 << std::endl;
1421 }
1422 // unfold the nuisance pdf if it is a prod pdf
1423 RooArgList pdfList;
1424 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(nuispdf.get());
1425 if (prod ) {
1426 pdfList.add(prod->pdfList());
1427 } else {
1428 // nothing to unfold - just use the pdf
1429 pdfList.add(*nuispdf);
1430 }
1431
1432 for (auto *cterm : static_range_cast<RooAbsPdf *>(pdfList)) {
1433 assert(dynamic_cast<RooAbsPdf *>(static_cast<RooAbsArg *>(cterm)) &&
1434 "AsimovUtils: a factor of the nuisance pdf is not a Pdf!");
1435
1436 if (!cterm->dependsOn(nuis)) continue; // dummy constraints
1437 // skip also the case of uniform components
1438 if (typeid(*cterm) == typeid(RooUniform)) continue;
1439
1440 std::unique_ptr<RooArgSet> cpars(cterm->getParameters(&gobs));
1441 std::unique_ptr<RooArgSet> cgobs(cterm->getObservables(&gobs));
1442 if (cgobs->size() > 1) {
1443 oocoutE(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1444 << " has multiple global observables -cannot generate - skip it" << std::endl;
1445 continue;
1446 }
1447 else if (cgobs->empty()) {
1448 oocoutW(nullptr, Generation)
1449 << "AsymptoticCalculator::MakeAsimovData: constraint term " << cterm->GetName()
1450 << " has no global observables - skip it" << std::endl;
1451 continue;
1452 }
1453 // the variable representing the global observable
1454 RooRealVar &rrv = dynamic_cast<RooRealVar &>(*cgobs->first());
1455
1456 // remove the constant parameters in cpars
1458 if (cpars->size() != 1) {
1459 oocoutE(nullptr, Generation)
1460 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1461 << cterm->GetName() << " has multiple floating params - cannot generate - skip it " << std::endl;
1462 continue;
1463 }
1464
1465 bool foundServer = false;
1466 // note : this will work only for this type of constraints
1467 // expressed as RooPoisson, RooGaussian, RooLognormal, RooGamma
1468 TClass * cClass = cterm->IsA();
1469 if (verbose > 2) std::cout << "Constraint " << cterm->GetName() << " of type " << cClass->GetName() << std::endl;
1473 TString className = (cClass) ? cClass->GetName() : "undefined";
1474 oocoutW(nullptr, Generation)
1475 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1476 << cterm->GetName() << " of type " << className
1477 << " is a non-supported type - result might be not correct " << std::endl;
1478 }
1479
1480 // in case of a Poisson constraint make sure the rounding is not set
1481 if (cClass == RooPoisson::Class() ) {
1482 RooPoisson * pois = static_cast<RooPoisson*>(cterm);
1483 assert(dynamic_cast<RooPoisson *>(cterm));
1484 pois->setNoRounding(true);
1485 }
1486
1487 // look at server of the constraint term and check if the global observable is part of the server
1488 RooAbsArg * arg = cterm->findServer(rrv);
1489 if (!arg) {
1490 // special case is for the Gamma where one might define the global observable n and you have a Gamma(b, n+1, ...._
1491 // 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
1492 // so in case of the Gamma ignore this test
1493 if ( cClass != RooGamma::Class() ) {
1494 oocoutE(nullptr, Generation)
1495 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1496 << cterm->GetName() << " has no direct dependence on global observable- cannot generate it " << std::endl;
1497 continue;
1498 }
1499 }
1500
1501 // loop on the server of the constraint term
1502 // need to treat the Gamma as a special case
1503 // the mode of the Gamma is (k-1)*theta where theta is the inverse of the rate parameter.
1504 // 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
1505 RooAbsReal * thetaGamma = nullptr;
1506 if ( cClass == RooGamma::Class() ) {
1507 for (RooAbsArg *a2 : cterm->servers()) {
1508 if (TString(a2->GetName()).Contains("theta") ) {
1509 thetaGamma = dynamic_cast<RooAbsReal*>(a2);
1510 break;
1511 }
1512 }
1513 if (thetaGamma == nullptr) {
1514 oocoutI(nullptr, Generation)
1515 << "AsymptoticCalculator::MakeAsimovData:constraint term "
1516 << cterm->GetName() << " is a Gamma distribution and no server named theta is found. Assume that the Gamma scale is 1 " << std::endl;
1517 }
1518 else if (verbose>2) {
1519 std::cout << "Gamma constraint has a scale " << thetaGamma->GetName() << " = " << thetaGamma->getVal() << std::endl;
1520 }
1521 }
1522 for (RooAbsArg *a2 : cterm->servers()) {
1523 RooAbsReal * rrv2 = dynamic_cast<RooAbsReal *>(a2);
1524 if (verbose > 2) std::cout << "Loop on constraint server term " << a2->GetName() << std::endl;
1525 if (rrv2 && rrv2->dependsOn(nuis) ) {
1526
1527
1528 // found server depending on nuisance
1529 if (foundServer) {
1530 oocoutE(nullptr,Generation) << "AsymptoticCalculator::MakeAsimovData:constraint term "
1531 << cterm->GetName() << " constraint term has more server depending on nuisance- cannot generate it " <<
1532 std::endl;
1533 foundServer = false;
1534 break;
1535 }
1536 if (thetaGamma && thetaGamma->getVal() > 0) {
1537 rrv.setVal( rrv2->getVal() / thetaGamma->getVal() );
1538 } else {
1539 rrv.setVal(rrv2->getVal());
1540 }
1541 foundServer = true;
1542
1543 if (verbose > 2) {
1544 std::cout << "setting global observable " << rrv.GetName() << " to value " << rrv.getVal()
1545 << " which comes from " << rrv2->GetName() << std::endl;
1546 }
1547 }
1548 }
1549
1550 if (!foundServer) {
1551 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;
1552 std::cerr << "Parameters: " << std::endl;
1553 cpars->Print("V");
1554 std::cerr << "Observables: " << std::endl;
1555 cgobs->Print("V");
1556 }
1557 }
1558
1559 // make a snapshot of global observables
1560 // needed this ?? (LM)
1561
1562 asimovGlobObs.removeAll();
1563 SetAllConstant(gobs, true);
1564 gobs.snapshot(asimovGlobObs);
1565
1566 // revert global observables to the data value
1567 gobs.assign(snapGlobalObsData);
1568
1569 if (verbose>0) {
1570 std::cout << "Generated Asimov data for global observables ";
1571 if (verbose == 1) gobs.Print();
1572 }
1573
1574 if (verbose > 1) {
1575 std::cout << "\nGlobal observables for data: " << std::endl;
1576 gobs.Print("V");
1577 std::cout << "\nGlobal observables for asimov: " << std::endl;
1578 asimovGlobObs.Print("V");
1579 }
1580
1581
1582 }
1583
1584 return asimov;
1585
1586}
#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 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
void Print(GNN_Data &d, std::string txt="")
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:77
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)
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
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:57
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:101
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:40
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:218
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
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:34
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)
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
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:226
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:139
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
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:379
MsgLevel
Verbosity level for RooMsgService::StreamConfig in RooMsgService.
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
std::string const & NLLOffsetMode()
Test what offsetting mode RooStats should use by default.
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
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:896
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:906
PaltFunction(double offset, double pval, int icase)