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